Loading...
Searching...
No Matches
PDRMesh.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-2016 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 PDRMesh
29
30Group
31 grpMeshAdvancedUtilities
32
33Description
34 Mesh and field preparation utility for PDR type simulations.
35
36 Reads
37 - cellSet giving blockedCells
38 - faceSets giving blockedFaces and the patch they should go into
39
40 NOTE: To avoid exposing wrong fields values faceSets should include
41 faces contained in the blockedCells cellset.
42
43 - coupledFaces reads coupledFacesSet to introduces mixed-coupled
44 duplicate baffles
45
46 Subsets out the blocked cells and splits the blockedFaces and updates
47 fields.
48
49 The face splitting is done by duplicating the faces. No duplication of
50 points for now (so checkMesh will show a lot of 'duplicate face' messages)
51
52\*---------------------------------------------------------------------------*/
53
54#include "fvMeshSubset.H"
55#include "argList.H"
56#include "cellSet.H"
57#include "BitOps.H"
58#include "IOobjectList.H"
59#include "volFields.H"
60#include "mapPolyMesh.H"
61#include "faceSet.H"
62#include "cellSet.H"
63#include "pointSet.H"
64#include "syncTools.H"
65#include "ReadFields.H"
66#include "polyTopoChange.H"
67#include "polyModifyFace.H"
68#include "polyAddFace.H"
69#include "regionSplit.H"
70#include "Tuple2.H"
71#include "cyclicFvPatch.H"
72
73using namespace Foam;
74
75// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76
77void modifyOrAddFace
78(
79 polyTopoChange& meshMod,
80 const face& f,
81 const label facei,
82 const label own,
83 const bool flipFaceFlux,
84 const label newPatchi,
85 const label zoneID,
86 const bool zoneFlip,
87
88 bitSet& modifiedFace
89)
90{
91 if (modifiedFace.set(facei))
92 {
93 // First usage of face. Modify.
94 meshMod.setAction
95 (
97 (
98 f, // modified face
99 facei, // label of face
100 own, // owner
101 -1, // neighbour
102 flipFaceFlux, // face flip
103 newPatchi, // patch for face
104 false, // remove from zone
105 zoneID, // zone for face
106 zoneFlip // face flip in zone
107 )
108 );
109 }
110 else
111 {
112 // Second or more usage of face. Add.
113 meshMod.setAction
114 (
116 (
117 f, // modified face
118 own, // owner
119 -1, // neighbour
120 -1, // master point
121 -1, // master edge
122 facei, // master face
123 flipFaceFlux, // face flip
124 newPatchi, // patch for face
125 zoneID, // zone for face
126 zoneFlip // face flip in zone
127 )
128 );
129 }
130}
131
132
133template<class Type>
135(
136 const fvMeshSubset& subsetter,
137 const IOobjectList& objects,
138 const label patchi,
139 const Type& exposedValue
140)
141{
143
144 const fvMesh& baseMesh = subsetter.baseMesh();
145
146 const UPtrList<const IOobject> fieldObjects
147 (
148 objects.csorted<GeoField>()
149 );
150
151 PtrList<GeoField> subFields(fieldObjects.size());
152
153 label nFields = 0;
154 for (const IOobject& io : fieldObjects)
155 {
156 if (!nFields)
157 {
158 Info<< "Subsetting " << GeoField::typeName << " (";
159 }
160 else
161 {
162 Info<< ' ';
163 }
164 Info<< " " << io.name() << endl;
165
166 // Read unregistered
168 GeoField origField(rio, baseMesh);
169
170 subFields.set(nFields, subsetter.interpolate(origField));
171 auto& subField = subFields[nFields];
172 ++nFields;
173
174
175 // Subsetting adds 'subset' prefix. Rename field to be like original.
176 subField.rename(io.name());
177 subField.writeOpt(IOobjectOption::AUTO_WRITE);
178
179
180 // Explicitly set exposed faces (in patchi) to exposedValue.
181 if (patchi >= 0)
182 {
183 fvPatchField<Type>& fld = subField.boundaryFieldRef()[patchi];
184
185 const label newStart = fld.patch().patch().start();
186 const label oldPatchi = subsetter.patchMap()[patchi];
187
188 if (oldPatchi == -1)
189 {
190 // New patch. Reset whole value.
191 fld = exposedValue;
192 }
193 else
194 {
195 // Reset faces that originate from different patch
196 // or internal faces.
197
198 const fvPatchField<Type>& origPfld =
199 origField.boundaryField()[oldPatchi];
200
201 const label oldSize = origPfld.size();
202 const label oldStart = origPfld.patch().patch().start();
203
204 forAll(fld, j)
205 {
206 const label oldFacei = subsetter.faceMap()[newStart+j];
207
208 if (oldFacei < oldStart || oldFacei >= oldStart+oldSize)
209 {
210 fld[j] = exposedValue;
211 }
212 }
213 }
214 }
215 }
216
217 if (nFields)
218 {
219 Info<< ')' << nl;
220 }
221
222 return subFields;
223}
224
225
226template<class Type>
228(
229 const fvMeshSubset& subsetter,
230 const IOobjectList& objects,
231 const label patchi,
232 const Type& exposedValue
233)
234{
236
237 const fvMesh& baseMesh = subsetter.baseMesh();
238
239 const UPtrList<const IOobject> fieldObjects
240 (
241 objects.csorted<GeoField>()
242 );
243
244 PtrList<GeoField> subFields(fieldObjects.size());
245
246 label nFields = 0;
247 for (const IOobject& io : fieldObjects)
248 {
249 if (!nFields)
250 {
251 Info<< "Subsetting " << GeoField::typeName << " (";
252 }
253 else
254 {
255 Info<< ' ';
256 }
257 Info<< io.name();
258
259 // Read unregistered
261 GeoField origField(rio, baseMesh);
262
263 subFields.set(nFields, subsetter.interpolate(origField));
264 auto& subField = subFields[nFields];
265 ++nFields;
266
267 // Subsetting adds 'subset' prefix. Rename field to be like original.
268 subField.rename(io.name());
269 subField.writeOpt(IOobjectOption::AUTO_WRITE);
270
271
272 // Explicitly set exposed faces (in patchi) to exposedValue.
273 if (patchi >= 0)
274 {
275 fvsPatchField<Type>& fld = subField.boundaryFieldRef()[patchi];
276
277 const label newStart = fld.patch().patch().start();
278 const label oldPatchi = subsetter.patchMap()[patchi];
279
280 if (oldPatchi == -1)
281 {
282 // New patch. Reset whole value.
283 fld = exposedValue;
284 }
285 else
286 {
287 // Reset faces that originate from different patch
288 // or internal faces.
289
290 const fvsPatchField<Type>& origPfld =
291 origField.boundaryField()[oldPatchi];
292
293 const label oldSize = origPfld.size();
294 const label oldStart = origPfld.patch().patch().start();
295
296 forAll(fld, j)
297 {
298 const label oldFacei = subsetter.faceMap()[newStart+j];
299
300 if (oldFacei < oldStart || oldFacei >= oldStart+oldSize)
301 {
302 fld[j] = exposedValue;
303 }
304 }
305 }
306 }
307 }
308
309 if (nFields)
310 {
311 Info<< ')' << nl;
312 }
313
314 return subFields;
315}
316
317
318// Faces introduced into zero-sized patches don't get a value at all.
319// This is hack to set them to an initial value.
320template<class GeoField>
321void initCreatedPatches
322(
323 const fvMesh& mesh,
324 const mapPolyMesh& map,
325 const typename GeoField::value_type initValue
326)
327{
328 for (const GeoField& field : mesh.objectRegistry::csorted<GeoField>())
329 {
330 auto& fieldBf = const_cast<GeoField&>(field).boundaryFieldRef();
331
332 forAll(fieldBf, patchi)
333 {
334 if (map.oldPatchSizes()[patchi] == 0)
335 {
336 // Not mapped.
337 fieldBf[patchi] = initValue;
338
339 if (fieldBf[patchi].fixesValue())
340 {
341 fieldBf[patchi] == initValue;
342 }
343 }
344 }
345 }
346}
347
348
349template<class TopoSet>
350void subsetTopoSets
351(
352 const fvMesh& mesh,
353 const IOobjectList& objects,
354 const labelList& map,
355 const fvMesh& subMesh,
356 PtrList<TopoSet>& subSets
357)
358{
359 // Read original sets
360 PtrList<TopoSet> sets;
361 ReadFields<TopoSet>(objects, sets);
362
363 subSets.resize_null(sets.size());
364
365 forAll(sets, seti)
366 {
367 const TopoSet& set = sets[seti];
368
369 Info<< "Subsetting " << set.type() << ' ' << set.name() << endl;
370
372 subset.reserve(Foam::min(set.size(), map.size()));
373
374 // Map the data
375 forAll(map, i)
376 {
377 if (set.contains(map[i]))
378 {
379 subset.insert(i);
380 }
381 }
382
383 subSets.set
384 (
385 seti,
386 new TopoSet
387 (
388 subMesh,
389 set.name(),
390 std::move(subset),
392 )
393 );
394 }
395}
396
397
398void createCoupledBaffles
399(
400 fvMesh& mesh,
401 const labelList& coupledWantedPatch,
402 polyTopoChange& meshMod,
403 bitSet& modifiedFace
404)
405{
406 const faceZoneMesh& faceZones = mesh.faceZones();
407
408 forAll(coupledWantedPatch, facei)
409 {
410 if (coupledWantedPatch[facei] != -1)
411 {
412 const face& f = mesh.faces()[facei];
413 label zoneID = faceZones.whichZone(facei);
414 bool zoneFlip = false;
415
416 if (zoneID >= 0)
417 {
418 const faceZone& fZone = faceZones[zoneID];
419 zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
420 }
421
422 // Use owner side of face
423 modifyOrAddFace
424 (
425 meshMod,
426 f, // modified face
427 facei, // label of face
428 mesh.faceOwner()[facei], // owner
429 false, // face flip
430 coupledWantedPatch[facei], // patch for face
431 zoneID, // zone for face
432 zoneFlip, // face flip in zone
433 modifiedFace // modify or add status
434 );
435
436 if (mesh.isInternalFace(facei))
437 {
438 label zoneID = faceZones.whichZone(facei);
439 bool zoneFlip = false;
440
441 if (zoneID >= 0)
442 {
443 const faceZone& fZone = faceZones[zoneID];
444 zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
445 }
446 // Use neighbour side of face
447 modifyOrAddFace
448 (
449 meshMod,
450 f.reverseFace(), // modified face
451 facei, // label of face
452 mesh.faceNeighbour()[facei],// owner
453 false, // face flip
454 coupledWantedPatch[facei], // patch for face
455 zoneID, // zone for face
456 zoneFlip, // face flip in zone
457 modifiedFace // modify or add status
458 );
459 }
460 }
461 }
462}
463
464
465void createCyclicCoupledBaffles
466(
467 fvMesh& mesh,
468 const labelList& cyclicMasterPatch,
469 const labelList& cyclicSlavePatch,
470 polyTopoChange& meshMod,
471 bitSet& modifiedFace
472)
473{
474 const faceZoneMesh& faceZones = mesh.faceZones();
475
476 forAll(cyclicMasterPatch, facei)
477 {
478 if (cyclicMasterPatch[facei] != -1)
479 {
480 const face& f = mesh.faces()[facei];
481
482 label zoneID = faceZones.whichZone(facei);
483 bool zoneFlip = false;
484
485 if (zoneID >= 0)
486 {
487 const faceZone& fZone = faceZones[zoneID];
488 zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
489 }
490
491 modifyOrAddFace
492 (
493 meshMod,
494 f.reverseFace(), // modified face
495 facei, // label of face
496 mesh.faceNeighbour()[facei], // owner
497 false, // face flip
498 cyclicMasterPatch[facei], // patch for face
499 zoneID, // zone for face
500 zoneFlip, // face flip in zone
501 modifiedFace // modify or add
502 );
503 }
504 }
505
506 forAll(cyclicSlavePatch, facei)
507 {
508 if (cyclicSlavePatch[facei] != -1)
509 {
510 const face& f = mesh.faces()[facei];
511 if (mesh.isInternalFace(facei))
512 {
513 label zoneID = faceZones.whichZone(facei);
514 bool zoneFlip = false;
515
516 if (zoneID >= 0)
517 {
518 const faceZone& fZone = faceZones[zoneID];
519 zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
520 }
521 // Use owner side of face
522 modifyOrAddFace
523 (
524 meshMod,
525 f, // modified face
526 facei, // label of face
527 mesh.faceOwner()[facei], // owner
528 false, // face flip
529 cyclicSlavePatch[facei], // patch for face
530 zoneID, // zone for face
531 zoneFlip, // face flip in zone
532 modifiedFace // modify or add status
533 );
534 }
535 }
536 }
537}
538
539
540void createBaffles
541(
542 fvMesh& mesh,
543 const labelList& wantedPatch,
544 polyTopoChange& meshMod
545)
546{
547 const faceZoneMesh& faceZones = mesh.faceZones();
548 Info << "faceZone:createBaffle " << faceZones << endl;
549 forAll(wantedPatch, facei)
550 {
551 if (wantedPatch[facei] != -1)
552 {
553 const face& f = mesh.faces()[facei];
554
555 label zoneID = faceZones.whichZone(facei);
556 bool zoneFlip = false;
557
558 if (zoneID >= 0)
559 {
560 const faceZone& fZone = faceZones[zoneID];
561 zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
562 }
563
564 meshMod.setAction
565 (
567 (
568 f, // modified face
569 facei, // label of face
570 mesh.faceOwner()[facei], // owner
571 -1, // neighbour
572 false, // face flip
573 wantedPatch[facei], // patch for face
574 false, // remove from zone
575 zoneID, // zone for face
576 zoneFlip // face flip in zone
577 )
578 );
579
580 if (mesh.isInternalFace(facei))
581 {
582 label zoneID = faceZones.whichZone(facei);
583 bool zoneFlip = false;
584
585 if (zoneID >= 0)
586 {
587 const faceZone& fZone = faceZones[zoneID];
588 zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
589 }
590
591 meshMod.setAction
592 (
594 (
595 f.reverseFace(), // modified face
596 mesh.faceNeighbour()[facei],// owner
597 -1, // neighbour
598 -1, // masterPointID
599 -1, // masterEdgeID
600 facei, // masterFaceID,
601 false, // face flip
602 wantedPatch[facei], // patch for face
603 zoneID, // zone for face
604 zoneFlip // face flip in zone
605 )
606 );
607 }
608 }
609 }
610}
611
612
613// Wrapper around find patch. Also makes sure same patch in parallel.
614label findPatch(const polyBoundaryMesh& patches, const word& patchName)
615{
616 const label patchi = patches.findPatchID(patchName);
617
618 if (patchi == -1)
619 {
621 << "Illegal patch " << patchName
622 << nl << "Valid patches are " << patches.names()
623 << exit(FatalError);
624 }
625
626 // Check same patch for all procs
627 {
628 const label newPatchi = returnReduce(patchi, minOp<label>());
629
630 if (newPatchi != patchi)
631 {
633 << "Patch " << patchName
634 << " should have the same patch index on all processors." << nl
635 << "On my processor it has index " << patchi
636 << " ; on some other processor it has index " << newPatchi
637 << exit(FatalError);
638 }
639 }
640 return patchi;
641}
642
643
644// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
645
646int main(int argc, char *argv[])
647{
649 (
650 "Mesh and field preparation utility for PDR type simulations."
651 );
652 #include "addOverwriteOption.H"
653
654 argList::noFunctionObjects(); // Never use function objects
655
656 #include "setRootCase.H"
657 #include "createTime.H"
658 #include "createNamedMesh.H"
659
660 // Read control dictionary
661 // ~~~~~~~~~~~~~~~~~~~~~~~
662
664 (
666 (
667 "PDRMeshDict",
668 runTime.system(),
669 mesh,
672 )
673 );
674
675 // Per faceSet the patch to put the baffles into
676 const List<Pair<word>> setsAndPatches(dict.lookup("blockedFaces"));
677
678 // Per faceSet the patch to put the coupled baffles into
679 DynamicList<FixedList<word, 3>> coupledAndPatches(10);
680
681 const dictionary& functionDicts = dict.subDict("coupledFaces");
682
683 for (const entry& dEntry : functionDicts)
684 {
685 if (!dEntry.isDict()) // Safety
686 {
687 continue;
688 }
689
690 const word& key = dEntry.keyword();
691 const dictionary& dict = dEntry.dict();
692
693 const word cyclicName = dict.get<word>("cyclicMasterPatch");
694 const word wallName = dict.get<word>("wallPatch");
695 FixedList<word, 3> nameAndType;
696 nameAndType[0] = key;
697 nameAndType[1] = wallName;
698 nameAndType[2] = cyclicName;
699 coupledAndPatches.append(nameAndType);
700 }
701
702 forAll(setsAndPatches, setI)
703 {
704 Info<< "Faces in faceSet " << setsAndPatches[setI][0]
705 << " become baffles in patch " << setsAndPatches[setI][1]
706 << endl;
707 }
708
709 forAll(coupledAndPatches, setI)
710 {
711 Info<< "Faces in faceSet " << coupledAndPatches[setI][0]
712 << " become coupled baffles in patch " << coupledAndPatches[setI][1]
713 << endl;
714 }
715
716 // All exposed faces that are not explicitly marked to be put into a patch
717 const word defaultPatch(dict.get<word>("defaultPatch"));
718
719 Info<< "Faces that get exposed become boundary faces in patch "
720 << defaultPatch << endl;
721
722 const word blockedSetName(dict.get<word>("blockedCells"));
723
724 Info<< "Reading blocked cells from cellSet " << blockedSetName
725 << endl;
726
727 const bool overwrite = args.found("overwrite");
728
729
730 // Read faceSets, lookup patches
731 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
732
733 // Check that face sets don't have coincident faces
734 labelList wantedPatch(mesh.nFaces(), -1);
735 forAll(setsAndPatches, setI)
736 {
737 faceSet fSet(mesh, setsAndPatches[setI][0]);
738
739 label patchi = findPatch
740 (
741 mesh.boundaryMesh(),
742 setsAndPatches[setI][1]
743 );
744
745 for (const label facei : fSet)
746 {
747 if (wantedPatch[facei] != -1)
748 {
750 << "Face " << facei
751 << " is in faceSet " << setsAndPatches[setI][0]
752 << " destined for patch " << setsAndPatches[setI][1]
753 << " but also in patch " << wantedPatch[facei]
754 << exit(FatalError);
755 }
756 wantedPatch[facei] = patchi;
757 }
758 }
759
760 // Per face the patch for coupled baffle or -1.
761 labelList coupledWantedPatch(mesh.nFaces(), -1);
762 labelList cyclicWantedPatch_half0(mesh.nFaces(), -1);
763 labelList cyclicWantedPatch_half1(mesh.nFaces(), -1);
764
765 forAll(coupledAndPatches, setI)
766 {
767 const polyBoundaryMesh& patches = mesh.boundaryMesh();
768 const label cyclicId =
769 findPatch(patches, coupledAndPatches[setI][2]);
770
771 const label cyclicSlaveId = findPatch
772 (
773 patches,
775 (
776 mesh.boundary()[cyclicId]
777 ).neighbFvPatch().name()
778 );
779
780 faceSet fSet(mesh, coupledAndPatches[setI][0]);
781 label patchi = findPatch(patches, coupledAndPatches[setI][1]);
782
783 for (const label facei : fSet)
784 {
785 if (coupledWantedPatch[facei] != -1)
786 {
788 << "Face " << facei
789 << " is in faceSet " << coupledAndPatches[setI][0]
790 << " destined for patch " << coupledAndPatches[setI][1]
791 << " but also in patch " << coupledWantedPatch[facei]
792 << exit(FatalError);
793 }
794
795 coupledWantedPatch[facei] = patchi;
796 cyclicWantedPatch_half0[facei] = cyclicId;
797 cyclicWantedPatch_half1[facei] = cyclicSlaveId;
798 }
799 }
800
801 // Exposed faces patch
802 label defaultPatchi = findPatch(mesh.boundaryMesh(), defaultPatch);
803
804
805 //
806 // Removing blockedCells (like subsetMesh)
807 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
808 //
809
810 // Mesh subsetting engine
811 fvMeshSubset subsetter
812 (
813 mesh,
815 (
816 mesh.nCells(),
817 cellSet(mesh, blockedSetName), // Blocked cells as labelHashSet
818 false // on=false: invert logic => retain the unblocked cells
819 ),
820 defaultPatchi,
821 true
822 );
823
824
825 // Subset wantedPatch. Note that might also include boundary faces
826 // that have been exposed by subsetting.
827 wantedPatch = IndirectList<label>(wantedPatch, subsetter.faceMap())();
828
829 coupledWantedPatch = IndirectList<label>
830 (
831 coupledWantedPatch,
832 subsetter.faceMap()
833 )();
834
835 cyclicWantedPatch_half0 = IndirectList<label>
836 (
837 cyclicWantedPatch_half0,
838 subsetter.faceMap()
839 )();
840
841 cyclicWantedPatch_half1 = IndirectList<label>
842 (
843 cyclicWantedPatch_half1,
844 subsetter.faceMap()
845 )();
846
847 // Read all fields in time and constant directories
848 IOobjectList objects(mesh, runTime.timeName());
849 {
850 IOobjectList timeObjects(mesh, mesh.facesInstance());
851
852 // Transfer specific types
853 forAllIters(timeObjects, iter)
854 {
855 autoPtr<IOobject> objPtr(timeObjects.remove(iter));
856 const auto& obj = *objPtr;
857
858 if
859 (
860 Foam::fieldTypes::is_volume(obj.headerClassName())
861 || Foam::fieldTypes::is_surface(obj.headerClassName())
862 )
863 {
864 objects.add(objPtr);
865 }
866 }
867 }
868
869 // Read vol fields and subset.
870 PtrList<volScalarField> scalarFlds
871 (
872 subsetVolFields<scalar>
873 (
874 subsetter,
875 objects,
876 defaultPatchi,
877 scalar(Zero)
878 )
879 );
880
881 PtrList<volVectorField> vectorFlds
882 (
883 subsetVolFields<vector>
884 (
885 subsetter,
886 objects,
887 defaultPatchi,
888 vector(Zero)
889 )
890 );
891
893 (
894 subsetVolFields<sphericalTensor>
895 (
896 subsetter,
897 objects,
898 defaultPatchi,
900 )
901 );
902
903 PtrList<volSymmTensorField> symmTensorFlds
904 (
905 subsetVolFields<symmTensor>
906 (
907 subsetter,
908 objects,
909 defaultPatchi,
911 )
912 );
913
914 PtrList<volTensorField> tensorFlds
915 (
916 subsetVolFields<tensor>
917 (
918 subsetter,
919 objects,
920 defaultPatchi,
921 tensor(Zero)
922 )
923 );
924
925 // Read surface fields and subset.
926 PtrList<surfaceScalarField> surfScalarFlds
927 (
928 subsetSurfaceFields<scalar>
929 (
930 subsetter,
931 objects,
932 defaultPatchi,
933 scalar(Zero)
934 )
935 );
936
937 PtrList<surfaceVectorField> surfVectorFlds
938 (
939 subsetSurfaceFields<vector>
940 (
941 subsetter,
942 objects,
943 defaultPatchi,
944 vector(Zero)
945 )
946 );
947
948 PtrList<surfaceSphericalTensorField> surfSphericalTensorFlds
949 (
950 subsetSurfaceFields<sphericalTensor>
951 (
952 subsetter,
953 objects,
954 defaultPatchi,
956 )
957 );
958
959 PtrList<surfaceSymmTensorField> surfSymmTensorFlds
960 (
961 subsetSurfaceFields<symmTensor>
962 (
963 subsetter,
964 objects,
965 defaultPatchi,
967 )
968 );
969
970 PtrList<surfaceTensorField> surfTensorFlds
971 (
972 subsetSurfaceFields<tensor>
973 (
974 subsetter,
975 objects,
976 defaultPatchi,
977 tensor(Zero)
978 )
979 );
980
981
982 // Set handling
983 PtrList<cellSet> cellSets;
984 PtrList<faceSet> faceSets;
985 PtrList<pointSet> pointSets;
986 {
987 IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
988 subsetTopoSets
989 (
990 mesh,
991 objects,
992 subsetter.cellMap(),
993 subsetter.subMesh(),
994 cellSets
995 );
996 subsetTopoSets
997 (
998 mesh,
999 objects,
1000 subsetter.faceMap(),
1001 subsetter.subMesh(),
1002 faceSets
1003 );
1004 subsetTopoSets
1005 (
1006 mesh,
1007 objects,
1008 subsetter.pointMap(),
1009 subsetter.subMesh(),
1010 pointSets
1011 );
1012 }
1013
1014
1015 if (!overwrite)
1016 {
1017 ++runTime;
1018 }
1019
1020 Info<< "Writing mesh without blockedCells to time "
1021 << runTime.value() << endl;
1022
1023 subsetter.subMesh().write();
1024
1025
1026 //
1027 // Splitting blockedFaces
1028 // ~~~~~~~~~~~~~~~~~~~~~~
1029 //
1030
1031 // Synchronize wantedPatch across coupled patches.
1033 (
1034 subsetter.subMesh(),
1035 wantedPatch,
1037 );
1038
1039 // Synchronize coupledWantedPatch across coupled patches.
1041 (
1042 subsetter.subMesh(),
1043 coupledWantedPatch,
1045 );
1046
1047 // Synchronize cyclicWantedPatch across coupled patches.
1049 (
1050 subsetter.subMesh(),
1051 cyclicWantedPatch_half0,
1053 );
1054
1055 // Synchronize cyclicWantedPatch across coupled patches.
1057 (
1058 subsetter.subMesh(),
1059 cyclicWantedPatch_half1,
1061 );
1062
1063 // Topochange container
1064 polyTopoChange meshMod(subsetter.subMesh());
1065
1066
1067 // Whether first use of face (modify) or consecutive (add)
1068 bitSet modifiedFace(mesh.nFaces());
1069
1070 // Create coupled wall-side baffles
1071 createCoupledBaffles
1072 (
1073 subsetter.subMesh(),
1074 coupledWantedPatch,
1075 meshMod,
1076 modifiedFace
1077 );
1078
1079 // Create coupled master/slave cyclic baffles
1080 createCyclicCoupledBaffles
1081 (
1082 subsetter.subMesh(),
1083 cyclicWantedPatch_half0,
1084 cyclicWantedPatch_half1,
1085 meshMod,
1086 modifiedFace
1087 );
1088
1089 // Create wall baffles
1090 createBaffles
1091 (
1092 subsetter.subMesh(),
1093 wantedPatch,
1094 meshMod
1095 );
1096
1097 if (!overwrite)
1098 {
1099 ++runTime;
1100 }
1101
1102 // Change the mesh. Change points directly (no inflation).
1103 autoPtr<mapPolyMesh> mapPtr =
1104 meshMod.changeMesh(subsetter.subMesh(), false);
1105 mapPolyMesh& map = *mapPtr;
1106
1107 // Update fields
1108 subsetter.subMesh().updateMesh(map);
1109
1110 // Fix faces that get mapped to zero-sized patches (these don't get any
1111 // value)
1112 initCreatedPatches<volScalarField>
1113 (
1114 subsetter.subMesh(),
1115 map,
1116 Zero
1117 );
1118 initCreatedPatches<volVectorField>
1119 (
1120 subsetter.subMesh(),
1121 map,
1122 Zero
1123 );
1124 initCreatedPatches<volSphericalTensorField>
1125 (
1126 subsetter.subMesh(),
1127 map,
1128 Zero
1129 );
1130 initCreatedPatches<volSymmTensorField>
1131 (
1132 subsetter.subMesh(),
1133 map,
1134 Zero
1135 );
1136 initCreatedPatches<volTensorField>
1137 (
1138 subsetter.subMesh(),
1139 map,
1140 Zero
1141 );
1142
1143 initCreatedPatches<surfaceScalarField>
1144 (
1145 subsetter.subMesh(),
1146 map,
1147 Zero
1148 );
1149 initCreatedPatches<surfaceVectorField>
1150 (
1151 subsetter.subMesh(),
1152 map,
1153 Zero
1154 );
1155 initCreatedPatches<surfaceSphericalTensorField>
1156 (
1157 subsetter.subMesh(),
1158 map,
1159 Zero
1160 );
1161 initCreatedPatches<surfaceSymmTensorField>
1162 (
1163 subsetter.subMesh(),
1164 map,
1165 Zero
1166 );
1167 initCreatedPatches<surfaceTensorField>
1168 (
1169 subsetter.subMesh(),
1170 map,
1171 Zero
1172 );
1173
1174 // Update numbering of topoSets
1175 topoSet::updateMesh(subsetter.subMesh().facesInstance(), map, cellSets);
1176 topoSet::updateMesh(subsetter.subMesh().facesInstance(), map, faceSets);
1177 topoSet::updateMesh(subsetter.subMesh().facesInstance(), map, pointSets);
1178
1179
1180 // Move mesh (since morphing might not do this)
1181 if (map.hasMotionPoints())
1182 {
1183 subsetter.subMesh().movePoints(map.preMotionPoints());
1184 }
1185
1186 Info<< "Writing mesh with split blockedFaces to time " << runTime.value()
1187 << endl;
1188
1189 subsetter.subMesh().write();
1190
1191
1192 //
1193 // Removing inaccessible regions
1194 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1195 //
1196
1197 // Determine connected regions. regionSplit is the labelList with the
1198 // region per cell.
1199 regionSplit cellRegion(subsetter.subMesh());
1200
1201 if (cellRegion.nRegions() > 1)
1202 {
1204 << "Removing blocked faces and cells created "
1205 << cellRegion.nRegions()
1206 << " regions that are not connected via a face." << nl
1207 << " This is not supported in solvers." << nl
1208 << " Use" << nl << nl
1209 << " splitMeshRegions <root> <case> -largestOnly" << nl << nl
1210 << " to extract a single region of the mesh." << nl
1211 << " This mesh will be written to a new timedirectory"
1212 << " so might have to be moved back to constant/" << nl
1213 << endl;
1214
1215 const word startFrom(runTime.controlDict().get<word>("startFrom"));
1216
1217 if (startFrom != "latestTime")
1218 {
1220 << "To run splitMeshRegions please set your"
1221 << " startFrom entry to latestTime" << endl;
1222 }
1223 }
1224
1225 Info<< "\nEnd\n" << endl;
1226
1227 return 0;
1228}
1229
1230
1231// ************************************************************************* //
Field reading functions for post-processing utilities.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
Generic GeometricField class.
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,...
UPtrList< const IOobject > csorted() const
The sorted list of IOobjects with headerClassName == Type::typeName.
bool add(std::unique_ptr< IOobject > &&objectPtr)
Move insert IOobject into the list.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A List with indirect addressing.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
void resize_null(const label newLen)
Set the addressed list to the given size, deleting all existing entries. Afterwards the list contains...
Definition PtrListI.H:113
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition ZoneMesh.C:410
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
Definition argList.C:562
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
A collection of cell labels.
Definition cellSet.H:50
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
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 face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
const fvMesh & baseMesh() const noexcept
Original mesh.
const labelList & faceMap() const
Return face map.
const labelList & cellMap() const
Return cell map.
static tmp< DimensionedField< Type, volMesh > > interpolate(const DimensionedField< Type, volMesh > &, const fvMesh &sMesh, const labelUList &cellMap)
Map volume internal (dimensioned) field.
const fvMesh & subMesh() const
Return reference to subset mesh.
const labelList & patchMap() const
Return patch map.
const labelList & pointMap() const
Return point map.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition fvMesh.C:960
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition fvMesh.C:1068
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition fvMesh.C:884
const fvPatch & patch() const noexcept
Return the patch.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
const polyPatch & patch() const noexcept
Return the polyPatch.
Definition fvPatch.H:202
const fvPatch & patch() const noexcept
Return the patch.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const pointField & preMotionPoints() const noexcept
Pre-motion point positions.
bool hasMotionPoints() const noexcept
Has valid preMotionPoints?
const labelList & oldPatchSizes() const noexcept
Return list of the old patch sizes.
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,...
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition polyMesh.C:844
Class describing modification of a face.
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition polyPatch.H:446
Direct mesh changes based on v1.3 polyTopoChange syntax.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
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.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const bool parRun=UPstream::parRun())
Synchronize values on all mesh faces.
Definition syncTools.H:465
Tensor of scalars, i.e. Tensor<scalar>.
virtual void updateMesh(const mapPolyMesh &morphMap)
Update any stored data for new labels. Not implemented.
Definition topoSet.C:687
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
rDeltaTY field()
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
engineTime & runTime
Required Classes.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
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
bitSet create(const label n, const labelHashSet &locations, const bool on=true)
Create a bitSet with length n with the specified on locations.
Definition BitOps.C:233
bool is_volume(const word &clsName)
Test if the class name appears to be a volume field.
Definition volFields.C:165
bool is_surface(const word &clsName)
Test if the class name appears to be a surface field.
constexpr auto key(const Type &t) noexcept
Helper function to return the enum value.
Namespace for OpenFOAM.
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 & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
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
messageStream Info
Information stream (stdout output on master, null elsewhere).
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
List< T > subset(const BoolListType &select, const UList< T > &input, const bool invert=false)
Extract elements of the input list when select is true.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars, i.e. SphericalTensor<scalar>.
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition symmTensor.H:55
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
dictionary dict
Foam::argList args(argc, argv)
Us boundaryFieldRef().evaluateCoupled< coupledFaPatch >()
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition stdFoam.H:214