Loading...
Searching...
No Matches
stitchMesh.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) 2017-2024 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 stitchMesh
29
30Group
31 grpMeshManipulationUtilities
32
33Description
34 'Stitches' a mesh.
35
36 Takes a mesh and two patches and merges the faces on the two patches
37 (if geometrically possible) so the faces become internal.
38
39 Can do
40 - 'perfect' match: faces and points on patches align exactly. Order might
41 be different though.
42 - 'integral' match: where the surfaces on both patches exactly
43 match but the individual faces not
44 - 'partial' match: where the non-overlapping part of the surface remains
45 in the respective patch.
46
47 Note : Is just a front-end to perfectInterface/slidingInterface.
48
49 Comparable to running a meshModifier of the form
50 (if masterPatch is called "M" and slavePatch "S"):
51
52 \verbatim
53 couple
54 {
55 type slidingInterface;
56 masterFaceZoneName MSMasterZone
57 slaveFaceZoneName MSSlaveZone
58 cutPointZoneName MSCutPointZone
59 cutFaceZoneName MSCutFaceZone
60 masterPatchName M;
61 slavePatchName S;
62 typeOfMatch partial or integral
63 }
64 \endverbatim
65
66
67\*---------------------------------------------------------------------------*/
68
69#include "fvCFD.H"
70#include "polyTopoChanger.H"
71#include "mapPolyMesh.H"
72#include "slidingInterface.H"
73#include "perfectInterface.H"
74#include "IOobjectList.H"
75#include "ReadFields.H"
76
77// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78
79// Checks whether patch present and non-zero
80bool checkPatch(const polyBoundaryMesh& bMesh, const word& name)
81{
82 const label patchi = bMesh.findPatchID(name);
83
84 if (patchi == -1)
85 {
86 Info<< "No patch " << name << " in mesh" << nl
87 << "Known patches: " << bMesh.names() << endl;
88
89 return false;
90 }
91
92 if (bMesh[patchi].empty())
93 {
94 Info<< "Patch " << name << " has zero size" << nl;
95
96 return false;
97 }
98
99 return true;
100}
101
102
103int main(int argc, char *argv[])
104{
105 argList::addNote
106 (
107 "Merge the faces on specified patches (if geometrically possible)"
108 " so that the\n"
109 "faces become internal.\n"
110 "This utility can be called without arguments (uses stitchMeshDict)"
111 " or with\n"
112 "two arguments (master/slave patch names)."
113 );
114
115 argList::noParallel();
116 argList::noFunctionObjects(); // Never use function objects
117
118 #include "addOverwriteOption.H"
119 #include "addRegionOption.H"
120
121 argList::addOption("dict", "file", "Alternative stitchMeshDict");
122
123 argList::addBoolOption
124 (
125 "integral",
126 "Couple integral master/slave patches (2 argument mode: default)"
127 );
128 argList::addBoolOption
129 (
130 "partial",
131 "Couple partially overlapping master/slave patches (2 argument mode)"
132 );
133 argList::addBoolOption
134 (
135 "perfect",
136 "Couple perfectly aligned master/slave patches (2 argument mode)"
137 );
138 argList::addBoolOption
139 (
140 "intermediate",
141 "Write intermediate stages, not just the final result"
142 );
143 argList::addOption
144 (
145 "toleranceDict",
146 "file",
147 "Dictionary file with tolerances"
148 );
149
150 // The arguments are optional (non-mandatory) when using dictionary mode
151 argList::noMandatoryArgs();
152 argList::addArgument
153 (
154 "master",
155 "The master patch name (non-dictionary mode)"
156 );
157 argList::addArgument
158 (
159 "slave",
160 "The slave patch name (non-dictionary mode)"
161 );
162
163 #include "setRootCase.H"
164
165 // We now handle checking args and general sanity etc.
166 const bool useCommandArgs = (args.size() > 1);
167
168 if (useCommandArgs)
169 {
170 if (args.found("dict"))
171 {
173 << "Cannot specify both dictionary and command-line arguments"
174 << nl
175 << endl;
176 }
177
178 // If we have arguments - we require all arguments!
179 if (!args.check(true, false))
180 {
182 }
183 }
184 else
185 {
186 // Carp about inapplicable options
187
188 if (args.found("integral"))
189 {
191 << "Only specify -integral with command-line arguments"
192 << endl;
193 }
194
195 if (args.found("partial"))
196 {
198 << "Only specify -partial with command-line arguments"
199 << endl;
200 }
201
202 if (args.found("perfect"))
203 {
205 << "Only specify -perfect with command-line arguments"
206 << endl;
207 }
208 }
209
210 #include "createTime.H"
211 #include "createNamedMesh.H"
212
213 const word oldInstance = mesh.pointsInstance();
214
215 const bool intermediate = args.found("intermediate");
216 const bool overwrite = args.found("overwrite");
217
218 const word dictName("stitchMeshDict");
219
220 // A validated input dictionary
221 dictionary validatedDict;
222
223 if (useCommandArgs)
224 {
225 // Command argument driven:
226 const int integralCover = args.found("integral");
227 const int partialCover = args.found("partial");
228 const int perfectCover = args.found("perfect");
229
230 if ((integralCover + partialCover + perfectCover) > 1)
231 {
233 << "Can only specify one of -integral | -partial | -perfect."
234 << nl
235 << "Use perfect match option if the patches perfectly align"
236 << " (both vertex positions and face centres)" << endl
237 << exit(FatalError);
238 }
239
240 // Patch names
241 const word masterPatchName(args[1]);
242 const word slavePatchName(args[2]);
243
244 // Patch names
245 Info<< " " << masterPatchName
246 << " / " << slavePatchName << nl;
247
248 // Bail out if either patch has problems
249 if
250 (
251 !checkPatch(mesh.boundaryMesh(), masterPatchName)
252 || !checkPatch(mesh.boundaryMesh(), slavePatchName)
253 )
254 {
256 << "Cannot continue"
257 << exit(FatalError);
258
259 return 1;
260 }
261
262 // Input was validated
264
265 if (perfectCover)
266 {
267 dict.add("match", word("perfect"));
268 }
269 else if (partialCover)
270 {
271 dict.add
272 (
273 "match",
274 slidingInterface::typeOfMatchNames[slidingInterface::PARTIAL]
275 );
276 }
277 else
278 {
279 dict.add
280 (
281 "match",
282 slidingInterface::typeOfMatchNames[slidingInterface::INTEGRAL]
283 );
284 }
285
286 // Patch names
287 dict.add("master", masterPatchName);
288 dict.add("slave", slavePatchName);
289
290 validatedDict.add("stitchMesh", dict);
291 }
292 else
293 {
294 // dictionary-driven:
295
297
298 Info<< "Reading " << dictIO.name() << flush;
299
300 IOdictionary stitchDict(dictIO);
301
302 Info<< " with " << stitchDict.size() << " entries" << nl;
303
304 // Suppress duplicate names
305 wordHashSet requestedPatches;
306
307 for (const entry& dEntry : stitchDict)
308 {
309 if (!dEntry.isDict())
310 {
311 Info<< "Ignoring non-dictionary entry: "
312 << dEntry.keyword() << nl;
313 continue;
314 }
315
316 const keyType& key = dEntry.keyword();
317 const dictionary& dict = dEntry.dict();
318
319 // Match type
320 word matchName;
321 if (dict.readIfPresent("match", matchName))
322 {
323 if
324 (
325 matchName != "perfect"
326 && !slidingInterface::typeOfMatchNames.found(matchName)
327 )
328 {
329 Info<< "Error: unknown match type - " << matchName
330 << " should be one of "
331 << slidingInterface::typeOfMatchNames.toc() << nl;
332 continue;
333 }
334 }
335
336 // Patch names
337 const word masterPatchName(dict.get<word>("master"));
338 const word slavePatchName(dict.get<word>("slave"));
339
340 // Patch names
341 Info<< " " << masterPatchName
342 << " / " << slavePatchName << nl;
343
344 if (!requestedPatches.insert(masterPatchName))
345 {
346 Info<< "Error: patch specified multiple times - "
347 << masterPatchName << nl;
348 continue;
349 }
350
351 if (!requestedPatches.insert(slavePatchName))
352 {
353 Info<< "Error: patch specified multiple times - "
354 << slavePatchName << nl;
355
356 requestedPatches.erase(masterPatchName);
357 continue;
358 }
359
360 // Bail out if either patch has problems
361 if
362 (
363 !checkPatch(mesh.boundaryMesh(), masterPatchName)
364 || !checkPatch(mesh.boundaryMesh(), slavePatchName)
365 )
366 {
367 requestedPatches.erase(masterPatchName);
368 requestedPatches.erase(slavePatchName);
369 continue;
370 }
371
372 // Input was validated
373
374 validatedDict.add(key, dict);
375 }
376 }
377
378 const label nActions = validatedDict.size();
379
380 Info<< nl << nActions << " validated actions" << endl;
381
382 if (!nActions)
383 {
384 Info<<"\nStopping" << nl << endl;
385 return 1;
386 }
387
388
389 // ------------------------------------------
390 // This is where the real work begins
391
392 // set up the tolerances for the sliding mesh
393 dictionary slidingTolerances;
394 if (args.found("toleranceDict"))
395 {
396 IOdictionary toleranceFile
397 (
398 IOobject
399 (
400 args["toleranceDict"],
401 runTime.constant(),
402 mesh,
403 IOobject::MUST_READ_IF_MODIFIED,
404 IOobject::NO_WRITE
405 )
406 );
407 slidingTolerances += toleranceFile;
408 }
409
410
411 // Search for list of objects for this time
412 IOobjectList objects(mesh, runTime.timeName());
413
414 // Read all current fvFields so they will get mapped
415 Info<< "Reading all current volfields" << endl;
416 PtrList<volScalarField> volScalarFields;
417 ReadFields(mesh, objects, volScalarFields);
418
419 PtrList<volVectorField> volVectorFields;
420 ReadFields(mesh, objects, volVectorFields);
421
422 PtrList<volSphericalTensorField> volSphericalTensorFields;
423 ReadFields(mesh, objects, volSphericalTensorFields);
424
425 PtrList<volSymmTensorField> volSymmTensorFields;
426 ReadFields(mesh, objects, volSymmTensorFields);
427
428 PtrList<volTensorField> volTensorFields;
429 ReadFields(mesh, objects, volTensorFields);
430
431 //- Uncomment if you want to interpolate surface fields (usually a bad idea)
432 //Info<< "Reading all current surfaceFields" << endl;
433 //PtrList<surfaceScalarField> surfaceScalarFields;
434 //ReadFields(mesh, objects, surfaceScalarFields);
435 //
436 //PtrList<surfaceVectorField> surfaceVectorFields;
437 //ReadFields(mesh, objects, surfaceVectorFields);
438 //
439 //PtrList<surfaceTensorField> surfaceTensorFields;
440 //ReadFields(mesh, objects, surfaceTensorFields);
441
442 // More precision (for points data)
443 IOstream::minPrecision(10);
444
445 polyTopoChanger stitcher(mesh, IOobject::NO_READ);
446
447 // Step through the topology changes
448 label actioni = 0;
449 for (const entry& dEntry : validatedDict)
450 {
451 const dictionary& dict = dEntry.dict();
452
453 // Match type
454 bool perfect = false;
455 slidingInterface::typeOfMatch matchType = slidingInterface::PARTIAL;
456
457 word matchName;
458 if (dict.readIfPresent("match", matchName))
459 {
460 if (matchName == "perfect")
461 {
462 perfect = true;
463 }
464 else
465 {
466 matchType = slidingInterface::typeOfMatchNames[matchName];
467 }
468 }
469
470 // Patch names
471 const word masterPatchName(dict.get<word>("master"));
472 const word slavePatchName(dict.get<word>("slave"));
473
474 // Zone names
475 const word mergePatchName(masterPatchName + slavePatchName);
476 const word cutZoneName(mergePatchName + "CutFaceZone");
477
478 Info<< nl << "========================================" << nl;
479
480 // Information messages
481 if (perfect)
482 {
483 Info<< "Coupling PERFECTLY aligned patches "
484 << masterPatchName << " / " << slavePatchName << nl << nl
485 << "Resulting (internal) faces in faceZone "
486 << cutZoneName << nl << nl
487 << "The patch vertices and face centres must align within a"
488 << " tolerance relative to the minimum edge length on the patch"
489 << nl << endl;
490 }
491 else if (matchType == slidingInterface::INTEGRAL)
492 {
493 Info<< "Coupling INTEGRALLY matching of patches "
494 << masterPatchName << " / " << slavePatchName << nl << nl
495 << "Resulting (internal) faces in faceZone "
496 << cutZoneName << nl << nl
497 << "The overall area covered by both patches should be"
498 << " identical!" << endl
499 << "If this is not the case use partial"
500 << nl << endl;
501 }
502 else
503 {
504 Info<< "Coupling PARTIALLY overlapping patches "
505 << masterPatchName << " / " << slavePatchName << nl
506 << "Resulting internal faces in faceZone "
507 << cutZoneName << nl
508 << "Uncovered faces remain in their patch"
509 << nl << endl;
510 }
511
512
513 // Master/slave patches
514 const polyPatch& masterPatch = mesh.boundaryMesh()[masterPatchName];
515 const polyPatch& slavePatch = mesh.boundaryMesh()[slavePatchName];
516
517 mesh.pointZones().clearAddressing();
518 mesh.faceZones().clearAddressing();
519 mesh.cellZones().clearAddressing();
520
521 // Lists of master and slave faces:
522 labelList faceIds;
523
524 // Markup master face ids
525 faceIds.resize_nocopy(masterPatch.size());
526 Foam::identity(faceIds, masterPatch.start());
527
528 stitcher.clear();
529 stitcher.setSize(1);
530
531 if (perfect)
532 {
533 // Add new (empty) zone for resulting internal faces
534 mesh.faceZones()
535 (
536 cutZoneName,
537 true // verbose
538 ).resetAddressing(std::move(faceIds), false);
539
540
541 // Add the perfect interface mesh modifier
542 stitcher.set
543 (
544 0,
545 new perfectInterface
546 (
547 "couple" + Foam::name(actioni),
548 0,
549 stitcher,
550 cutZoneName,
551 masterPatchName,
552 slavePatchName
553 )
554 );
555 }
556 else
557 {
558 mesh.pointZones()
559 (
560 mergePatchName + "CutPointZone",
561 true // verbose
562 ) = labelList();
563
564 mesh.faceZones()
565 (
566 mergePatchName + "Side0Zone",
567 true // verbose
568 ).resetAddressing(std::move(faceIds), false);
569
570 // Markup slave face ids
571 faceIds.resize_nocopy(slavePatch.size());
572 Foam::identity(faceIds, slavePatch.start());
573
574 mesh.faceZones()
575 (
576 mergePatchName + "Side1Zone",
577 true // verbose
578 ).resetAddressing(std::move(faceIds), false);
579
580 // Add empty zone for cut faces
581 mesh.faceZones()
582 (
583 cutZoneName,
584 true // verbose
585 ).resetAddressing(labelList(), false);
586
587
588 // Add the sliding interface mesh modifier
589 stitcher.set
590 (
591 0,
592 new slidingInterface
593 (
594 "couple" + Foam::name(actioni),
595 0,
596 stitcher,
597 mergePatchName + "Side0Zone",
598 mergePatchName + "Side1Zone",
599 mergePatchName + "CutPointZone",
600 cutZoneName,
601 masterPatchName,
602 slavePatchName,
603 matchType, // integral or partial
604 true // couple/decouple mode
605 )
606 );
607
608 static_cast<slidingInterface&>(stitcher[0]).setTolerances
609 (
610 slidingTolerances,
611 true
612 );
613 }
614
615 ++actioni;
616
617 // Advance time for intermediate results or only on final
618 if (!overwrite && (intermediate || actioni == nActions))
619 {
620 ++runTime;
621 }
622
623 // Execute all polyMeshModifiers
624 autoPtr<mapPolyMesh> morphMap = stitcher.changeMesh(true);
625
626 mesh.movePoints(morphMap->preMotionPoints());
627
628 // Write mesh
629 if (overwrite)
630 {
631 mesh.setInstance(oldInstance);
632 stitcher.instance() = oldInstance;
633 }
634
635 if (intermediate || actioni == nActions)
636 {
637 Info<< nl << "Writing polyMesh to time "
638 << runTime.timeName() << endl;
639
640 // Bypass runTime write (since only writes at writeTime)
641 if
642 (
643 !runTime.objectRegistry::writeObject
644 (
645 IOstreamOption
646 (
647 runTime.writeFormat(),
648 runTime.writeCompression()
649 ),
650 true
651 )
652 )
653 {
655 << "Failed writing polyMesh."
656 << exit(FatalError);
657 }
658
659 mesh.faceZones().write();
660 mesh.pointZones().write();
661 mesh.cellZones().write();
662
663 // Write fields
664 runTime.write();
665 }
666 }
667
668 Info<< "\nEnd\n" << endl;
669
670 return 0;
671}
672
673
674// ************************************************************************* //
Field reading functions for post-processing utilities.
bool found
Required Classes.
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
void exit(const int errNo=1)
Exit : can be called for any error to exit program.
Definition error.C:352
dynamicFvMesh & mesh
engineTime & runTime
Required Classes.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const word dictName("faMeshDefinition")
auto & name
constexpr auto key(const Type &t) noexcept
Helper function to return the enum value.
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.
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition HashSet.H:80
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
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...
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
PrimitivePatch< List< face >, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points).
Definition bMesh.H:39
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
Ostream & flush(Ostream &os)
Flush stream.
Definition Ostream.H:509
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
bool checkPatch(const bool allGeometry, const std::string &name, const polyMesh &mesh, const PatchType &pp, const labelUList &meshEdges, labelHashSet *pointSetPtr=nullptr, labelHashSet *badEdgesPtr=nullptr)
dictionary dict
IOobject dictIO
Foam::argList args(argc, argv)