1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985
|
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Application
createBaffles
Description
Makes internal faces into boundary faces. Does not duplicate points, unlike
mergeOrSplitBaffles.
Note: if any coupled patch face is selected for baffling the opposite
member has to be selected for baffling as well.
- if the patch already exists will not override it nor its fields
- if the patch does not exist it will be created together with 'calculated'
patchfields unless the field is mentioned in the patchFields section.
- any 0-sized patches (since faces have been moved out) will get removed
\*---------------------------------------------------------------------------*/
#include "argList.H"
#include "Time.H"
#include "polyTopoChange.H"
#include "polyModifyFace.H"
#include "polyAddFace.H"
#include "ReadFields.H"
#include "volFields.H"
#include "surfaceFields.H"
#include "fvMeshMapper.H"
#include "faceSelection.H"
#include "fvMeshTools.H"
using namespace Foam;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
label addPatch
(
fvMesh& mesh,
const word& patchName,
const word& groupName,
const dictionary& patchDict
)
{
const polyBoundaryMesh& pbm = mesh.boundaryMesh();
if (pbm.findPatchID(patchName) == -1)
{
autoPtr<polyPatch> ppPtr
(
polyPatch::New
(
patchName,
patchDict,
0,
pbm
)
);
polyPatch& pp = ppPtr();
if (!groupName.empty() && !pp.inGroup(groupName))
{
pp.inGroups().append(groupName);
}
// Add patch, create calculated everywhere
fvMeshTools::addPatch
(
mesh,
pp,
dictionary(), // do not set specialised patchFields
calculatedFvPatchField<scalar>::typeName,
true // parallel sync'ed addition
);
}
else
{
Info<< "Patch '" << patchName
<< "' already exists. Only "
<< "moving patch faces - type will remain the same"
<< endl;
}
return pbm.findPatchID(patchName);
}
// Filter out the empty patches.
void filterPatches(fvMesh& mesh, const HashSet<word>& addedPatchNames)
{
// Remove any zero-sized ones. Assumes
// - processor patches are already only there if needed
// - all other patches are available on all processors
// - but coupled ones might still be needed, even if zero-size
// (e.g. processorCyclic)
// See also logic in createPatch.
const polyBoundaryMesh& pbm = mesh.boundaryMesh();
labelList oldToNew(pbm.size(), -1);
label newPatchi = 0;
forAll(pbm, patchi)
{
const polyPatch& pp = pbm[patchi];
if (!isA<processorPolyPatch>(pp))
{
if
(
isA<coupledPolyPatch>(pp)
|| returnReduce(pp.size(), sumOp<label>())
|| addedPatchNames.found(pp.name())
)
{
// Coupled (and unknown size) or uncoupled and used
oldToNew[patchi] = newPatchi++;
}
}
}
forAll(pbm, patchi)
{
const polyPatch& pp = pbm[patchi];
if (isA<processorPolyPatch>(pp))
{
oldToNew[patchi] = newPatchi++;
}
}
const label nKeepPatches = newPatchi;
// Shuffle unused ones to end
if (nKeepPatches != pbm.size())
{
Info<< endl
<< "Removing zero-sized patches:" << endl << incrIndent;
forAll(oldToNew, patchi)
{
if (oldToNew[patchi] == -1)
{
Info<< indent << pbm[patchi].name()
<< " type " << pbm[patchi].type()
<< " at position " << patchi << endl;
oldToNew[patchi] = newPatchi++;
}
}
Info<< decrIndent;
fvMeshTools::reorderPatches(mesh, oldToNew, nKeepPatches, true);
Info<< endl;
}
}
void modifyOrAddFace
(
polyTopoChange& meshMod,
const face& f,
const label facei,
const label own,
const bool flipFaceFlux,
const label newPatchi,
const label zoneID,
const bool zoneFlip,
PackedBoolList& modifiedFace
)
{
if (!modifiedFace[facei])
{
// First usage of face. Modify.
meshMod.setAction
(
polyModifyFace
(
f, // modified face
facei, // label of face
own, // owner
-1, // neighbour
flipFaceFlux, // face flip
newPatchi, // patch for face
false, // remove from zone
zoneID, // zone for face
zoneFlip // face flip in zone
)
);
modifiedFace[facei] = 1;
}
else
{
// Second or more usage of face. Add.
meshMod.setAction
(
polyAddFace
(
f, // modified face
own, // owner
-1, // neighbour
-1, // master point
-1, // master edge
facei, // master face
flipFaceFlux, // face flip
newPatchi, // patch for face
zoneID, // zone for face
zoneFlip // face flip in zone
)
);
}
}
// Create faces for fZone faces. Usually newMasterPatches, newSlavePatches
// only size one but can be more for duplicate baffle sets
void createFaces
(
const bool internalFacesOnly,
const fvMesh& mesh,
const faceZone& fZone,
const labelList& newMasterPatches,
const labelList& newSlavePatches,
polyTopoChange& meshMod,
PackedBoolList& modifiedFace,
label& nModified
)
{
const polyBoundaryMesh& pbm = mesh.boundaryMesh();
forAll(newMasterPatches, i)
{
// Pass 1. Do selected side of zone
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
{
label zoneFacei = fZone.whichFace(facei);
if (zoneFacei != -1)
{
if (!fZone.flipMap()[zoneFacei])
{
// Use owner side of face
modifyOrAddFace
(
meshMod,
mesh.faces()[facei], // modified face
facei, // label of face
mesh.faceOwner()[facei],// owner
false, // face flip
newMasterPatches[i], // patch for face
fZone.index(), // zone for face
false, // face flip in zone
modifiedFace // modify or add status
);
}
else
{
// Use neighbour side of face.
// To keep faceZone pointing out of original neighbour
// we don't need to set faceFlip since that cell
// now becomes the owner
modifyOrAddFace
(
meshMod,
mesh.faces()[facei].reverseFace(), // modified face
facei, // label of face
mesh.faceNeighbour()[facei],// owner
true, // face flip
newMasterPatches[i], // patch for face
fZone.index(), // zone for face
false, // face flip in zone
modifiedFace // modify or add status
);
}
nModified++;
}
}
// Pass 2. Do other side of zone
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
{
label zoneFacei = fZone.whichFace(facei);
if (zoneFacei != -1)
{
if (!fZone.flipMap()[zoneFacei])
{
// Use neighbour side of face
modifyOrAddFace
(
meshMod,
mesh.faces()[facei].reverseFace(), // modified face
facei, // label of face
mesh.faceNeighbour()[facei], // owner
true, // face flip
newSlavePatches[i], // patch for face
fZone.index(), // zone for face
true, // face flip in zone
modifiedFace // modify or add
);
}
else
{
// Use owner side of face
modifyOrAddFace
(
meshMod,
mesh.faces()[facei], // modified face
facei, // label of face
mesh.faceOwner()[facei],// owner
false, // face flip
newSlavePatches[i], // patch for face
fZone.index(), // zone for face
true, // face flip in zone
modifiedFace // modify or add status
);
}
}
}
// Modify any boundary faces
// ~~~~~~~~~~~~~~~~~~~~~~~~~
// Normal boundary:
// - move to new patch. Might already be back-to-back baffle
// you want to add cyclic to. Do warn though.
//
// Processor boundary:
// - do not move to cyclic
// - add normal patches though.
// For warning once per patch.
labelHashSet patchWarned;
forAll(pbm, patchi)
{
const polyPatch& pp = pbm[patchi];
label newPatchi = newMasterPatches[i];
if (pp.coupled() && pbm[newPatchi].coupled())
{
// Do not allow coupled faces to be moved to different
// coupled patches.
}
else if (pp.coupled() || !internalFacesOnly)
{
forAll(pp, i)
{
label facei = pp.start()+i;
label zoneFacei = fZone.whichFace(facei);
if (zoneFacei != -1)
{
if (patchWarned.insert(patchi))
{
WarningInFunction
<< "Found boundary face (in patch "
<< pp.name()
<< ") in faceZone " << fZone.name()
<< " to convert to baffle patch "
<< pbm[newPatchi].name()
<< endl
<< " Run with -internalFacesOnly option"
<< " if you don't wish to convert"
<< " boundary faces." << endl;
}
modifyOrAddFace
(
meshMod,
mesh.faces()[facei], // modified face
facei, // label of face
mesh.faceOwner()[facei], // owner
false, // face flip
newPatchi, // patch for face
fZone.index(), // zone for face
fZone.flipMap()[zoneFacei], // face flip in zone
modifiedFace // modify or add
);
nModified++;
}
}
}
}
}
}
int main(int argc, char *argv[])
{
argList::addNote
(
"Makes internal faces into boundary faces.\n"
"Does not duplicate points."
);
#include "addDictOption.H"
#include "addOverwriteOption.H"
#include "addDictOption.H"
#include "addRegionOption.H"
#include "setRootCase.H"
#include "createTime.H"
runTime.functionObjects().off();
#include "createNamedMesh.H"
const bool overwrite = args.optionFound("overwrite");
const word oldInstance = mesh.pointsInstance();
const word dictName("createBafflesDict");
#include "setSystemMeshDictionaryIO.H"
Switch internalFacesOnly(false);
Switch noFields(false);
PtrList<faceSelection> selectors;
{
Info<< "Reading baffle criteria from " << dictName << nl << endl;
IOdictionary dict(dictIO);
dict.lookup("internalFacesOnly") >> internalFacesOnly;
noFields = dict.lookupOrDefault("noFields", false);
const dictionary& selectionsDict = dict.subDict("baffles");
label n = 0;
forAllConstIter(dictionary, selectionsDict, iter)
{
if (iter().isDict())
{
n++;
}
}
selectors.setSize(n);
n = 0;
forAllConstIter(dictionary, selectionsDict, iter)
{
if (iter().isDict())
{
selectors.set
(
n++,
faceSelection::New(iter().keyword(), mesh, iter().dict())
);
}
}
}
if (internalFacesOnly)
{
Info<< "Not converting faces on non-coupled patches." << nl << endl;
}
// Read objects in time directory
IOobjectList objects(mesh, runTime.timeName());
// Read vol fields.
Info<< "Reading geometric fields" << nl << endl;
PtrList<volScalarField> vsFlds;
if (!noFields) ReadFields(mesh, objects, vsFlds);
PtrList<volVectorField> vvFlds;
if (!noFields) ReadFields(mesh, objects, vvFlds);
PtrList<volSphericalTensorField> vstFlds;
if (!noFields) ReadFields(mesh, objects, vstFlds);
PtrList<volSymmTensorField> vsymtFlds;
if (!noFields) ReadFields(mesh, objects, vsymtFlds);
PtrList<volTensorField> vtFlds;
if (!noFields) ReadFields(mesh, objects, vtFlds);
// Read surface fields.
PtrList<surfaceScalarField> ssFlds;
if (!noFields) ReadFields(mesh, objects, ssFlds);
PtrList<surfaceVectorField> svFlds;
if (!noFields) ReadFields(mesh, objects, svFlds);
PtrList<surfaceSphericalTensorField> sstFlds;
if (!noFields) ReadFields(mesh, objects, sstFlds);
PtrList<surfaceSymmTensorField> ssymtFlds;
if (!noFields) ReadFields(mesh, objects, ssymtFlds);
PtrList<surfaceTensorField> stFlds;
if (!noFields) ReadFields(mesh, objects, stFlds);
// Creating (if necessary) faceZones
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
forAll(selectors, selectorI)
{
const word& name = selectors[selectorI].name();
if (mesh.faceZones().findZoneID(name) == -1)
{
mesh.faceZones().clearAddressing();
label sz = mesh.faceZones().size();
labelList addr(0);
boolList flip(0);
mesh.faceZones().setSize(sz+1);
mesh.faceZones().set
(
sz,
new faceZone(name, addr, flip, sz, mesh.faceZones())
);
}
}
// Select faces
// ~~~~~~~~~~~~
//- Per face zoneID it is in and flip status.
labelList faceToZoneID(mesh.nFaces(), -1);
boolList faceToFlip(mesh.nFaces(), false);
forAll(selectors, selectorI)
{
const word& name = selectors[selectorI].name();
label zoneID = mesh.faceZones().findZoneID(name);
selectors[selectorI].select(zoneID, faceToZoneID, faceToFlip);
}
// Add faces to faceZones
labelList nFaces(mesh.faceZones().size(), 0);
forAll(faceToZoneID, facei)
{
label zoneID = faceToZoneID[facei];
if (zoneID != -1)
{
nFaces[zoneID]++;
}
}
forAll(selectors, selectorI)
{
const word& name = selectors[selectorI].name();
label zoneID = mesh.faceZones().findZoneID(name);
label& n = nFaces[zoneID];
labelList addr(n);
boolList flip(n);
n = 0;
forAll(faceToZoneID, facei)
{
label zone = faceToZoneID[facei];
if (zone == zoneID)
{
addr[n] = facei;
flip[n] = faceToFlip[facei];
n++;
}
}
Info<< "Created zone " << name
<< " at index " << zoneID
<< " with " << n << " faces" << endl;
mesh.faceZones().set
(
zoneID,
new faceZone(name, addr, flip, zoneID, mesh.faceZones())
);
}
// Count patches to add
// ~~~~~~~~~~~~~~~~~~~~
HashSet<word> bafflePatches;
{
forAll(selectors, selectorI)
{
const dictionary& dict = selectors[selectorI].dict();
if (dict.found("patches"))
{
const dictionary& patchSources = dict.subDict("patches");
forAllConstIter(dictionary, patchSources, iter)
{
const word patchName(iter().dict()["name"]);
bafflePatches.insert(patchName);
}
}
else
{
const word masterName = selectors[selectorI].name() + "_master";
bafflePatches.insert(masterName);
const word slaveName = selectors[selectorI].name() + "_slave";
bafflePatches.insert(slaveName);
}
}
}
// Create baffles
// ~~~~~~~~~~~~~~
// Is done in multiple steps
// - create patches with 'calculated' patchFields
// - move faces into these patches
// - change the patchFields to the wanted type
// This order is done so e.g. fixedJump works:
// - you cannot create patchfields at the same time as patches since
// they do an evaluate upon construction
// - you want to create the patchField only after you have faces
// so you don't get the 'create-from-nothing' mapping problem.
// Pass 1: add patches
// ~~~~~~~~~~~~~~~~~~~
//HashSet<word> addedPatches;
{
const polyBoundaryMesh& pbm = mesh.boundaryMesh();
forAll(selectors, selectorI)
{
const dictionary& dict = selectors[selectorI].dict();
const word& groupName = selectors[selectorI].name();
if (dict.found("patches"))
{
const dictionary& patchSources = dict.subDict("patches");
forAllConstIter(dictionary, patchSources, iter)
{
const word patchName(iter().dict()["name"]);
if (pbm.findPatchID(patchName) == -1)
{
dictionary patchDict = iter().dict();
patchDict.set("nFaces", 0);
patchDict.set("startFace", 0);
// Note: do not set coupleGroup if constructed from
// baffles so you have freedom specifying it
// yourself.
//patchDict.set("coupleGroup", groupName);
addPatch(mesh, patchName, groupName, patchDict);
}
else
{
Info<< "Patch '" << patchName
<< "' already exists. Only "
<< "moving patch faces - type will remain the same"
<< endl;
}
}
}
else
{
const dictionary& patchSource = dict.subDict("patchPairs");
const word masterName = groupName + "_master";
const word slaveName = groupName + "_slave";
word groupNameMaster = groupName;
word groupNameSlave = groupName;
dictionary patchDictMaster(patchSource);
patchDictMaster.set("nFaces", 0);
patchDictMaster.set("startFace", 0);
patchDictMaster.set("coupleGroup", groupName);
dictionary patchDictSlave(patchDictMaster);
// Note: This is added for the particular case where we want
// master and slave in different groupNames
// (ie 3D thermal baffles)
Switch sameGroup
(
patchSource.lookupOrDefault("sameGroup", true)
);
if (!sameGroup)
{
groupNameMaster = groupName + "Group_master";
groupNameSlave = groupName + "Group_slave";
patchDictMaster.set("coupleGroup", groupNameMaster);
patchDictSlave.set("coupleGroup", groupNameSlave);
}
addPatch(mesh, masterName, groupNameMaster, patchDictMaster);
addPatch(mesh, slaveName, groupNameSlave, patchDictSlave);
}
}
}
// Make sure patches and zoneFaces are synchronised across couples
mesh.boundaryMesh().checkParallelSync(true);
mesh.faceZones().checkParallelSync(true);
// Mesh change container
polyTopoChange meshMod(mesh);
const polyBoundaryMesh& pbm = mesh.boundaryMesh();
// Do the actual changes. Note:
// - loop in incrementing face order (not necessary if faceZone ordered).
// Preserves any existing ordering on patch faces.
// - two passes, do non-flip faces first and flip faces second. This
// guarantees that when e.g. creating a cyclic all faces from one
// side come first and faces from the other side next.
// Whether first use of face (modify) or consecutive (add)
PackedBoolList modifiedFace(mesh.nFaces());
label nModified = 0;
forAll(selectors, selectorI)
{
const word& name = selectors[selectorI].name();
label zoneID = mesh.faceZones().findZoneID(name);
const faceZone& fZone = mesh.faceZones()[zoneID];
const dictionary& dict = selectors[selectorI].dict();
DynamicList<label> newMasterPatches;
DynamicList<label> newSlavePatches;
if (dict.found("patches"))
{
const dictionary& patchSources = dict.subDict("patches");
bool master = true;
forAllConstIter(dictionary, patchSources, iter)
{
const word patchName(iter().dict()["name"]);
label patchi = pbm.findPatchID(patchName);
if (master)
{
newMasterPatches.append(patchi);
}
else
{
newSlavePatches.append(patchi);
}
master = !master;
}
}
else
{
const word masterName = selectors[selectorI].name() + "_master";
newMasterPatches.append(pbm.findPatchID(masterName));
const word slaveName = selectors[selectorI].name() + "_slave";
newSlavePatches.append(pbm.findPatchID(slaveName));
}
createFaces
(
internalFacesOnly,
mesh,
fZone,
newMasterPatches,
newSlavePatches,
meshMod,
modifiedFace,
nModified
);
}
Info<< "Converted " << returnReduce(nModified, sumOp<label>())
<< " faces into boundary faces in patches "
<< bafflePatches.sortedToc() << nl << endl;
if (!overwrite)
{
runTime++;
}
// Change the mesh. Change points directly (no inflation).
autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
// Update fields
mesh.updateMesh(map);
// Correct boundary faces mapped-out-of-nothing.
// This is just a hack to correct the value field.
{
fvMeshMapper mapper(mesh, map);
bool hasWarned = false;
forAllConstIter(HashSet<word>, bafflePatches, iter)
{
label patchi = mesh.boundaryMesh().findPatchID(iter.key());
const fvPatchMapper& pm = mapper.boundaryMap()[patchi];
if (pm.sizeBeforeMapping() == 0)
{
if (!hasWarned)
{
hasWarned = true;
WarningInFunction
<< "Setting field on boundary faces to zero." << endl
<< "You might have to edit these fields." << endl;
}
fvMeshTools::zeroPatchFields(mesh, patchi);
}
}
}
// Pass 2: change patchFields
// ~~~~~~~~~~~~~~~~~~~~~~~~~~
{
const polyBoundaryMesh& pbm = mesh.boundaryMesh();
forAll(selectors, selectorI)
{
const dictionary& dict = selectors[selectorI].dict();
if (dict.found("patches"))
{
const dictionary& patchSources = dict.subDict("patches");
forAllConstIter(dictionary, patchSources, iter)
{
const word patchName(iter().dict()["name"]);
label patchi = pbm.findPatchID(patchName);
if (iter().dict().found("patchFields"))
{
const dictionary& patchFieldsDict =
iter().dict().subDict
(
"patchFields"
);
fvMeshTools::setPatchFields
(
mesh,
patchi,
patchFieldsDict
);
}
}
}
else
{
const dictionary& patchSource = dict.subDict("patchPairs");
Switch sameGroup
(
patchSource.lookupOrDefault("sameGroup", true)
);
const word& groupName = selectors[selectorI].name();
if (patchSource.found("patchFields"))
{
dictionary patchFieldsDict = patchSource.subDict
(
"patchFields"
);
if (sameGroup)
{
// Add coupleGroup to all entries
forAllIter(dictionary, patchFieldsDict, iter)
{
if (iter().isDict())
{
dictionary& dict = iter().dict();
dict.set("coupleGroup", groupName);
}
}
const labelList& patchIDs =
pbm.groupPatchIDs()[groupName];
forAll(patchIDs, i)
{
fvMeshTools::setPatchFields
(
mesh,
patchIDs[i],
patchFieldsDict
);
}
}
else
{
const word masterPatchName(groupName + "_master");
const word slavePatchName(groupName + "_slave");
label patchiMaster = pbm.findPatchID(masterPatchName);
label patchiSlave = pbm.findPatchID(slavePatchName);
fvMeshTools::setPatchFields
(
mesh,
patchiMaster,
patchFieldsDict
);
fvMeshTools::setPatchFields
(
mesh,
patchiSlave,
patchFieldsDict
);
}
}
}
}
}
// Move mesh (since morphing might not do this)
if (map().hasMotionPoints())
{
mesh.movePoints(map().preMotionPoints());
}
// Remove any now zero-sized patches
filterPatches(mesh, bafflePatches);
if (overwrite)
{
mesh.setInstance(oldInstance);
}
Info<< "Writing mesh to " << runTime.timeName() << endl;
mesh.write();
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //
|