faMeshPatches.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) 2017 Wikki Ltd
9 Copyright (C) 2021 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
27\*---------------------------------------------------------------------------*/
28
29#include "faMesh.H"
30#include "faPatchData.H"
31#include "processorPolyPatch.H"
32#include "processorFaPatch.H"
33#include "foamVtkLineWriter.H"
34
35// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
36
38(
39 faPatchList& plist,
40 const bool validBoundary
41)
42{
43 if (!boundary().empty())
44 {
46 << "boundary already exists"
47 << abort(FatalError);
48 }
49
50 globalMeshDataPtr_.reset(nullptr);
51
52 boundary_.transfer(plist);
53
54 setPrimitiveMeshData();
55
56 if (validBoundary)
57 {
58 boundary_.checkDefinition();
59 }
60}
61
62
64(
65 const List<faPatch*>& p,
66 const bool validBoundary
67)
68{
69 // Acquire ownership of the pointers
70 faPatchList plist(const_cast<List<faPatch*>&>(p));
71
72 addFaPatches(plist, validBoundary);
73}
74
75
76Foam::faPatchList Foam::faMesh::createOnePatch
77(
78 const word& patchName,
79 const word& patchType
80) const
81{
82 dictionary onePatchDict;
83 if (!patchName.empty())
84 {
85 onePatchDict.add("name", patchName);
86 }
87 if (!patchType.empty())
88 {
89 onePatchDict.add("type", patchType);
90 }
91
92 return createPatchList
93 (
95 word::null, // Name for empty patch placeholder
96 &onePatchDict // Definitions for defaultPatch
97 );
98}
99
100
101Foam::faPatchList Foam::faMesh::createPatchList
102(
103 const dictionary& bndDict,
104 const word& emptyPatchName,
105 const dictionary* defaultPatchDefinition
106) const
107{
108 const polyBoundaryMesh& pbm = mesh().boundaryMesh();
109
110 // Transcribe into patch definitions
111 DynamicList<faPatchData> faPatchDefs(bndDict.size() + 8);
112 for (const entry& dEntry : bndDict)
113 {
114 if (!dEntry.isDict())
115 {
117 << "Not a dictionary entry: " << dEntry.name() << nl;
118 continue;
119 }
120 const dictionary& patchDict = dEntry.dict();
121
122 // Add entry
123 faPatchDefs.append(faPatchData());
124
125 auto& patchDef = faPatchDefs.last();
126 patchDef.name_ = dEntry.keyword();
127 patchDef.type_ = patchDict.get<word>("type");
128
129 word patchName;
130
131 // Optional: ownerPolyPatch
132 if (patchDict.readIfPresent("ownerPolyPatch", patchName))
133 {
134 patchDef.ownerPolyPatchId_ = pbm.findPatchID(patchName);
135 if (patchDef.ownerPolyPatchId_ < 0)
136 {
138 << "ownerPolyPatch " << patchName << " not found"
139 << exit(FatalError);
140 }
141 }
142
143 // Mandatory: neighbourPolyPatch
144 patchDict.readEntry("neighbourPolyPatch", patchName);
145 {
146 patchDef.neighPolyPatchId_ = pbm.findPatchID(patchName);
147 if (patchDef.neighPolyPatchId_ < 0)
148 {
150 << "neighbourPolyPatch " << patchName << " not found"
151 << exit(FatalError);
152 }
153 }
154 }
155
156 // Additional empty placeholder patch?
157 if (!emptyPatchName.empty())
158 {
159 faPatchDefs.append(faPatchData());
160
161 auto& patchDef = faPatchDefs.last();
162 patchDef.name_ = emptyPatchName;
163 patchDef.type_ = "empty";
164 }
165
166 label nWarnUndefinedPatch(5);
167
168 // Placeholder for any undefined edges
169 const label undefPatchIndex = faPatchDefs.size();
170 {
171 faPatchDefs.append(faPatchData());
172
173 auto& patchDef = faPatchDefs.last();
174 patchDef.name_ = "undefined";
175 patchDef.type_ = "patch";
176
177 if (defaultPatchDefinition)
178 {
179 if
180 (
181 (*defaultPatchDefinition).readIfPresent("name", patchDef.name_)
182 )
183 {
184 // Suppress warnings if defaultPatch name was specified
185 // - probably means we want to use this mechanism
186 nWarnUndefinedPatch = 0;
187 }
188 (*defaultPatchDefinition).readIfPresent("type", patchDef.type_);
189 }
190 }
191
192 // ----------------------------------------------------------------------
193
194 // Get edge connections (in globally consistent ordering)
195 List<Pair<patchTuple>> bndEdgeConnections
196 (
197 this->getBoundaryEdgeConnections()
198 );
199
200 // Update accordingly
201 this->setBoundaryConnections(bndEdgeConnections);
202
203
204 // Lookup of patchDef for each connection. Initially all -1 (unvisited)
205 labelList patchDefLookup(bndEdgeConnections.size(), -1);
206
207 Map<labelHashSet> procConnections;
208 labelHashSet patchDefsUsed;
209
210 label nBadEdges(0);
211 labelHashSet badEdges(2*bndEdgeConnections.size());
212
213 forAll(bndEdgeConnections, connecti)
214 {
215 const Pair<patchTuple>& connection = bndEdgeConnections[connecti];
216 const auto& a = connection.first();
217 const auto& b = connection.second();
218
219 edge patchPair;
220
221 if (a.is_finiteArea())
222 {
223 if (b.is_finiteArea())
224 {
225 // A processor-processor connection
226 if (a.procNo() == b.procNo())
227 {
229 << "Processor-processor addressing error:" << nl
230 << "Both connections have the same processor: "
231 << a.procNo() << nl
232 << abort(FatalError);
233 }
234 else if (a.is_localProc())
235 {
236 procConnections(b.procNo()).insert(connecti);
237 }
238 else
239 {
240 procConnections(a.procNo()).insert(connecti);
241 }
242
243 continue;
244 }
245 else if (a.is_localProc())
246 {
247 patchPair.first() = a.realPatchi();
248 patchPair.second() = b.realPatchi();
249 }
250 }
251 else if (b.is_finiteArea() && b.is_localProc())
252 {
253 patchPair.first() = b.realPatchi();
254 patchPair.second() = a.realPatchi();
255 }
256
257
258 // Find first definition with a matching neighbour and
259 // possibly with a matching owner.
260
261 label bestPatchDefi = -1;
262
263 const label nPatchDefs = (patchPair.valid() ? faPatchDefs.size() : 0);
264
265 for (label patchDefi = 0; patchDefi < nPatchDefs; ++patchDefi)
266 {
267 const int match = faPatchDefs[patchDefi].matchPatchPair(patchPair);
268 if (match == 3)
269 {
270 // Exact match (owner/neighbour) - done!
271 bestPatchDefi = patchDefi;
272 break;
273 }
274 else if (match == 2 && bestPatchDefi < 0)
275 {
276 // Match (neighbour) - keep looking for exact match
277 bestPatchDefi = patchDefi;
278 }
279 }
280
281 if (bestPatchDefi < 0)
282 {
283 bestPatchDefi = undefPatchIndex; // Missed?
284 }
285
286 patchDefLookup[connecti] = bestPatchDefi;
287 patchDefsUsed.insert(bestPatchDefi);
288 }
289
290 // Remove undefPatchIndex if not actually needed
291 if (!returnReduce(patchDefsUsed.found(undefPatchIndex), orOp<bool>()))
292 {
293 faPatchDefs.remove(undefPatchIndex);
294 }
295 else
296 {
297 patchDefsUsed.insert(undefPatchIndex); // Parallel consistency
298
299 // Report locations of undefined edges
300
301 badEdges.clear();
302 forAll(patchDefLookup, connecti)
303 {
304 if (patchDefLookup[connecti] == undefPatchIndex)
305 {
306 const auto& connection = bndEdgeConnections[connecti];
307
308 const auto& a = connection.first();
309 const auto& b = connection.second();
310
311 if (a.is_localProc() && a.is_finiteArea())
312 {
313 badEdges.insert(a.patchEdgei());
314 }
315 else if (b.is_localProc() && b.is_finiteArea())
316 {
317 badEdges.insert(b.patchEdgei());
318 }
319
320 if (badEdges.size() <= nWarnUndefinedPatch)
321 {
322 Pout<< "Undefined connection: "
323 << "(patch:" << a.realPatchi()
324 << " face:" << a.meshFacei()
325 << ") and (patch:" << b.realPatchi()
326 << " face:" << b.meshFacei() << ") connects: "
327 << pbm[a.realPatchi()].name() << " to "
328 << pbm[b.realPatchi()].name() << nl;
329 }
330 }
331 }
332
333 if ((nBadEdges = returnReduce(badEdges.size(), sumOp<label>())) > 0)
334 {
335 // Report directly as Info, not InfoInFunction
336 // since it can also be an expected result when
337 // nWarnUndefinedPatch == 0
338 Info<< "Had " << nBadEdges << '/'
339 << returnReduce(patch().nBoundaryEdges(), sumOp<label>())
340 << " undefined edge connections, added to defaultPatch: "
341 << faPatchDefs[undefPatchIndex].name_ << nl;
342
343 if (nWarnUndefinedPatch)
344 {
345 edgeList dumpEdges(patch().edges(), badEdges.sortedToc());
346
347 vtk::lineWriter writer
348 (
349 patch().localPoints(),
350 dumpEdges,
351 fileName
352 (
353 mesh().time().globalPath()
354 / ("faMesh-construct.undefEdges")
355 )
356 );
357
359
360 // CellData
363
365 << "(debug) wrote " << writer.output().name() << nl;
366 }
367 }
368 }
369
370 // Create processor-processor definitions
371 Map<label> procToDefLookup(2*procConnections.size());
372 {
373 faPatchDefs.reserve(faPatchDefs.size() + procConnections.size());
374
375 for (const label otherProci : procConnections.sortedToc())
376 {
377 const label patchDefi = faPatchDefs.size();
378 procToDefLookup.insert(otherProci, patchDefi);
379
380 // Add entry
381 faPatchDefs.append(faPatchData());
382 auto& patchDef = faPatchDefs.last();
383
384 patchDef.assign_coupled(Pstream::myProcNo(), otherProci);
385 }
386 }
387
388
389 // Extract which edges map to which patch
390 // and set edgeLabels for each faPatch
391
392 DynamicList<label> selectEdges(bndEdgeConnections.size());
393 label nOffProcessorEdges = 0;
394
395 for (const label patchDefi : patchDefsUsed.sortedToc())
396 {
397 auto& patchDef = faPatchDefs[patchDefi];
398 selectEdges.clear();
399
400 // Find the corresponding entries
401 // and extract the patchEdgeId
402
403 forAll(patchDefLookup, connecti)
404 {
405 if (patchDefLookup[connecti] == patchDefi)
406 {
407 const auto& a = bndEdgeConnections[connecti].first();
408 const auto& b = bndEdgeConnections[connecti].second();
409
410 if (a.is_localProc() && a.is_finiteArea())
411 {
412 selectEdges.append(a.patchEdgei());
413 }
414 else if (b.is_localProc() && b.is_finiteArea())
415 {
416 selectEdges.append(b.patchEdgei());
417 }
418 else
419 {
421 << "Error in programming logic" << nl
422 << abort(FatalError);
423 }
424
425 if (a.is_localProc() != b.is_localProc())
426 {
427 ++nOffProcessorEdges;
428 }
429
430 patchDefLookup[connecti] = -2; // Mark as handled
431 }
432 }
433
434 // Remove any cosmetic sorting artifacts from off-processor
435 // termination by doing using a regular sort here.
436
437 Foam::sort(selectEdges);
438 patchDef.edgeLabels_ = selectEdges;
439 }
440
441 if (debug)
442 {
443 Pout<< "Had " << nOffProcessorEdges
444 << " patch edges connected off-processor" << endl;
445
447 << "Total "
448 << returnReduce(nOffProcessorEdges, sumOp<label>())
449 << " patch edges connected off-processor" << endl;
450 }
451
452
453 // Processor patches
454 for (const label otherProci : procToDefLookup.sortedToc())
455 {
456 const label patchDefi = procToDefLookup[otherProci];
457
458 auto& patchDef = faPatchDefs[patchDefi];
459 selectEdges.clear();
460
461 // Find the corresponding entries
462 // and extract the patchEdgeId
463
464 for (const label connecti : procConnections(otherProci).sortedToc())
465 {
466 const auto& connection = bndEdgeConnections[connecti];
467 const auto& a = connection.first();
468 const auto& b = connection.second();
469
470 if (a.is_localProc())
471 {
472 selectEdges.append(a.patchEdgei());
473 }
474 else if (b.is_localProc())
475 {
476 selectEdges.append(b.patchEdgei());
477 }
478 else
479 {
481 << "Error in programming logic" << nl
482 << abort(FatalError);
483 }
484
485 patchDefLookup[connecti] = -2; // Mark as handled
486 }
487
488 // The edge order is guaranteed to be consistent from the original
489 // getBoundaryEdgeConnections() - sorted by proc/patch/edge
490
491 patchDef.edgeLabels_ = selectEdges;
492 }
493
494
495 // Now convert list of definitions to list of patches
496
497 label nPatches = 0;
498 faPatchList newPatches(faPatchDefs.size());
499
500 for (faPatchData& patchDef : faPatchDefs)
501 {
502 newPatches.set
503 (
504 nPatches,
506 (
507 patchDef.name(), // name
508 patchDef.dict(false), // withEdgeLabels == false
509 nPatches, // index
510 boundary()
511 )
512 );
513
514 // Transfer edge labels
515 newPatches[nPatches].resetEdges(std::move(patchDef.edgeLabels_));
516 ++nPatches;
517 }
518
519
520 if (debug > 1)
521 {
522 label nTotal = 0;
523 Pout<< "Created new finiteArea patches:" << nl;
524 for (const faPatch& p : newPatches)
525 {
526 Pout<< " size: " << p.size()
527 << " name:" << p.name()
528 << " type:" << p.type() << nl;
529 nTotal += p.size();
530 }
531
532 Pout<< "addressed: " << nTotal
533 << '/' << patch().nBoundaryEdges() << " edges" << endl;
534 }
535
536 return newPatches;
537}
538
539
540// ************************************************************************* //
vtk::internalMeshWriter writer(topoMesh, topoCells, vtk::formatType::INLINE_ASCII, runTime.path()/"blockTopology")
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
void transfer(PtrList< T > &list)
Transfer into this list and annul the argument list.
Definition: PtrListI.H:269
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
T & first()
Return the first element of the list.
Definition: UListI.H:202
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:640
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:394
bool checkDefinition(const bool report=false) const
Check boundary definition.
const faBoundaryMesh & boundary() const noexcept
Return constant reference to boundary mesh.
Definition: faMeshI.H:38
void addFaPatches(faPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: faMeshPatches.C:38
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition: fileNameI.H:199
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
int myProcNo() const noexcept
Return processor number.
const fileName & output() const noexcept
The current output file name.
bool writeProcIDs()
Write processor ids as CellData. This is no-op in serial.
virtual bool beginCellData(label nFields=0)
Begin CellData output section for specified number of fields.
virtual bool writeGeometry()
Write mesh topology.
A class for handling words, derived from Foam::string.
Definition: word.H:68
static const word null
An empty word.
Definition: word.H:80
volScalarField & p
faceListList boundary
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const label nPatches
#define WarningInFunction
Report a warning using Foam::Warning.
#define InfoInFunction
Report an information message using Foam::Info.
List< label > sortedToc(const UList< bool > &bools)
Return the (sorted) values corresponding to 'true' entries.
Definition: BitOps.C:203
const std::string patch
OpenFOAM patch number as a std::string.
bool match(const UList< wordRe > &patterns, const std::string &text)
Return true if text matches one of the regular expressions.
Definition: stringOps.H:76
List< label > labelList
A List of labels.
Definition: List.H:66
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:342
errorManip< error > abort(error &err)
Definition: errorManip.H:144
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
error FatalError
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
PtrList< faPatch > faPatchList
Store lists of faPatch as a PtrList.
Definition: faPatch.H:66
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & b
Definition: createFields.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333