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-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "emptyFaPatch.H"
32 #include "ignoreFaPatch.H"
33 #include "processorFaPatch.H"
34 #include "processorPolyPatch.H"
35 #include "foamVtkLineWriter.H"
36 
37 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 
42 // Write edges in VTK format
43 template<class PatchType>
44 static void vtkWritePatchEdges
45 (
46  const PatchType& p,
47  const labelList& selectEdges,
48  const fileName& outputPath,
49  const word& outputName
50 )
51 {
52  edgeList dumpEdges(p.edges(), selectEdges);
53 
55  (
56  p.localPoints(),
57  dumpEdges,
58  outputPath/outputName
59  );
60 
62 
63  // CellData
66  writer.close();
67 }
68 
69 } // End namespace Foam
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
75 (
76  faPatchList& plist,
77  const bool validBoundary
78 )
79 {
80  if (!boundary().empty())
81  {
83  << "boundary already exists"
84  << abort(FatalError);
85  }
86 
87  globalMeshDataPtr_.reset(nullptr);
88 
89  boundary_.transfer(plist);
90 
91  setPrimitiveMeshData();
92 
93  if (validBoundary)
94  {
95  boundary_.checkDefinition();
96  }
97 }
98 
99 
101 (
102  const List<faPatch*>& p,
103  const bool validBoundary
104 )
105 {
106  // Acquire ownership of the pointers
107  faPatchList plist(const_cast<List<faPatch*>&>(p));
108 
109  addFaPatches(plist, validBoundary);
110 }
111 
112 
113 Foam::faPatchList Foam::faMesh::createOnePatch
114 (
115  const word& patchName,
116  const word& patchType
117 ) const
118 {
119  dictionary onePatchDict;
120  if (!patchName.empty())
121  {
122  onePatchDict.add("name", patchName);
123  }
124  if (!patchType.empty())
125  {
126  onePatchDict.add("type", patchType);
127  }
128 
129  return createPatchList
130  (
132  word::null, // Name for empty patch placeholder
133  &onePatchDict // Definitions for defaultPatch
134  );
135 }
136 
137 
138 Foam::faPatchList Foam::faMesh::createPatchList
139 (
140  const dictionary& bndDict,
141  const word& emptyPatchName,
142  const dictionary* defaultPatchDefinition
143 ) const
144 {
145  const polyBoundaryMesh& pbm = mesh().boundaryMesh();
146 
147  // Transcribe into patch definitions
148  DynamicList<faPatchData> faPatchDefs(bndDict.size() + 8);
149  for (const entry& dEntry : bndDict)
150  {
151  if (!dEntry.isDict())
152  {
154  << "Not a dictionary entry: " << dEntry.name() << nl;
155  continue;
156  }
157  const dictionary& patchDict = dEntry.dict();
158 
159  // Add entry
160  auto& patchDef = faPatchDefs.emplace_back();
161  patchDef.name_ = dEntry.keyword();
162  patchDef.type_ = patchDict.get<word>("type");
163 
164  word patchName;
165 
166  // Optional: ownerPolyPatch
167  if (patchDict.readIfPresent("ownerPolyPatch", patchName))
168  {
169  patchDef.ownerPolyPatchId_ = pbm.findPatchID(patchName);
170  if (patchDef.ownerPolyPatchId_ < 0)
171  {
173  << "ownerPolyPatch " << patchName << " not found"
174  << exit(FatalError);
175  }
176  }
177 
178  // Mandatory: neighbourPolyPatch
179  patchDict.readEntry("neighbourPolyPatch", patchName);
180  {
181  patchDef.neighPolyPatchId_ = pbm.findPatchID(patchName);
182  if (patchDef.neighPolyPatchId_ < 0)
183  {
185  << "neighbourPolyPatch " << patchName << " not found"
186  << exit(FatalError);
187  }
188  }
189  }
190 
191  // Additional empty placeholder patch?
192  if (!emptyPatchName.empty())
193  {
194  auto& patchDef = faPatchDefs.emplace_back();
195  patchDef.name_ = emptyPatchName;
196  patchDef.type_ = emptyFaPatch::typeName_();
197  }
198 
199  label nWarnUndefinedPatch(5);
200 
201  // Placeholder for any undefined edges
202  const label undefPatchIndex = faPatchDefs.size();
203  {
204  auto& patchDef = faPatchDefs.emplace_back();
205  patchDef.name_ = "undefined";
206  patchDef.type_ = "patch";
207 
208  if (defaultPatchDefinition)
209  {
210  if
211  (
212  (*defaultPatchDefinition).readIfPresent("name", patchDef.name_)
213  )
214  {
215  // Suppress warnings if defaultPatch name was specified
216  // - probably means we want to use this mechanism
217  nWarnUndefinedPatch = 0;
218  }
219  (*defaultPatchDefinition).readIfPresent("type", patchDef.type_);
220  }
221  }
222 
223  // Placeholder for any undefined edges
224  const label ignorePatchIndex = faPatchDefs.size();
225  {
226  auto& patchDef = faPatchDefs.emplace_back();
227  patchDef.name_ = "_ignore_edges_";
228  patchDef.type_ = ignoreFaPatch::typeName_();
229  }
230 
231  // ----------------------------------------------------------------------
232 
233  // Get edge connections (in globally consistent ordering)
234  List<Pair<patchTuple>> bndEdgeConnections
235  (
236  this->getBoundaryEdgeConnections()
237  );
238 
239  // Update accordingly
240  this->setBoundaryConnections(bndEdgeConnections);
241 
242 
243  // Lookup of patchDef for each connection. Initially all -1 (unvisited)
244  labelList patchDefLookup(bndEdgeConnections.size(), -1);
245 
246  Map<labelHashSet> procConnections;
247  labelHashSet patchDefsUsed;
248 
249  labelHashSet badEdges(2*bndEdgeConnections.size());
250 
251  forAll(bndEdgeConnections, connecti)
252  {
253  const Pair<patchTuple>& connection = bndEdgeConnections[connecti];
254  const auto& a = connection.first();
255  const auto& b = connection.second();
256 
257  edge patchPair;
258 
259  if (!a.valid() || !b.valid())
260  {
261  // If either is invalid, mark as an 'ignore' edge
262  patchDefLookup[connecti] = ignorePatchIndex;
263  patchDefsUsed.insert(ignorePatchIndex);
264  continue;
265  }
266  else if (a.is_finiteArea())
267  {
268  if (b.is_finiteArea())
269  {
270  // Expecting an inter-processor connection
271 
272  if (a.procNo() == b.procNo())
273  {
274  // An intra-processor connection (should not be possible)
276  << "Processor-processor addressing error:" << nl
277  << "Both connections have the same processor: "
278  << a.procNo() << nl
279  << "Connecting patches "
280  << a.realPatchi() << " and " << b.realPatchi() << nl
281  << abort(FatalError);
282  }
283  else if (a.is_localProc())
284  {
285  procConnections(b.procNo()).insert(connecti);
286  }
287  else
288  {
289  procConnections(a.procNo()).insert(connecti);
290  }
291 
292  continue;
293  }
294  else if (a.is_localProc())
295  {
296  patchPair.first() = a.realPatchi();
297  patchPair.second() = b.realPatchi();
298  }
299  }
300  else if (b.is_finiteArea() && b.is_localProc())
301  {
302  patchPair.first() = b.realPatchi();
303  patchPair.second() = a.realPatchi();
304  }
305 
306 
307  // Find first definition with a matching neighbour and
308  // possibly with a matching owner.
309 
310  label bestPatchDefi = -1;
311 
312  const label nPatchDefs = (patchPair.valid() ? faPatchDefs.size() : 0);
313 
314  for (label patchDefi = 0; patchDefi < nPatchDefs; ++patchDefi)
315  {
316  const int match = faPatchDefs[patchDefi].matchPatchPair(patchPair);
317  if (match == 3)
318  {
319  // Exact match (owner/neighbour) - done!
320  bestPatchDefi = patchDefi;
321  break;
322  }
323  else if (match == 2 && bestPatchDefi < 0)
324  {
325  // Match (neighbour) - keep looking for exact match
326  bestPatchDefi = patchDefi;
327  }
328  }
329 
330  if (bestPatchDefi < 0)
331  {
332  bestPatchDefi = undefPatchIndex; // Missed?
333  }
334 
335  patchDefLookup[connecti] = bestPatchDefi;
336  patchDefsUsed.insert(bestPatchDefi);
337  }
338 
339 
340  bool reportBadEdges = false;
341 
342  // Skip undefPatchIndex if not actually needed anywhere
343  if (!returnReduceOr(patchDefsUsed.found(undefPatchIndex)))
344  {
345  faPatchDefs[undefPatchIndex].clear();
346  }
347  else
348  {
349  patchDefsUsed.insert(undefPatchIndex); // Parallel consistency
350  reportBadEdges = true;
351  }
352 
353  // Skip ignorePatchIndex if not actually needed anywhere
354  if (!returnReduceOr(patchDefsUsed.found(ignorePatchIndex)))
355  {
356  faPatchDefs[ignorePatchIndex].clear();
357  }
358  else
359  {
360  patchDefsUsed.insert(ignorePatchIndex); // Parallel consistency
361  reportBadEdges = true;
362  }
363 
364  // Report locations of undefined edges
365  if (reportBadEdges)
366  {
367  badEdges.clear();
368  forAll(patchDefLookup, connecti)
369  {
370  if (patchDefLookup[connecti] == undefPatchIndex)
371  {
372  const auto& connection = bndEdgeConnections[connecti];
373 
374  const auto& a = connection.first();
375  const auto& b = connection.second();
376 
377  if (a.is_localProc() && a.is_finiteArea())
378  {
379  badEdges.insert(a.patchEdgei());
380  }
381  else if (b.is_localProc() && b.is_finiteArea())
382  {
383  badEdges.insert(b.patchEdgei());
384  }
385 
386  if (badEdges.size() <= nWarnUndefinedPatch)
387  {
388  Pout<< "Undefined connection: "
389  << "(patch:" << a.realPatchi()
390  << " face:" << a.meshFacei()
391  << " proc:" << a.procNo()
392  << ") and (patch:" << b.realPatchi()
393  << " face:" << b.meshFacei()
394  << " proc:" << b.procNo()
395  << ") patch:"
396  <<
397  (
398  a.realPatchi() >= 0
399  ? pbm[a.realPatchi()].name()
400  : word::null
401  )
402  << " and patch:"
403  <<
404  (
405  b.realPatchi() >= 0
406  ? pbm[b.realPatchi()].name()
407  : word::null
408  )
409  << nl;
410  }
411  }
412  }
413 
415  {
416  // Report directly as Info, not InfoInFunction
417  // since it can also be an expected result when
418  // nWarnUndefinedPatch == 0
419  Info<< nl
420  << "Had "
421  << returnReduce(badEdges.size(), sumOp<label>()) << '/'
422  << returnReduce(patch().nBoundaryEdges(), sumOp<label>())
423  << " undefined edge connections, added to defaultPatch: "
424  << faPatchDefs[undefPatchIndex].name_ << nl << nl
425  << "==> Could indicate a non-manifold patch geometry" << nl
426  << nl;
427 
428  if (nWarnUndefinedPatch)
429  {
430  labelList selectEdges(badEdges.sortedToc());
431  word outputName("faMesh-construct.undefEdges");
432 
434  (
435  patch(),
436  selectEdges,
437  thisDb().time().globalPath(),
438  outputName
439  );
440 
442  << "(debug) wrote " << outputName << nl;
443  }
444  }
445  }
446 
447  // Report locations of undefined edges
448  if (reportBadEdges)
449  {
450  badEdges.clear();
451  forAll(patchDefLookup, connecti)
452  {
453  if (patchDefLookup[connecti] == ignorePatchIndex)
454  {
455  const auto& connection = bndEdgeConnections[connecti];
456 
457  const auto& a = connection.first();
458  const auto& b = connection.second();
459 
460  if (a.is_localProc() && a.is_finiteArea())
461  {
462  badEdges.insert(a.patchEdgei());
463  }
464  else if (b.is_localProc() && b.is_finiteArea())
465  {
466  badEdges.insert(b.patchEdgei());
467  }
468 
469  if (badEdges.size() <= nWarnUndefinedPatch)
470  {
471  Pout<< "Illegal connection: "
472  << "(patch:" << a.realPatchi()
473  << " face:" << a.meshFacei()
474  << " proc:" << a.procNo()
475  << ") and (patch:" << b.realPatchi()
476  << " face:" << b.meshFacei()
477  << " proc:" << b.procNo()
478  << ") patch:"
479  <<
480  (
481  a.realPatchi() >= 0
482  ? pbm[a.realPatchi()].name()
483  : word::null
484  )
485  << " and patch:"
486  <<
487  (
488  b.realPatchi() >= 0
489  ? pbm[b.realPatchi()].name()
490  : word::null
491  )
492  << nl;
493  }
494  }
495  }
496 
498  {
499  // Report directly as Info, not InfoInFunction
500  // since it can also be an expected result when
501  // nWarnUndefinedPatch == 0
502  Info<< nl
503  << "Had "
504  << returnReduce(badEdges.size(), sumOp<label>()) << '/'
505  << returnReduce(patch().nBoundaryEdges(), sumOp<label>())
506  << " illegal edge connections, added to "
507  << faPatchDefs[ignorePatchIndex].name_ << nl << nl
508  << "==> Could indicate a non-manifold patch geometry" << nl
509  << nl;
510 
511  if (nWarnUndefinedPatch)
512  {
513  labelList selectEdges(badEdges.sortedToc());
514  word outputName("faMesh-construct.ignoreEdges");
515 
517  (
518  patch(),
519  selectEdges,
520  thisDb().time().globalPath(),
521  outputName
522  );
523 
525  << "(debug) wrote " << outputName << nl;
526  }
527  }
528  }
529 
530  // Create processor-processor definitions
531  Map<label> procToDefLookup(2*procConnections.size());
532  {
533  faPatchDefs.reserve(faPatchDefs.size() + procConnections.size());
534 
535  for (const label otherProci : procConnections.sortedToc())
536  {
537  const label patchDefi = faPatchDefs.size();
538  procToDefLookup.insert(otherProci, patchDefi);
539 
540  // Add entry
541  auto& patchDef = faPatchDefs.emplace_back();
542 
543  patchDef.assign_coupled(UPstream::myProcNo(), otherProci);
544  }
545  }
546 
547 
548  // Extract which edges map to which patch
549  // and set edgeLabels for each faPatch
550 
551  DynamicList<label> selectEdges(bndEdgeConnections.size());
552  label nOffProcessorEdges = 0;
553 
554  for (const label patchDefi : patchDefsUsed.sortedToc())
555  {
556  auto& patchDef = faPatchDefs[patchDefi];
557  selectEdges.clear();
558 
559  // Find the corresponding entries
560  // and extract the patchEdgeId
561 
562  forAll(patchDefLookup, connecti)
563  {
564  if (patchDefLookup[connecti] == patchDefi)
565  {
566  const auto& a = bndEdgeConnections[connecti].first();
567  const auto& b = bndEdgeConnections[connecti].second();
568 
569  if (a.is_localProc() && a.is_finiteArea())
570  {
571  selectEdges.push_back(a.patchEdgei());
572  }
573  else if (b.is_localProc() && b.is_finiteArea())
574  {
575  selectEdges.push_back(b.patchEdgei());
576  }
577  else if (a.valid() && b.valid())
578  {
580  << "Error in programming logic" << nl
581  << abort(FatalError);
582  }
583 
584  if (a.is_localProc() != b.is_localProc())
585  {
586  ++nOffProcessorEdges;
587  }
588 
589  patchDefLookup[connecti] = -2; // Mark as handled
590  }
591  }
592 
593  // Remove any cosmetic sorting artifacts from off-processor
594  // termination by doing using a regular sort here.
595 
596  Foam::sort(selectEdges);
597  patchDef.edgeLabels_ = selectEdges;
598  }
599 
600  if (debug)
601  {
602  Pout<< "Had " << nOffProcessorEdges
603  << " patch edges connected off-processor" << endl;
604 
606  << "Total "
607  << returnReduce(nOffProcessorEdges, sumOp<label>())
608  << " patch edges connected off-processor" << endl;
609  }
610 
611 
612  // Processor patches
613  for (const label otherProci : procToDefLookup.sortedToc())
614  {
615  const label patchDefi = procToDefLookup[otherProci];
616 
617  auto& patchDef = faPatchDefs[patchDefi];
618  selectEdges.clear();
619 
620  // Find the corresponding entries
621  // and extract the patchEdgeId
622 
623  for (const label connecti : procConnections(otherProci).sortedToc())
624  {
625  const auto& connection = bndEdgeConnections[connecti];
626  const auto& a = connection.first();
627  const auto& b = connection.second();
628 
629  if (a.is_localProc())
630  {
631  selectEdges.push_back(a.patchEdgei());
632  }
633  else if (b.is_localProc())
634  {
635  selectEdges.push_back(b.patchEdgei());
636  }
637  else
638  {
640  << "Error in programming logic" << nl
641  << abort(FatalError);
642  }
643 
644  patchDefLookup[connecti] = -2; // Mark as handled
645  }
646 
647  // The edge order is guaranteed to be consistent from the original
648  // getBoundaryEdgeConnections() - sorted by proc/patch/edge
649 
650  patchDef.edgeLabels_ = selectEdges;
651  }
652 
653 
654  // Now convert list of definitions to list of patches
655 
656  label nPatches = 0;
657  faPatchList newPatches(faPatchDefs.size());
658 
659  for (faPatchData& patchDef : faPatchDefs)
660  {
661  if (!patchDef.good())
662  {
663  continue;
664  }
665 
666  newPatches.set
667  (
668  nPatches,
670  (
671  patchDef.name(), // name
672  patchDef.dict(false), // withEdgeLabels == false
673  nPatches, // index
674  boundary()
675  )
676  );
677 
678  // Transfer edge labels
679  newPatches[nPatches].resetEdges(std::move(patchDef.edgeLabels_));
680  ++nPatches;
681  }
682  newPatches.resize(nPatches);
683 
684  if (debug > 1)
685  {
686  label nTotal = 0;
687  Pout<< "Created new finiteArea patches:" << nl;
688  for (const faPatch& p : newPatches)
689  {
690  Pout<< " size: " << p.size()
691  << " name:" << p.name()
692  << " type:" << p.type() << nl;
693  nTotal += p.size();
694  }
695 
696  Pout<< "addressed: " << nTotal
697  << '/' << patch().nBoundaryEdges() << " edges" << endl;
698  }
699 
700  return newPatches;
701 }
702 
703 
704 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:396
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
PtrList< faPatch > faPatchList
Store lists of faPatch as a PtrList.
Definition: faPatch.H:59
faceListList boundary
bool match(const UList< wordRe > &patterns, const std::string &text)
Return true if text matches one of the regular expressions.
Definition: stringOps.H:77
const polyBoundaryMesh & pbm
static autoPtr< faPatch > New(const word &name, const dictionary &dict, const label index, const faBoundaryMesh &bm)
Return pointer to a new patch created on freestore from dictionary.
Definition: faPatchNew.C:28
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
bool checkDefinition(const bool report=false) const
Check boundary definition.
void addFaPatches(faPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: faMeshPatches.C:68
A class for handling file names.
Definition: fileName.H:72
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/"finiteArea-edgesCentres"))
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
srcOptions insert("case", fileName(rootDirSource/caseDirSource))
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
bool writeProcIDs()
Write processor ids for each line as CellData or for each point as PointData, depending on isPointDat...
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:625
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:232
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
labelHashSet badEdges(pp.nEdges()/20)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
label size() const noexcept
The number of elements in table.
Definition: HashTable.H:342
word outputName("finiteArea-edges.obj")
dynamicFvMesh & mesh
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:296
void clear()
Remove all entries from table.
Definition: HashTable.C:725
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:474
static const word null
An empty word.
Definition: word.H:84
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Write edge/points (optionally with fields) as a vtp file or a legacy vtk file.
static void vtkWritePatchEdges(const PatchType &p, const labelList &selectEdges, const fileName &outputPath, const word &outputName)
Definition: faMeshPatches.C:38
virtual bool beginCellData(label nFields=0)
Begin CellData output section for specified number of fields.
int debug
Static debugging option.
const faBoundaryMesh & boundary() const noexcept
Return constant reference to boundary mesh.
Definition: faMeshI.H:31
#define WarningInFunction
Report a warning using Foam::Warning.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:139
void close()
End the file contents and close the file after writing.
const std::string patch
OpenFOAM patch number as a std::string.
List< label > sortedToc(const UList< bool > &bools)
Return the (sorted) values corresponding to &#39;true&#39; entries.
Definition: BitOps.C:195
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual bool writeGeometry()
Write patch topology.
void transfer(PtrList< T > &list)
Transfer into this list and annul the argument list.
Definition: PtrListI.H:272
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Namespace for OpenFOAM.
#define InfoInFunction
Report an information message using Foam::Info.