Loading...
Searching...
No Matches
createPatch.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Application
28 createPatch
29
30Group
31 grpMeshManipulationUtilities
32
33Description
34 Create patches out of selected boundary faces, which are either
35 from existing patches or from a faceSet.
36
37 More specifically it:
38 - creates new patches (from selected boundary faces).
39 Synchronise faces on coupled patches.
40 - synchronises points on coupled boundaries
41 - remove patches with 0 faces in them
42
43\*---------------------------------------------------------------------------*/
44
45#include "argList.H"
46#include "Time.H"
47#include "OFstream.H"
48#include "meshTools.H"
49#include "faceSet.H"
50#include "IOPtrList.H"
51#include "cyclicPolyPatch.H"
52#include "syncTools.H"
53#include "polyTopoChange.H"
54#include "polyModifyFace.H"
55#include "polyAddFace.H"
56#include "wordRes.H"
57#include "processorMeshes.H"
58#include "IOdictionary.H"
59#include "regionProperties.H"
60#include "faceAreaWeightAMI2D.H"
61#include "fvMeshTools.H"
62#include "ReadFields.H"
63
64using namespace Foam;
65
66// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67
68namespace Foam
69{
71}
72
73
74word patchName(const word& name, const fvMesh& mesh0, const fvMesh& mesh1)
75{
76 word pName(name != "none" ? name : word::null);
77 pName += mesh0.name();
78 pName += "_to_";
79 pName += mesh1.name();
80 return pName;
81}
82
83
84void matchPatchFaces
85(
86 const word& entryName,
87 const word& AMIMethod,
88 const dictionary& AMIDict,
90
91 const label meshi,
92 const label nSourcei,
93 const label sourcei,
94 const labelList& patchesi,
95
96 const label meshj,
97 const label nSourcej,
98 const label sourcej,
99 const labelList& patchesj,
100
101 DynamicList<labelList>& interfaceMesh0,
102 DynamicList<label>& interfaceSource0,
103 DynamicList<labelList>& interfacePatch0,
104 DynamicList<wordList>& interfaceNames0,
105 DynamicList<List<DynamicList<label>>>& interfaceFaces0,
106
107 DynamicList<labelList>& interfaceMesh1,
108 DynamicList<label>& interfaceSource1,
109 DynamicList<labelList>& interfacePatch1,
110 DynamicList<wordList>& interfaceNames1,
111 DynamicList<List<DynamicList<label>>>& interfaceFaces1
112)
113{
114 // Now we have:
115 // - meshi, sourcei, patchesi
116 // - meshj, sourcej, patchesj
117
118
119 // Attempt to match patches
120 forAll(patchesi, i)
121 {
122 const auto& ppi = meshes[meshi].boundaryMesh()[patchesi[i]];
123
124 forAll(patchesj, j)
125 {
126 const auto& ppj = meshes[meshj].boundaryMesh()[patchesj[j]];
127
128 // Use AMI to try and find matches
129 auto AMPtr(AMIInterpolation::New(AMIMethod, AMIDict));
130
131 AMPtr->calculate(ppi, ppj, nullptr);
132 if
133 (
134 gAverage(AMPtr->tgtWeightsSum()) > SMALL
135 || gAverage(AMPtr->srcWeightsSum()) > SMALL
136 )
137 {
138 const label inti = interfaceMesh0.size();
139
140 Info<< "Introducing interface " << inti << " between"
141 << " mesh " << meshes[meshi].name()
142 //<< " source:" << sourcei
143 << " patch " << ppi.name()
144 << " and mesh " << meshes[meshj].name()
145 //<< " source:" << sourcej
146 << " patch " << ppj.name()
147 << endl;
148
149 // Mesh 0
150 //~~~~~~~
151
152 auto& intMesh0 = interfaceMesh0.emplace_back(nSourcei, -1);
153 intMesh0[sourcei] = meshi;
154
155 interfaceSource0.push_back(sourcei);
156
157 auto& intPatch0 = interfacePatch0.emplace_back(nSourcei, -1);
158 intPatch0[sourcei] = ppi.index();
159
160 auto& intNames0 = interfaceNames0.emplace_back(nSourcei);
161 intNames0[sourcei] =
162 patchName(entryName, meshes[meshi], meshes[meshj]);
163
164
165 // Mesh 1
166 //~~~~~~~
167
168 auto& intMesh1 = interfaceMesh1.emplace_back(nSourcej, -1);
169 intMesh1[sourcej] = meshj;
170
171 interfaceSource1.push_back(sourcej);
172
173 auto& intPatch1 = interfacePatch1.emplace_back(nSourcej, -1);
174 intPatch1[sourcej] = ppj.index();
175
176 auto& intNames1 = interfaceNames1.emplace_back(nSourcej);
177 intNames1[sourcej] =
178 patchName(entryName, meshes[meshj], meshes[meshi]);
179
180 auto& intFaces0 = interfaceFaces0.emplace_back(nSourcei);
181 DynamicList<label>& faces0 = intFaces0[sourcei];
182 faces0.setCapacity(ppi.size());
183
184 auto& intFaces1 = interfaceFaces1.emplace_back(nSourcej);
185 DynamicList<label>& faces1 = intFaces1[sourcej];
186 faces1.setCapacity(ppj.size());
187
188
189 // Mark any mesh faces:
190 // - that have a weight > 0.5.
191 // - contribute to these faces (even if < 0.5)
192
193 faces0.clear();
194 faces1.clear();
195
196 // Marked as donor of face with high weight
197 scalarField targetMask;
198 {
199 const scalarField& weights = AMPtr->srcWeightsSum();
200 scalarField mask(weights.size(), Zero);
201 forAll(weights, facei)
202 {
203 const scalar sum = weights[facei];
204 if (sum > 0.5)
205 {
206 mask[facei] = sum;
207 }
208 }
209
210 // Push source field mask to target
211 targetMask = AMPtr->interpolateToTarget(mask);
212 }
213
214 scalarField sourceMask;
215 {
216 const scalarField& weights = AMPtr->tgtWeightsSum();
217 scalarField mask(weights.size(), Zero);
218 forAll(weights, facei)
219 {
220 const scalar sum = weights[facei];
221 if (sum > 0.5)
222 {
223 mask[facei] = sum;
224 }
225 }
226
227 // Push target mask back to source
228 sourceMask = AMPtr->interpolateToSource(mask);
229 }
230
231 {
232 const scalarField& weights = AMPtr->srcWeightsSum();
233 forAll(weights, facei)
234 {
235 if (weights[facei] > 0.5 || sourceMask[facei] > SMALL)
236 {
237 faces0.push_back(ppi.start()+facei);
238 }
239 }
240 }
241 {
242 const scalarField& weights = AMPtr->tgtWeightsSum();
243 forAll(weights, facei)
244 {
245 if (weights[facei] > 0.5 || targetMask[facei] > SMALL)
246 {
247 faces1.push_back(ppj.start()+facei);
248 }
249 }
250 }
251
252 faces0.shrink();
253 faces1.shrink();
254 }
255 }
256 }
257}
258
259
260void matchPatchFaces
261(
262 const bool includeOwn,
263 const PtrList<fvMesh>& meshes,
264 const List<DynamicList<label>>& interRegionSources,
265 const List<wordList>& patchNames,
266 const labelListListList& matchPatchIDs,
267 const List<wordList>& AMIMethods,
268 List<PtrList<dictionary>> patchInfoDicts,
269
270 DynamicList<labelList>& interfaceMesh0,
271 DynamicList<label>& interfaceSource0,
272 DynamicList<labelList>& interfacePatch0,
273 DynamicList<List<DynamicList<label>>>& interfaceFaces0,
274 DynamicList<wordList>& interfaceNames0,
275
276 DynamicList<labelList>& interfaceMesh1,
277 DynamicList<label>& interfaceSource1,
278 DynamicList<labelList>& interfacePatch1,
279 DynamicList<List<DynamicList<label>>>& interfaceFaces1,
280 DynamicList<wordList>& interfaceNames1
281)
282{
283 // Add approximate matches
284 forAll(meshes, meshi)
285 {
286 const labelListList& allPatchesi = matchPatchIDs[meshi];
287
288 const label meshjStart(includeOwn ? meshi : meshi+1);
289
290 for (label meshj = meshjStart; meshj < meshes.size(); meshj++)
291 {
292 const labelListList& allPatchesj = matchPatchIDs[meshj];
293
294 for (const label sourcei : interRegionSources[meshi])
295 {
296 const labelList& patchesi = allPatchesi[sourcei];
297
298 for (const label sourcej : interRegionSources[meshj])
299 {
300 const labelList& patchesj = allPatchesj[sourcej];
301
302 // Now we have:
303 // - meshi, sourcei, patchesi
304 // - meshj, sourcej, patchesj
305
306 matchPatchFaces
307 (
308 patchNames[meshi][sourcei],
309 (
310 AMIMethods[meshi][sourcei].size()
311 ? AMIMethods[meshi][sourcei]
312 : faceAreaWeightAMI2D::typeName
313 ),
314 patchInfoDicts[meshi][sourcei],
315 meshes,
316
317 meshi,
318 allPatchesi.size(), // nSourcei,
319 sourcei,
320 patchesi,
321
322 meshj,
323 allPatchesj.size(), // nSourcej,
324 sourcej,
325 patchesj,
326
327 interfaceMesh0,
328 interfaceSource0,
329 interfacePatch0,
330 interfaceNames0,
331 interfaceFaces0,
332
333 interfaceMesh1,
334 interfaceSource1,
335 interfacePatch1,
336 interfaceNames1,
337 interfaceFaces1
338 );
339 }
340 }
341 }
342 }
343}
344
345
346void changePatchID
347(
348 const bool modify,
349 const fvMesh& mesh,
350 const label faceID,
351 const label patchID,
352 polyTopoChange& meshMod
353)
354{
355 const label zoneID = mesh.faceZones().whichZone(faceID);
356
357 bool zoneFlip = false;
358
359 if (zoneID >= 0)
360 {
361 const faceZone& fZone = mesh.faceZones()[zoneID];
362
363 zoneFlip = fZone.flipMap()[fZone.whichFace(faceID)];
364 }
365
366 if (modify)
367 {
368 meshMod.setAction
369 (
371 (
372 mesh.faces()[faceID], // face
373 faceID, // face ID
374 mesh.faceOwner()[faceID], // owner
375 -1, // neighbour
376 false, // flip flux
377 patchID, // patch ID
378 false, // remove from zone
379 zoneID, // zone ID
380 zoneFlip // zone flip
381 )
382 );
383 }
384 else
385 {
386 meshMod.setAction
387 (
389 (
390 mesh.faces()[faceID], // modified face
391 mesh.faceOwner()[faceID], // owner
392 -1, // neighbour
393 -1, // master point
394 -1, // master edge
395 faceID, // master face
396 false, // face flip
397 patchID, // patch for face
398 zoneID, // zone for face
399 zoneFlip // face flip in zone
400 )
401 );
402 }
403}
404
405
406void changePatchID
407(
408 const fvMesh& mesh,
409 const labelList& faceLabels,
410 const label patchID,
411 bitSet& isRepatchedBoundary,
412 polyTopoChange& meshMod
413)
414{
415 for (const label facei : faceLabels)
416 {
417 if (mesh.isInternalFace(facei))
418 {
420 << "Face " << facei
421 << " is not an external face of the mesh." << endl
422 << "This application can only repatch"
423 << " existing boundary faces." << exit(FatalError);
424 }
425
426 const bool isFirst =
427 isRepatchedBoundary.set(facei-mesh.nInternalFaces());
428 if (!isFirst)
429 {
430 static label nWarnings = 0;
431 if (nWarnings == 0)
432 {
433 const label newPatchi = meshMod.region()[facei];
434 //FatalErrorInFunction
436 << "Face " << facei
437 << " at " << mesh.faceCentres()[facei]
438 << " marked for patch " << patchID
439 << " name " << mesh.boundaryMesh()[patchID].name()
440 << " is already marked for patch " << newPatchi
441 << " name " << mesh.boundaryMesh()[newPatchi].name()
442 << ". Creating duplicate face. Suppressing further warnings"
443 //<< exit(FatalError);
444 << endl;
445 }
446 nWarnings++;
447 }
448
449 changePatchID(isFirst, mesh, facei, patchID, meshMod);
450 }
451}
452
453
454// Dump for all patches the current match
455void dumpCyclicMatch(const fileName& prefix, const polyMesh& mesh)
456{
457 for (const polyPatch& pp : mesh.boundaryMesh())
458 {
460
461 if (cpp && cpp->owner())
462 {
463 const auto& cycPatch = *cpp;
464 const auto& nbrPatch = cycPatch.neighbPatch();
465
466 // Dump patches
467 {
468 OFstream str(prefix+cycPatch.name()+".obj");
469 Pout<< "Dumping " << cycPatch.name()
470 << " faces to " << str.name() << endl;
472 (
473 str,
474 cycPatch,
475 cycPatch.points()
476 );
477 }
478
479 {
480 OFstream str(prefix+nbrPatch.name()+".obj");
481 Pout<< "Dumping " << nbrPatch.name()
482 << " faces to " << str.name() << endl;
484 (
485 str,
486 nbrPatch,
487 nbrPatch.points()
488 );
489 }
490
491
492 // Lines between corresponding face centres
493 OFstream str(prefix+cycPatch.name()+nbrPatch.name()+"_match.obj");
494 label vertI = 0;
495
496 Pout<< "Dumping cyclic match as lines between face centres to "
497 << str.name() << endl;
498
499 forAll(cycPatch, facei)
500 {
501 const point& fc0 = mesh.faceCentres()[cycPatch.start()+facei];
502 meshTools::writeOBJ(str, fc0);
503 vertI++;
504 const point& fc1 = mesh.faceCentres()[nbrPatch.start()+facei];
505 meshTools::writeOBJ(str, fc1);
506 vertI++;
507
508 str<< "l " << vertI-1 << ' ' << vertI << nl;
509 }
510 }
511 }
512}
513
514
515void separateList
516(
517 const vectorField& separation,
519)
520{
521 if (separation.size() == 1)
522 {
523 // Single value for all.
524
525 forAll(field, i)
526 {
527 field[i] += separation[0];
528 }
529 }
530 else if (separation.size() == field.size())
531 {
532 forAll(field, i)
533 {
534 field[i] += separation[i];
535 }
536 }
537 else
538 {
540 << "Sizes of field and transformation not equal. field:"
541 << field.size() << " transformation:" << separation.size()
542 << abort(FatalError);
543 }
544}
545
546
547// Synchronise points on both sides of coupled boundaries.
548template<class CombineOp>
549void syncPoints
550(
551 const polyMesh& mesh,
553 const CombineOp& cop,
554 const point& nullValue
555)
556{
557 if (points.size() != mesh.nPoints())
558 {
560 << "Number of values " << points.size()
561 << " is not equal to the number of points in the mesh "
562 << mesh.nPoints() << abort(FatalError);
563 }
564
565 const polyBoundaryMesh& patches = mesh.boundaryMesh();
566
567 // Is there any coupled patch with transformation?
568 bool hasTransformation = false;
569
570 if (UPstream::parRun())
571 {
572 const labelList& procPatches = mesh.globalData().processorPatches();
573
574 // Send
575 for (const label patchi : procPatches)
576 {
577 const polyPatch& pp = patches[patchi];
578 const auto& procPatch = refCast<const processorPolyPatch>(pp);
579
580 if (pp.nPoints() && procPatch.owner())
581 {
582 // Get data per patchPoint in neighbouring point numbers.
583 pointField patchInfo(procPatch.nPoints(), nullValue);
584
585 const labelList& meshPts = procPatch.meshPoints();
586 const labelList& nbrPts = procPatch.neighbPoints();
587
588 forAll(nbrPts, pointi)
589 {
590 label nbrPointi = nbrPts[pointi];
591 if (nbrPointi >= 0 && nbrPointi < patchInfo.size())
592 {
593 patchInfo[nbrPointi] = points[meshPts[pointi]];
594 }
595 }
596
597 OPstream toNbr
598 (
600 procPatch.neighbProcNo()
601 );
602 toNbr << patchInfo;
603 }
604 }
605
606
607 // Receive and set.
608
609 for (const label patchi : procPatches)
610 {
611 const polyPatch& pp = patches[patchi];
612 const auto& procPatch = refCast<const processorPolyPatch>(pp);
613
614 if (pp.nPoints() && !procPatch.owner())
615 {
616 // We do not know the number of points on the other side
617 // so cannot use UIPstream::read
618
619 pointField nbrPatchInfo;
620 IPstream::recv(nbrPatchInfo, procPatch.neighbProcNo());
621
622 // Null any value which is not on neighbouring processor
623 nbrPatchInfo.setSize(procPatch.nPoints(), nullValue);
624
625 if (!procPatch.parallel())
626 {
627 hasTransformation = true;
628 transformList(procPatch.forwardT(), nbrPatchInfo);
629 }
630 else if (procPatch.separated())
631 {
632 hasTransformation = true;
633 separateList(-procPatch.separation(), nbrPatchInfo);
634 }
635
636 const labelList& meshPts = procPatch.meshPoints();
637
638 forAll(meshPts, pointi)
639 {
640 label meshPointi = meshPts[pointi];
641 points[meshPointi] = nbrPatchInfo[pointi];
642 }
643 }
644 }
645 }
646
647 // Do the cyclics.
648 for (const polyPatch& pp : patches)
649 {
651
652 if (cpp && cpp->owner())
653 {
654 // Owner does all.
655
656 const auto& cycPatch = *cpp;
657 const auto& nbrPatch = cycPatch.neighbPatch();
658
659 const edgeList& coupledPoints = cycPatch.coupledPoints();
660 const labelList& meshPts = cycPatch.meshPoints();
661 const labelList& nbrMeshPts = nbrPatch.meshPoints();
662
663 pointField half0Values(coupledPoints.size());
664
665 forAll(coupledPoints, i)
666 {
667 const edge& e = coupledPoints[i];
668 label point0 = meshPts[e[0]];
669 half0Values[i] = points[point0];
670 }
671
672 if (!cycPatch.parallel())
673 {
674 hasTransformation = true;
675 transformList(cycPatch.reverseT(), half0Values);
676 }
677 else if (cycPatch.separated())
678 {
679 hasTransformation = true;
680 separateList(cycPatch.separation(), half0Values);
681 }
682
683 forAll(coupledPoints, i)
684 {
685 const edge& e = coupledPoints[i];
686 label point1 = nbrMeshPts[e[1]];
687 points[point1] = half0Values[i];
688 }
689 }
690 }
691
692 //- Note: hasTransformation is only used for warning messages so
693 // reduction not strictly necessary.
694 //Pstream::reduceOr(hasTransformation);
695
696 // Synchronize multiple shared points.
697 const globalMeshData& pd = mesh.globalData();
698
699 if (pd.nGlobalPoints() > 0)
700 {
701 if (hasTransformation)
702 {
704 << "There are decomposed cyclics in this mesh with"
705 << " transformations." << endl
706 << "This is not supported. The result will be incorrect"
707 << endl;
708 }
709
710
711 // Values on shared points.
712 pointField sharedPts(pd.nGlobalPoints(), nullValue);
713
715 {
716 label meshPointi = pd.sharedPointLabels()[i];
717 // Fill my entries in the shared points
718 sharedPts[pd.sharedPointAddr()[i]] = points[meshPointi];
719 }
720
721 // Combine - globally consistent
722 Pstream::listCombineReduce(sharedPts, cop);
723
724 // Now we will all have the same information. Merge it back with
725 // my local information.
727 {
728 label meshPointi = pd.sharedPointLabels()[i];
729 points[meshPointi] = sharedPts[pd.sharedPointAddr()[i]];
730 }
731 }
732}
733
734
735
736int main(int argc, char *argv[])
737{
739 (
740 "Create patches out of selected boundary faces, which are either"
741 " from existing patches or from a faceSet"
742 );
743
744 #include "addOverwriteOption.H"
745 #include "addAllRegionOptions.H"
746
747 argList::addOption("dict", "file", "Alternative createPatchDict");
749 (
750 "writeObj",
751 "Write obj files showing the cyclic matching process"
752 );
753
754 argList::noFunctionObjects(); // Never use function objects
755
756 #include "setRootCase.H"
757 #include "createTime.H"
758 #include "getAllRegionOptions.H"
759
760 const bool overwrite = args.found("overwrite");
761
762 #include "createNamedMeshes.H"
763
764 const bool writeObj = args.found("writeObj");
765
766
767 // Read dictionaries and extract various
768 wordList oldInstances(meshes.size());
769 PtrList<dictionary> dicts(meshes.size());
770 List<wordList> patchNames(meshes.size());
771 List<labelListList> matchPatchIDs(meshes.size());
772 List<wordList> matchMethods(meshes.size());
773 List<PtrList<dictionary>> patchInfoDicts(meshes.size());
774
775 List<DynamicList<label>> interRegionSources(meshes.size());
776 forAll(meshes, meshi)
777 {
778 fvMesh& mesh = meshes[meshi];
779 const polyBoundaryMesh& patches = mesh.boundaryMesh();
780
781 // If running parallel check same patches everywhere
782 patches.checkParallelSync(true);
783
784 oldInstances[meshi] = mesh.pointsInstance();
785
786 const word dictName("createPatchDict");
788 Info<< "Reading " << dictIO.instance()/dictIO.name() << nl << endl;
789
790 dicts.set(meshi, new IOdictionary(dictIO));
791
792 const auto& dict = dicts[meshi];
793 PtrList<dictionary> patchSources(dict.lookup("patches"));
794 patchNames[meshi].setSize(patchSources.size());
795 matchPatchIDs[meshi].setSize(patchSources.size());
796 matchMethods[meshi].setSize(patchSources.size());
797 patchInfoDicts[meshi].setSize(patchSources.size());
798
799 // Read patch construct info from dictionary
800 forAll(patchSources, sourcei)
801 {
802 const auto& pDict = patchSources[sourcei];
803 patchNames[meshi][sourcei] = pDict.getOrDefault<word>
804 (
805 "name",
808 );
809
810 patchInfoDicts[meshi].set
811 (
812 sourcei,
813 new dictionary(pDict.subDict("patchInfo"))
814 );
815 dictionary& patchDict = patchInfoDicts[meshi][sourcei];
816 if (patchDict.found("AMIMethod"))
817 {
818 matchMethods[meshi][sourcei] = patchDict.get<word>("AMIMethod");
819 // Disable full matching since we're trying to use AMIMethod to
820 // find out actual overlap
821 patchDict.add("requireMatch", false);
822 }
823
824 wordRes matchNames;
825 if (pDict.readIfPresent("patches", matchNames, keyType::LITERAL))
826 {
827 matchPatchIDs[meshi][sourcei] =
828 patches.patchSet(matchNames).sortedToc();
829 }
830
831 if (pDict.get<word>("constructFrom") == "autoPatch")
832 {
833 interRegionSources[meshi].push_back(sourcei);
834 }
835 }
836 }
837
838
839
840 // Determine matching faces. This is an interface between two
841 // sets of faces:
842 // - originating mesh
843 // - originating patch
844 // - originating (mesh)facelabels
845 // It matches all mesh against each other. Lower numbered mesh gets
846 // postfix 0, higher numbered mesh postfix 1.
847
848 // Per interface, per patchSource:
849 // 1. the lower numbered mesh
850 DynamicList<labelList> interfaceMesh0;
851 // 1b. the source index (i.e. the patch dictionary)
852 DynamicList<label> interfaceSource0;
853 // 2. the patch on the interfaceMesh0
854 DynamicList<labelList> interfacePatch0;
855 // 3. the facelabels on the interfaceMesh0
856 DynamicList<List<DynamicList<label>>> interfaceFaces0;
857 // 4. generated interface name
858 DynamicList<wordList> interfaceNames0;
859
860 // Same for the higher numbered mesh
861 DynamicList<labelList> interfaceMesh1;
862 DynamicList<label> interfaceSource1;
863 DynamicList<labelList> interfacePatch1;
864 DynamicList<List<DynamicList<label>>> interfaceFaces1;
865 DynamicList<wordList> interfaceNames1;
866 {
867 // Whether to match to patches in own mesh
868 const bool includeOwn = (meshes.size() == 1);
869
870 matchPatchFaces
871 (
872 includeOwn,
873 meshes,
874 interRegionSources,
876 matchPatchIDs,
877 matchMethods,
878 patchInfoDicts,
879
880 interfaceMesh0,
881 interfaceSource0,
882 interfacePatch0,
883 interfaceFaces0,
884 interfaceNames0,
885
886 interfaceMesh1,
887 interfaceSource1,
888 interfacePatch1,
889 interfaceFaces1,
890 interfaceNames1
891 );
892 }
893
894
895
896 // Read fields
897 List<PtrList<volScalarField>> vsFlds(meshes.size());
898 List<PtrList<volVectorField>> vvFlds(meshes.size());
899 List<PtrList<volSphericalTensorField>> vstFlds(meshes.size());
900 List<PtrList<volSymmTensorField>> vsymtFlds(meshes.size());
901 List<PtrList<surfaceScalarField>> ssFlds(meshes.size());
902 List<PtrList<volTensorField>> vtFlds(meshes.size());
903 List<PtrList<surfaceVectorField>> svFlds(meshes.size());
904 List<PtrList<surfaceSphericalTensorField>> sstFlds(meshes.size());
905 List<PtrList<surfaceSymmTensorField>> ssymtFlds(meshes.size());
906 List<PtrList<surfaceTensorField>> stFlds(meshes.size());
907
908 {
909 forAll(meshes, meshi)
910 {
911 const fvMesh& mesh = meshes[meshi];
912
913 bool noFields = true;
914 for (const auto& d : patchInfoDicts[meshi])
915 {
916 if (d.found("patchFields"))
917 {
918 noFields = false;
919 }
920 }
921
922 if (!noFields)
923 {
924 // Read objects in time directory
925 IOobjectList objects(mesh, runTime.timeName());
926
927 // Read volume fields.
928 ReadFields(mesh, objects, vsFlds[meshi]);
929 ReadFields(mesh, objects, vvFlds[meshi]);
930 ReadFields(mesh, objects, vstFlds[meshi]);
931 ReadFields(mesh, objects, vsymtFlds[meshi]);
932 ReadFields(mesh, objects, vtFlds[meshi]);
933
934 // Read surface fields.
935 ReadFields(mesh, objects, ssFlds[meshi]);
936 ReadFields(mesh, objects, svFlds[meshi]);
937 ReadFields(mesh, objects, sstFlds[meshi]);
938 ReadFields(mesh, objects, ssymtFlds[meshi]);
939 ReadFields(mesh, objects, stFlds[meshi]);
940 }
941 }
942 }
943
944
945 // Maintain list of added patches so we exclude them from filtering
946 // later on
947 List<DynamicList<word>> allAddedPatches(meshes.size());
948
949 // Loop over all regions
950
951 forAll(meshes, meshi)
952 {
953 fvMesh& mesh = meshes[meshi];
954
955 Info<< "\n\nAdding patches to mesh " << mesh.name() << nl << endl;
956
957 const polyBoundaryMesh& patches = mesh.boundaryMesh();
958 const dictionary& dict = dicts[meshi];
959
960 if (writeObj)
961 {
962 dumpCyclicMatch("initial_", mesh);
963 }
964
965 // Read patch construct info from dictionary
966 PtrList<dictionary> patchSources(dict.lookup("patches"));
967
968
969 // 1. Add all new patches
970 // ~~~~~~~~~~~~~~~~~~~~~~
971
972 forAll(patchSources, sourcei)
973 {
974 const dictionary& dict = patchSources[sourcei];
975 const word sourceType(dict.get<word>("constructFrom"));
976 const word patchName
977 (
978 dict.getOrDefault<word>
979 (
980 "name",
983 )
984 );
985
986 dictionary patchDict(patchInfoDicts[meshi][sourcei]);
987 patchDict.set("nFaces", 0);
988 patchDict.set("startFace", 0); // Gets overwritten
989
990
991 if (sourceType == "autoPatch")
992 {
993 // Special : automatically create the necessary inter-region
994 // coupling patches
995
996 forAll(interfaceMesh0, inti)
997 {
998 const labelList& allMeshes0 = interfaceMesh0[inti];
999 const wordList& allNames0 = interfaceNames0[inti];
1000 const labelList& allMeshes1 = interfaceMesh1[inti];
1001 const wordList& allNames1 = interfaceNames1[inti];
1002
1003 if
1004 (
1005 interfaceSource0[inti] == sourcei
1006 && allMeshes0[sourcei] == meshi
1007 )
1008 {
1009 // Current mesh is mesh0. mesh1 is the remote mesh.
1010 const label sourcej = interfaceSource1[inti];
1011 const word& patchName = allNames0[sourcei];
1012 if (patches.findPatchID(patchName) == -1)
1013 {
1014 dictionary allDict(patchDict);
1015 const auto& mesh1 = meshes[allMeshes1[sourcej]];
1016 allDict.set("sampleRegion", mesh1.name());
1017 const auto& destPatch = allNames1[sourcej];
1018 allDict.set("samplePatch", destPatch);
1019 allDict.set("neighbourPatch", destPatch);
1020
1021 Info<< "Adding new patch " << patchName
1022 << " from " << allDict << endl;
1023
1024 autoPtr<polyPatch> ppPtr
1025 (
1027 (
1028 patchName,
1029 allDict,
1030 0, // overwritten
1031 patches
1032 )
1033 );
1035 (
1036 mesh,
1037 ppPtr(),
1038 patchDict.subOrEmptyDict("patchFields"),
1040 true
1041 );
1042 allAddedPatches[meshi].append(ppPtr->name());
1043 }
1044 }
1045 }
1046
1047 forAll(interfaceMesh1, inti)
1048 {
1049 const labelList& allMeshes0 = interfaceMesh0[inti];
1050 const wordList& allNames0 = interfaceNames0[inti];
1051 const labelList& allMeshes1 = interfaceMesh1[inti];
1052 const wordList& allNames1 = interfaceNames1[inti];
1053
1054 if
1055 (
1056 interfaceSource1[inti] == sourcei
1057 && allMeshes1[sourcei] == meshi
1058 )
1059 {
1060 // Current mesh is mesh1. mesh0 is the remote mesh.
1061
1062 const label sourcej = interfaceSource0[inti];
1063 const word& patchName = allNames1[sourcei];
1064 if (patches.findPatchID(patchName) == -1)
1065 {
1066 dictionary allDict(patchDict);
1067 const auto& destPatch = allNames0[sourcej];
1068 const auto& mesh0 = meshes[allMeshes0[sourcej]];
1069 allDict.set("sampleRegion", mesh0.name());
1070 allDict.set("samplePatch", destPatch);
1071 allDict.set("neighbourPatch", destPatch);
1072
1073 Info<< "Adding new patch " << patchName
1074 << " from " << allDict << endl;
1075
1076 autoPtr<polyPatch> ppPtr
1077 (
1079 (
1080 patchName,
1081 allDict,
1082 0, // overwritten
1083 patches
1084 )
1085 );
1087 (
1088 mesh,
1089 ppPtr(),
1090 patchDict.subOrEmptyDict("patchFields"),
1092 true
1093 );
1094 allAddedPatches[meshi].append(ppPtr->name());
1095 }
1096 }
1097 }
1098 }
1099 else
1100 {
1101 if (patches.findPatchID(patchName) == -1)
1102 {
1103 Info<< "Adding new patch " << patchName
1104 << " from " << patchDict << endl;
1105
1106 autoPtr<polyPatch> ppPtr
1107 (
1109 (
1110 patchName,
1111 patchDict,
1112 0, // overwritten
1113 patches
1114 )
1115 );
1117 (
1118 mesh,
1119 ppPtr(),
1120 patchDict.subOrEmptyDict("patchFields"),
1122 true
1123 );
1124 allAddedPatches[meshi].append(ppPtr->name());
1125 }
1126 }
1127 }
1128
1129 Info<< endl;
1130 }
1131
1132
1133 // 2. Repatch faces
1134 // ~~~~~~~~~~~~~~~~
1135
1136 forAll(meshes, meshi)
1137 {
1138 fvMesh& mesh = meshes[meshi];
1139
1140 Info<< "\n\nRepatching mesh " << mesh.name() << nl << endl;
1141
1142
1143 const polyBoundaryMesh& patches = mesh.boundaryMesh();
1144 const dictionary& dict = dicts[meshi];
1145
1146 // Whether to synchronise points
1147 const bool pointSync(dict.get<bool>("pointSync"));
1148
1149 // Read patch construct info from dictionary
1150 PtrList<dictionary> patchSources(dict.lookup("patches"));
1151
1152
1153 polyTopoChange meshMod(mesh);
1154
1155 // Mark all repatched faces. This makes sure that the faces to repatch
1156 // do not overlap
1157 bitSet isRepatchedBoundary(mesh.nBoundaryFaces());
1158
1159 forAll(patchSources, sourcei)
1160 {
1161 const dictionary& dict = patchSources[sourcei];
1162 const word patchName
1163 (
1164 dict.getOrDefault<word>
1165 (
1166 "name",
1167 word::null,
1169 )
1170 );
1171 const word sourceType(dict.get<word>("constructFrom"));
1172
1173 if (sourceType == "autoPatch")
1174 {
1175 forAll(interfaceMesh0, inti)
1176 {
1177 const labelList& allMeshes0 = interfaceMesh0[inti];
1178 const wordList& allNames0 = interfaceNames0[inti];
1179 const auto& allFaces0 = interfaceFaces0[inti];
1180
1181 const labelList& allMeshes1 = interfaceMesh1[inti];
1182 const wordList& allNames1 = interfaceNames1[inti];
1183 const auto& allFaces1 = interfaceFaces1[inti];
1184
1185 if
1186 (
1187 interfaceSource0[inti] == sourcei
1188 && allMeshes0[sourcei] == meshi
1189 )
1190 {
1191 // Current mesh is mesh0. mesh1 is the remote mesh.
1192
1193 const label destPatchi =
1194 patches.findPatchID(allNames0[sourcei], false);
1195
1196 //const auto& mesh1 =
1197 // meshes[allMeshes1[interfaceSource1[inti]]];
1198 //Pout<< "Matched mesh:" << mesh.name()
1199 // << " to mesh:" << mesh1.name()
1200 // << " through:" << allNames0[sourcei] << endl;
1201
1202 changePatchID
1203 (
1204 mesh,
1205 allFaces0[sourcei],
1206 destPatchi,
1207 isRepatchedBoundary,
1208 meshMod
1209 );
1210 }
1211 if
1212 (
1213 interfaceSource1[inti] == sourcei
1214 && allMeshes1[sourcei] == meshi
1215 )
1216 {
1217 // Current mesh is mesh1. mesh0 is the remote mesh.
1218
1219 const label destPatchi =
1220 patches.findPatchID(allNames1[sourcei], false);
1221
1222 //const auto& mesh0 =
1223 // meshes[allMeshes0[interfaceSource0[inti]]];
1224 //Pout<< "Matched mesh:" << mesh.name()
1225 // << " to mesh:" << mesh0.name()
1226 // << " through:" << allNames1[sourcei] << endl;
1227
1228 changePatchID
1229 (
1230 mesh,
1231 allFaces1[sourcei],
1232 destPatchi,
1233 isRepatchedBoundary,
1234 meshMod
1235 );
1236 }
1237 }
1238 }
1239 else if (sourceType == "patches")
1240 {
1241 const label destPatchi = patches.findPatchID(patchName, false);
1242 labelHashSet patchSources
1243 (
1244 patches.patchSet(dict.get<wordRes>("patches"))
1245 );
1246
1247 // Repatch faces of the patches.
1248 for (const label patchi : patchSources)
1249 {
1250 const polyPatch& pp = patches[patchi];
1251
1252 Info<< "Moving faces from patch " << pp.name()
1253 << " to patch " << destPatchi << endl;
1254
1255 changePatchID
1256 (
1257 mesh,
1258 pp.start()+identity(pp.size()),
1259 destPatchi,
1260 isRepatchedBoundary,
1261 meshMod
1262 );
1263 }
1264 }
1265 else if (sourceType == "set")
1266 {
1267 const label destPatchi = patches.findPatchID(patchName, false);
1268 const word setName(dict.get<word>("set"));
1269
1270 faceSet set(mesh, setName);
1271
1272 Info<< "Read " << returnReduce(set.size(), sumOp<label>())
1273 << " faces from faceSet " << set.name() << endl;
1274
1275 // Sort (since faceSet contains faces in arbitrary order)
1276 labelList faceLabels(set.sortedToc());
1277
1278 changePatchID
1279 (
1280 mesh,
1281 faceLabels,
1282 destPatchi,
1283 isRepatchedBoundary,
1284 meshMod
1285 );
1286 }
1287 else
1288 {
1290 << "Invalid source type " << sourceType << endl
1291 << "Valid source types are 'patches' 'set'"
1292 << exit(FatalError);
1293 }
1294 }
1295
1296
1297 // Change mesh, use inflation to reforce calculation of transformation
1298 // tensors.
1299 Info<< "Doing topology modification to order faces." << nl << endl;
1300 autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, true);
1301
1302 if (map().hasMotionPoints())
1303 {
1304 mesh.movePoints(map().preMotionPoints());
1305 }
1306 else
1307 {
1308 // Force calculation of transformation tensors
1309 mesh.movePoints(pointField(mesh.points()));
1310 }
1311
1312
1313 // Update fields
1314 mesh.updateMesh(map());
1315
1316
1317 // Update numbering pointing to meshi
1318 const auto& oldToNew = map().reverseFaceMap();
1319 forAll(interfaceMesh0, inti)
1320 {
1321 const labelList& allMeshes0 = interfaceMesh0[inti];
1322 forAll(allMeshes0, sourcei)
1323 {
1324 if (allMeshes0[sourcei] == meshi)
1325 {
1326 inplaceRenumber(oldToNew, interfaceFaces0[inti][sourcei]);
1327 }
1328 }
1329 }
1330 forAll(interfaceMesh1, inti)
1331 {
1332 const labelList& allMeshes1 = interfaceMesh1[inti];
1333 forAll(allMeshes1, sourcei)
1334 {
1335 if (allMeshes1[sourcei] == meshi)
1336 {
1337 inplaceRenumber(oldToNew, interfaceFaces1[inti][sourcei]);
1338 }
1339 }
1340 }
1341
1342
1343 if (writeObj)
1344 {
1345 dumpCyclicMatch("coupled_", mesh);
1346 }
1347
1348 // Synchronise points.
1349 if (!pointSync)
1350 {
1351 Info<< "Not synchronising points." << nl << endl;
1352 }
1353 else
1354 {
1355 Info<< "Synchronising points." << nl << endl;
1356
1357 // This is a bit tricky. Both normal and position might be out and
1358 // current separation also includes the normal
1359 // ( separation_ = (nf&(Cr - Cf))*nf ).
1360
1361 // For cyclic patches:
1362 // - for separated ones use user specified offset vector
1363
1364 forAll(mesh.boundaryMesh(), patchi)
1365 {
1366 const polyPatch& pp = mesh.boundaryMesh()[patchi];
1367
1368 if (pp.size() && isA<coupledPolyPatch>(pp))
1369 {
1370 const coupledPolyPatch& cpp =
1372
1373 if (cpp.separated())
1374 {
1375 Info<< "On coupled patch " << pp.name()
1376 << " separation[0] was "
1377 << cpp.separation()[0] << endl;
1378
1379 if (isA<cyclicPolyPatch>(pp) && pp.size())
1380 {
1381 const auto& cycpp =
1383
1384 if
1385 (
1386 cycpp.transform()
1388 )
1389 {
1390 // Force to wanted separation
1391 Info<< "On cyclic translation patch "
1392 << pp.name()
1393 << " forcing uniform separation of "
1394 << cycpp.separationVector() << endl;
1395 const_cast<vectorField&>(cpp.separation()) =
1396 pointField(1, cycpp.separationVector());
1397 }
1398 else
1399 {
1400 const auto& nbr = cycpp.neighbPatch();
1401 const_cast<vectorField&>(cpp.separation()) =
1403 (
1404 1,
1405 nbr[0].centre(mesh.points())
1406 - cycpp[0].centre(mesh.points())
1407 );
1408 }
1409 }
1410 Info<< "On coupled patch " << pp.name()
1411 << " forcing uniform separation of "
1412 << cpp.separation() << endl;
1413 }
1414 else if (!cpp.parallel())
1415 {
1416 Info<< "On coupled patch " << pp.name()
1417 << " forcing uniform rotation of "
1418 << cpp.forwardT()[0] << endl;
1419
1420 const_cast<tensorField&>
1421 (
1422 cpp.forwardT()
1423 ).setSize(1);
1424 const_cast<tensorField&>
1425 (
1426 cpp.reverseT()
1427 ).setSize(1);
1428
1429 Info<< "On coupled patch " << pp.name()
1430 << " forcing uniform rotation of "
1431 << cpp.forwardT() << endl;
1432 }
1433 }
1434 }
1435
1436 Info<< "Synchronising points." << endl;
1437
1438 pointField newPoints(mesh.points());
1439
1440 syncPoints
1441 (
1442 mesh,
1443 newPoints,
1445 point(GREAT, GREAT, GREAT)
1446 );
1447
1448 scalarField diff(mag(newPoints-mesh.points()));
1449 Info<< "Points changed by average:" << gAverage(diff)
1450 << " max:" << gMax(diff) << nl << endl;
1451
1452 mesh.movePoints(newPoints);
1453 }
1454
1455
1456 // 3. Remove zeros-sized patches
1457 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1458
1459 Info<< "Removing patches with no faces in them." << nl << endl;
1460 const wordList oldPatchNames(mesh.boundaryMesh().names());
1461 const wordList oldPatchTypes(mesh.boundaryMesh().types());
1462 fvMeshTools::removeEmptyPatches(mesh, allAddedPatches[meshi], true);
1463 forAll(oldPatchNames, patchi)
1464 {
1465 const word& pName = oldPatchNames[patchi];
1466 if (mesh.boundaryMesh().findPatchID(pName) == -1)
1467 {
1468 Info<< "Removed zero-sized patch " << pName
1469 << " type " << oldPatchTypes[patchi]
1470 << " at position " << patchi << endl;
1471 }
1472 }
1473
1474
1475 if (writeObj)
1476 {
1477 dumpCyclicMatch("final_", mesh);
1478 }
1479 }
1480
1481 if (!overwrite)
1482 {
1483 ++runTime;
1484 }
1485
1486 forAll(meshes, meshi)
1487 {
1488 fvMesh& mesh = meshes[meshi];
1489 mesh.setInstance(overwrite ? oldInstances[meshi] : runTime.timeName());
1490 }
1491
1492 // More precision (for points data)
1494
1495 // Write resulting mesh
1496 forAll(meshes, meshi)
1497 {
1498 fvMesh& mesh = meshes[meshi];
1499
1500 // Override bcs with explicitly provided info. Done late so there
1501 // are already patch faces
1502 forAll(patchInfoDicts[meshi], sourcei)
1503 {
1504 const dictionary& patchDict = patchInfoDicts[meshi][sourcei];
1505 const word& patchName = patchNames[meshi][sourcei];
1506 const label patchID = mesh.boundary().findPatchID(patchName);
1507 if (patchID != -1 && patchDict.found("patchFields"))
1508 {
1509 const dictionary& pfd = patchDict.subDict("patchFields");
1510 fvMeshTools::setPatchFields(mesh, patchID, pfd);
1511 }
1512 }
1513
1514 Info<< "\n\nWriting repatched mesh " << mesh.name()
1515 << " to " << runTime.timeName() << nl << endl;
1516 mesh.clearOut(); // remove meshPhi
1517 mesh.write();
1520 }
1521
1522 Info<< "End\n" << endl;
1523
1524 return 0;
1525}
1526
1527
1528// ************************************************************************* //
Field reading functions for post-processing utilities.
Required Classes.
labelList faceLabels(nFaceLabels)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
static autoPtr< AMIInterpolation > New(const word &modelName, const dictionary &dict, const bool reverseTarget=false)
Selector for dictionary.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
T & emplace_back(Args &&... args)
Construct an element at the end of the list, return reference to the new list element.
void push_back(const T &val)
Copy append an element to the end of this list.
DynamicList< T, SizeMin > & shrink()
Calls shrink_to_fit() and returns a reference to the DynamicList.
void setCapacity(const label len)
Alter the size of the underlying storage.
A PtrList of objects of type <T> with automated input and output.
Definition IOPtrList.H:53
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
static unsigned int minPrecision(unsigned int prec) noexcept
Set the minimum default precision.
Definition IOstream.H:459
static void recv(Type &value, const int fromProcNo, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm, IOstreamOption::streamFormat fmt=IOstreamOption::BINARY)
Receive and deserialize a value. Uses operator>> for de-serialization.
Definition IPstream.H:80
bool set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition List.H:469
void setSize(label n)
Alias for resize().
Definition List.H:536
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
Output inter-processor communications stream.
Definition OPstream.H:53
static void listCombineReduce(UList< T > &values, CombineOp cop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Forwards to Pstream::listReduce with an in-place cop.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
@ buffered
"buffered" : (MPI_Bsend, MPI_Recv)
Definition UPstream.H:82
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
Definition argList.C:562
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition argList.C:389
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition argList.C:400
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
virtual bool separated() const
Are the planes separated.
virtual bool parallel() const
Are the cyclic planes parallel.
virtual const vectorField & separation() const
If the planes are separated the separation vector.
virtual const tensorField & forwardT() const
Return face transformation tensor.
Cyclic plane patch.
virtual bool owner() const
Does this side own the patch ?
const cyclicPolyPatch & neighbPatch() const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Find and return a sub-dictionary as a copy, otherwise return an empty dictionary.
Definition dictionary.C:521
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition dictionary.C:625
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition dictionary.C:765
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
A list of face labels.
Definition faceSet.H:50
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
label whichFace(const label meshFaceID) const
The local index of the given mesh face, -1 if not in the zone.
Definition faceZone.H:394
const boolList & flipMap() const noexcept
Return face flip map.
Definition faceZone.H:389
A class for handling file names.
Definition fileName.H:75
static labelList removeEmptyPatches(fvMesh &, const bool validBoundary)
Remove zero sized patches. All but processor patches are.
static label addPatch(fvMesh &mesh, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add patch. Inserts patch before all processor patches.
Definition fvMeshTools.C:38
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const word & name() const
Return reference to name.
Definition fvMesh.H:387
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
const labelList & sharedPointAddr() const
Return addressing into the complete globally shared points list.
label nGlobalPoints() const
Return number of globally shared points.
const labelList & sharedPointLabels() const
Return indices of local points that are globally shared.
@ LITERAL
String literal.
Definition keyType.H:82
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition polyAddFace.H:50
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Class describing modification of a face.
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return pointer to a new patch created on freestore from components.
Direct mesh changes based on v1.3 polyTopoChange syntax.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
const DynamicList< label > & region() const
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const labelUList &patchMap, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
static void removeFiles(const polyMesh &mesh)
Helper: remove all procAddressing files from mesh instance.
static void removeFiles(const polyMesh &)
Helper: remove all sets files from mesh instance.
Definition topoSet.C:693
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
#define defineTemplateTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information for templates, useful.
Definition className.H:158
rDeltaTY field()
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
engineTime & runTime
Required Variables.
Foam::PtrList< Foam::fvMesh > meshes(regionNames.size())
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const word dictName("faMeshDefinition")
auto & name
const pointField & points
return returnReduce(nRefine-oldNRefine, sumOp< label >())
#define WarningInFunction
Report a warning using Foam::Warning.
void set(List< bool > &bools, const labelUList &locations)
Set the listed locations (assign 'true').
Definition BitOps.C:35
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition meshTools.C:196
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, const bool syncPar=true, const bool readOldTime=false)
Read Geometric fields of templated type.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< word > wordList
List of word.
Definition fileName.H:60
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values within a list.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
List< labelListList > labelListListList
List of labelListList.
Definition labelList.H:41
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition triad.C:373
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
void transformList(const tensor &rotTensor, UList< T > &field)
Inplace transform a list of elements.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
wordList patchNames(nPatches)
dictionary dict
IOobject dictIO
Foam::argList args(argc, argv)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299