cyclicACMIPolyPatch.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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2020 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 "cyclicACMIPolyPatch.H"
30 #include "SubField.H"
31 #include "Time.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(cyclicACMIPolyPatch, 0);
39 
40  addToRunTimeSelectionTable(polyPatch, cyclicACMIPolyPatch, word);
41  addToRunTimeSelectionTable(polyPatch, cyclicACMIPolyPatch, dictionary);
42 }
43 
44 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
45 
47 (
48  const word& name,
49  const scalarField& weightSum
50 ) const
51 {
52  label nUncovered = 0;
53  label nCovered = 0;
54  for (const scalar sum : weightSum)
55  {
56  if (sum < tolerance_)
57  {
58  ++nUncovered;
59  }
60  else if (sum > scalar(1) - tolerance_)
61  {
62  ++nCovered;
63  }
64  }
65  reduce(nUncovered, sumOp<label>());
66  reduce(nCovered, sumOp<label>());
67  label nTotal = returnReduce(weightSum.size(), sumOp<label>());
68 
69  Info<< "ACMI: Patch " << name << " uncovered/blended/covered = "
70  << nUncovered << ", " << nTotal-nUncovered-nCovered
71  << ", " << nCovered << endl;
72 }
73 
74 
76 {
78 }
79 
80 
82 {
83  if (!owner())
84  {
85  return;
86  }
87 
88  const polyPatch& nonOverlapPatch = this->nonOverlapPatch();
89 
90  DebugPout
91  << "cyclicACMIPolyPatch::resetAMI : recalculating weights"
92  << " for " << name() << " and " << nonOverlapPatch.name()
93  << endl;
94 
95  const polyMesh& mesh = boundaryMesh().mesh();
96 
97  if (!createAMIFaces_ && mesh.hasCellCentres())
98  {
99  DebugPout
100  << "cyclicACMIPolyPatch::resetAMI : clearing cellCentres"
101  << " for " << name() << " and " << nonOverlapPatch.name()
102  << endl;
103 
105  << "The mesh already has cellCentres calculated when"
106  << " resetting ACMI " << name() << "." << nl
107  << "This is a problem since ACMI adapts the face areas"
108  << " (to close cells) so this has" << nl
109  << "to be done before cell centre calculation." << nl
110  << "This can happen if e.g. the cyclicACMI is after"
111  << " any processor patches in the boundary." << endl;
112  const_cast<polyMesh&>(mesh).primitiveMesh::clearGeom();
113  }
114 
115 
116  // Trigger re-building of faceAreas
117  (void)mesh.faceAreas();
118 
119 
120  // Calculate the AMI using partial face-area-weighted. This leaves
121  // the weights as fractions of local areas (sum(weights) = 1 means
122  // face is fully covered)
124 
125  const AMIPatchToPatchInterpolation& AMI = this->AMI();
126 
127  // Output some statistics
128  reportCoverage("source", AMI.srcWeightsSum());
129  reportCoverage("target", AMI.tgtWeightsSum());
130 
131  // Set the mask fields
132  // Note:
133  // - assumes that the non-overlap patches are decomposed using the same
134  // decomposition as the coupled patches (per side)
135  srcMask_ = min(scalar(1), max(scalar(0), AMI.srcWeightsSum()));
136  tgtMask_ = min(scalar(1), max(scalar(0), AMI.tgtWeightsSum()));
137 
138  if (debug)
139  {
140  Pout<< "resetAMI" << endl;
141  {
142  const cyclicACMIPolyPatch& patch = *this;
143  Pout<< "patch:" << patch.name() << " size:" << patch.size()
144  << " non-overlap patch: " << patch.nonOverlapPatch().name()
145  << " size:" << patch.nonOverlapPatch().size()
146  << " mask size:" << patch.srcMask().size() << endl;
147  }
148  {
149  const cyclicACMIPolyPatch& patch = this->neighbPatch();
150  Pout<< "patch:" << patch.name() << " size:" << patch.size()
151  << " non-overlap patch: " << patch.nonOverlapPatch().name()
152  << " size:" << patch.nonOverlapPatch().size()
153  << " mask size:" << patch.neighbPatch().tgtMask().size()
154  << endl;
155  }
156  }
157 }
158 
159 
161 {
162  if (!owner() || !canResetAMI())
163  {
164  return;
165  }
166 
167  scalePatchFaceAreas(*this);
168  scalePatchFaceAreas(this->neighbPatch());
169 }
170 
171 
173 (
174  const cyclicACMIPolyPatch& acmipp
175 )
176 {
177  // Primitive patch face areas have been cleared/reset based on the raw
178  // points - need to reset to avoid double-accounting of face areas
179 
180  const scalar maxTol = scalar(1) - tolerance_;
181  const scalarField& mask = acmipp.mask();
182 
183  const polyPatch& nonOverlapPatch = acmipp.nonOverlapPatch();
184  vectorField::subField noSf = nonOverlapPatch.faceAreas();
185 
186  DebugPout
187  << "rescaling non-overlap patch areas for: " << nonOverlapPatch.name()
188  << endl;
189 
190 
191  if (mask.size() != noSf.size())
192  {
194  << "Inconsistent sizes for patch: " << acmipp.name()
195  << " - not manipulating patches" << nl
196  << " - size: " << size() << nl
197  << " - non-overlap patch size: " << noSf.size() << nl
198  << " - mask size: " << mask.size() << nl
199  << "This is OK for decomposition but should be considered fatal "
200  << "at run-time" << endl;
201 
202  return;
203  }
204 
205  forAll(noSf, facei)
206  {
207  const scalar w = min(maxTol, max(tolerance_, mask[facei]));
208  noSf[facei] *= scalar(1) - w;
209  }
210 
211  if (!createAMIFaces_)
212  {
213  // Note: for topological update (createAMIFaces_ = true)
214  // AMI coupled patch face areas are updated as part of the topological
215  // updates, e.g. by the calls to cyclicAMIPolyPatch's setTopology and
216  // initMovePoints
217  DebugPout
218  << "scaling coupled patch areas for: " << acmipp.name() << endl;
219 
220  // Scale the coupled patch face areas
221  vectorField::subField Sf = acmipp.faceAreas();
222 
223  forAll(Sf, facei)
224  {
225  Sf[facei] *= max(tolerance_, mask[facei]);
226  }
227 
228  // Re-normalise the weights since the effect of overlap is already
229  // accounted for in the area
230  auto& weights = const_cast<scalarListList&>(acmipp.weights());
231  auto& weightsSum = const_cast<scalarField&>(acmipp.weightsSum());
232  forAll(weights, i)
233  {
234  scalarList& wghts = weights[i];
235  if (wghts.size())
236  {
237  scalar& sum = weightsSum[i];
238 
239  forAll(wghts, j)
240  {
241  wghts[j] /= sum;
242  }
243  sum = 1.0;
244  }
245  }
246  }
247 }
248 
249 
251 {
252  DebugPout << "cyclicACMIPolyPatch::initGeometry : " << name() << endl;
253 
254  // Note: calculates transformation and triggers face centre calculation
256 
257  // Initialise the AMI early to make sure we adapt the face areas before the
258  // cell centre calculation gets triggered.
259  if (!createAMIFaces_ && canResetAMI())
260  {
261  resetAMI();
262  }
263 
264  scalePatchFaceAreas();
265 }
266 
267 
269 {
270  DebugPout << "cyclicACMIPolyPatch::calcGeometry : " << name() << endl;
271 
273 }
274 
275 
277 (
278  PstreamBuffers& pBufs,
279  const pointField& p
280 )
281 {
282  DebugPout<< "cyclicACMIPolyPatch::initMovePoints : " << name() << endl;
283 
284  // Note: calculates transformation and triggers face centre calculation
285  // - Note: resetAMI called by cyclicAMIPolyPatch::initMovePoints
287 
288  scalePatchFaceAreas();
289 }
290 
291 
293 (
294  PstreamBuffers& pBufs,
295  const pointField& p
296 )
297 {
298  DebugPout << "cyclicACMIPolyPatch::movePoints : " << name() << endl;
299 
300  // When topology is changing, this will scale the duplicate AMI faces
302 }
303 
304 
306 {
307  DebugPout << "cyclicACMIPolyPatch::initUpdateMesh : " << name() << endl;
308 
310 }
311 
312 
314 {
315  DebugPout << "cyclicACMIPolyPatch::updateMesh : " << name() << endl;
316 
318 }
319 
320 
322 {
323  DebugPout << "cyclicACMIPolyPatch::clearGeom : " << name() << endl;
324 
326 }
327 
328 
330 {
331  return srcMask_;
332 }
333 
334 
336 {
337  return tgtMask_;
338 }
339 
340 
341 // * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
342 
344 (
345  const word& name,
346  const label size,
347  const label start,
348  const label index,
349  const polyBoundaryMesh& bm,
350  const word& patchType,
351  const transformType transform,
352  const word& defaultAMIMethod
353 )
354 :
356  (
357  name,
358  size,
359  start,
360  index,
361  bm,
362  patchType,
363  transform,
364  defaultAMIMethod
365  ),
366  nonOverlapPatchName_(word::null),
367  nonOverlapPatchID_(-1),
368  srcMask_(),
369  tgtMask_()
370 {
371  AMIPtr_->setRequireMatch(false);
372 
373  // Non-overlapping patch might not be valid yet so cannot determine
374  // associated patchID
375 }
376 
377 
379 (
380  const word& name,
381  const dictionary& dict,
382  const label index,
383  const polyBoundaryMesh& bm,
384  const word& patchType,
385  const word& defaultAMIMethod
386 )
387 :
388  cyclicAMIPolyPatch(name, dict, index, bm, patchType, defaultAMIMethod),
389  nonOverlapPatchName_(dict.get<word>("nonOverlapPatch")),
390  nonOverlapPatchID_(-1),
391  srcMask_(),
392  tgtMask_()
393 {
394  AMIPtr_->setRequireMatch(false);
395 
396  if (nonOverlapPatchName_ == name)
397  {
399  << "Non-overlapping patch name " << nonOverlapPatchName_
400  << " cannot be the same as this patch " << name
401  << exit(FatalIOError);
402  }
403 
404  // Non-overlapping patch might not be valid yet so cannot determine
405  // associated patchID
406 }
407 
408 
410 (
411  const cyclicACMIPolyPatch& pp,
412  const polyBoundaryMesh& bm
413 )
414 :
415  cyclicAMIPolyPatch(pp, bm),
416  nonOverlapPatchName_(pp.nonOverlapPatchName_),
417  nonOverlapPatchID_(-1),
418  srcMask_(),
419  tgtMask_()
420 {
421  AMIPtr_->setRequireMatch(false);
422 
423  // Non-overlapping patch might not be valid yet so cannot determine
424  // associated patchID
425 }
426 
427 
429 (
430  const cyclicACMIPolyPatch& pp,
431  const polyBoundaryMesh& bm,
432  const label index,
433  const label newSize,
434  const label newStart,
435  const word& nbrPatchName,
436  const word& nonOverlapPatchName
437 )
438 :
439  cyclicAMIPolyPatch(pp, bm, index, newSize, newStart, nbrPatchName),
440  nonOverlapPatchName_(nonOverlapPatchName),
441  nonOverlapPatchID_(-1),
442  srcMask_(),
443  tgtMask_()
444 {
445  AMIPtr_->setRequireMatch(false);
446 
447  if (nonOverlapPatchName_ == name())
448  {
450  << "Non-overlapping patch name " << nonOverlapPatchName_
451  << " cannot be the same as this patch " << name()
452  << exit(FatalError);
453  }
454 
455  // Non-overlapping patch might not be valid yet so cannot determine
456  // associated patchID
457 }
458 
459 
461 (
462  const cyclicACMIPolyPatch& pp,
463  const polyBoundaryMesh& bm,
464  const label index,
465  const labelUList& mapAddressing,
466  const label newStart
467 )
468 :
469  cyclicAMIPolyPatch(pp, bm, index, mapAddressing, newStart),
470  nonOverlapPatchName_(pp.nonOverlapPatchName_),
471  nonOverlapPatchID_(-1),
472  srcMask_(),
473  tgtMask_()
474 {
475  AMIPtr_->setRequireMatch(false);
476 }
477 
478 
479 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
480 
482 {
483  const polyPatch& pp = this->boundaryMesh()[neighbPatchID()];
484  return refCast<const cyclicACMIPolyPatch>(pp);
485 }
486 
487 
489 {
490  if (nonOverlapPatchID_ == -1)
491  {
492  nonOverlapPatchID_ =
493  this->boundaryMesh().findPatchID(nonOverlapPatchName_);
494 
495  if (nonOverlapPatchID_ == -1)
496  {
498  << "Illegal non-overlapping patch name " << nonOverlapPatchName_
499  << nl << "Valid patch names are "
500  << this->boundaryMesh().names()
501  << exit(FatalError);
502  }
503 
504  if (nonOverlapPatchID_ < index())
505  {
507  << "Boundary ordering error: " << type()
508  << " patch must be defined prior to its non-overlapping patch"
509  << nl
510  << type() << " patch: " << name() << ", ID:" << index() << nl
511  << "Non-overlap patch: " << nonOverlapPatchName_
512  << ", ID:" << nonOverlapPatchID_ << nl
513  << exit(FatalError);
514  }
515 
516  const polyPatch& noPp = this->boundaryMesh()[nonOverlapPatchID_];
517 
518  bool ok = true;
519 
520  if (size() == noPp.size())
521  {
522  const scalarField magSf(mag(faceAreas()));
523  const scalarField noMagSf(mag(noPp.faceAreas()));
524 
525  forAll(magSf, facei)
526  {
527  scalar ratio = mag(magSf[facei]/(noMagSf[facei] + ROOTVSMALL));
528 
529  if (ratio - 1 > tolerance_)
530  {
531  ok = false;
532  break;
533  }
534  }
535  }
536  else
537  {
538  ok = false;
539  }
540 
541  if (!ok)
542  {
544  << "Inconsistent ACMI patches " << name() << " and "
545  << noPp.name() << ". Patches should have identical topology"
546  << exit(FatalError);
547  }
548  }
549 
550  return nonOverlapPatchID_;
551 }
552 
553 
555 (
556  PstreamBuffers& pBufs,
557  const primitivePatch& pp
558 ) const
559 {
561 }
562 
563 
565 (
566  PstreamBuffers& pBufs,
567  const primitivePatch& pp,
569  labelList& rotation
570 ) const
571 {
572  return cyclicAMIPolyPatch::order(pBufs, pp, faceMap, rotation);
573 }
574 
575 
577 {
579 
580  os.writeEntry("nonOverlapPatch", nonOverlapPatchName_);
581 }
582 
583 
584 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::cyclicAMIPolyPatch::order
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
Definition: cyclicAMIPolyPatch.C:1040
Foam::boundaryMesh::mesh
const bMesh & mesh() const
Definition: boundaryMesh.H:206
Foam::primitiveMesh::clearGeom
void clearGeom()
Clear geometry.
Definition: primitiveMeshClear.C:127
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::cyclicACMIPolyPatch::calcGeometry
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
Definition: cyclicACMIPolyPatch.C:268
Foam::cyclicACMIPolyPatch::movePoints
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
Definition: cyclicACMIPolyPatch.C:293
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeTopological.C:94
Foam::cyclicAMIPolyPatch::resetAMI
virtual void resetAMI() const
Reset the AMI interpolator, use current patch points.
Definition: cyclicAMIPolyPatch.C:323
cyclicACMIPolyPatch.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
SubField.H
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::cyclicAMIPolyPatch::initOrder
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Definition: cyclicAMIPolyPatch.C:1032
Foam::AMIInterpolation::tgtWeightsSum
const scalarField & tgtWeightsSum() const
Definition: AMIInterpolationI.H:213
Foam::cyclicAMIPolyPatch::initGeometry
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: cyclicAMIPolyPatch.C:436
Foam::PstreamBuffers
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Definition: PstreamBuffers.H:88
Foam::cyclicAMIPolyPatch::weights
const scalarListList & weights() const
Helper function to return the weights.
Definition: cyclicAMIPolyPatchI.H:61
Foam::AMIInterpolation::srcWeightsSum
const scalarField & srcWeightsSum() const
Definition: AMIInterpolationI.H:147
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::cyclicACMIPolyPatch::resetAMI
virtual void resetAMI() const
Reset the AMI interpolator, use current patch points.
Definition: cyclicACMIPolyPatch.C:75
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:519
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::cyclicAMIPolyPatch::initMovePoints
virtual void initMovePoints(PstreamBuffers &pBufs, const pointField &)
Initialise the patches for moving points.
Definition: cyclicAMIPolyPatch.C:460
Foam::cyclicAMIPolyPatch::updateMesh
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: cyclicAMIPolyPatch.C:532
Foam::SubField
SubField is a Field obtained as a section of another Field.
Definition: Field.H:64
Foam::cyclicACMIPolyPatch::updateMesh
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: cyclicACMIPolyPatch.C:313
Foam::cyclicACMIPolyPatch::write
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: cyclicACMIPolyPatch.C:576
Foam::cyclicAMIPolyPatch::weightsSum
const scalarField & weightsSum() const
Helper function to return the weights sum.
Definition: cyclicAMIPolyPatchI.H:72
Foam::cyclicAMIPolyPatch::calcGeometry
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
Definition: cyclicAMIPolyPatch.C:453
Foam::Field< scalar >
Foam::cyclicAMIPolyPatch::write
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: cyclicAMIPolyPatch.C:1104
Foam::polyPatch::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:307
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:67
Foam::cyclicACMIPolyPatch::srcMask
virtual const scalarField & srcMask() const
Return the mask/weighting for the source patch.
Definition: cyclicACMIPolyPatch.C:329
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::primitiveMesh::hasCellCentres
bool hasCellCentres() const
Definition: primitiveMeshI.H:186
Foam::cyclicAMIPolyPatch::movePoints
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
Definition: cyclicAMIPolyPatch.C:496
Foam::cyclicACMIPolyPatch::clearGeom
virtual void clearGeom()
Clear geometry.
Definition: cyclicACMIPolyPatch.C:321
Foam::cyclicACMIPolyPatch::order
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
Definition: cyclicACMIPolyPatch.C:565
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::cyclicACMIPolyPatch::initOrder
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
Definition: cyclicACMIPolyPatch.C:555
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::cyclicACMIPolyPatch::nonOverlapPatchID
virtual label nonOverlapPatchID() const
Non-overlapping patch ID.
Definition: cyclicACMIPolyPatch.C:488
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
reduce
reduce(hasMovingMesh, orOp< bool >())
Foam::PrimitivePatch::points
const Field< point_type > & points() const
Return reference to global points.
Definition: PrimitivePatch.H:305
Foam::cyclicACMIPolyPatch::cyclicACMIPolyPatch
cyclicACMIPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform=UNKNOWN, const word &defaultAMIMethod=partialFaceAreaWeightAMI::typeName)
Construct from (base coupled patch) components.
Definition: cyclicACMIPolyPatch.C:344
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
DebugPout
#define DebugPout
Report an information message using Foam::Pout.
Definition: messageStream.H:365
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::cyclicACMIPolyPatch::initMovePoints
virtual void initMovePoints(PstreamBuffers &pBufs, const pointField &)
Initialise the patches for moving points.
Definition: cyclicACMIPolyPatch.C:277
Foam::cyclicACMIPolyPatch::scalePatchFaceAreas
virtual void scalePatchFaceAreas()
Scale patch face areas to maintain physical area.
Definition: cyclicACMIPolyPatch.C:160
Time.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::cyclicACMIPolyPatch::neighbPatch
virtual const cyclicACMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
Definition: cyclicACMIPolyPatch.C:481
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::cyclicAMIPolyPatch::initUpdateMesh
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: cyclicAMIPolyPatch.C:519
Foam::cyclicAMIPolyPatch::clearGeom
virtual void clearGeom()
Clear geometry.
Definition: cyclicAMIPolyPatch.C:541
Foam::List< scalarList >
Foam::polyPatch::faceAreas
const vectorField::subField faceAreas() const
Return face normals.
Definition: polyPatch.C:319
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::AMIInterpolation
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
Definition: AMIInterpolation.H:79
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:103
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
Foam::cyclicACMIPolyPatch::nonOverlapPatch
const polyPatch & nonOverlapPatch() const
Return a const reference to the non-overlapping patch.
Definition: cyclicACMIPolyPatchI.H:36
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:62
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:392
Foam::cyclicACMIPolyPatch::tgtMask
virtual const scalarField & tgtMask() const
Return the mask/weighting for the target patch.
Definition: cyclicACMIPolyPatch.C:335
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::cyclicACMIPolyPatch
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI)
Definition: cyclicACMIPolyPatch.H:54
Foam::coupledPolyPatch::transformType
transformType
Definition: coupledPolyPatch.H:59
Foam::cyclicACMIPolyPatch::mask
const scalarField & mask() const
Mask field where 1 = overlap(coupled), 0 = no-overlap.
Definition: cyclicACMIPolyPatchI.H:54
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::patchIdentifier::name
const word & name() const
The patch name.
Definition: patchIdentifier.H:134
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:298
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:85
Foam::cyclicACMIPolyPatch::initGeometry
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: cyclicACMIPolyPatch.C:250
Foam::cyclicACMIPolyPatch::reportCoverage
void reportCoverage(const word &name, const scalarField &weightSum) const
Definition: cyclicACMIPolyPatch.C:47
Foam::cyclicACMIPolyPatch::initUpdateMesh
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: cyclicACMIPolyPatch.C:305
Foam::primitiveMesh::faceAreas
const vectorField & faceAreas() const
Definition: primitiveMeshFaceCentresAndAreas.C:155
Foam::cyclicAMIPolyPatch
Cyclic patch for Arbitrary Mesh Interface (AMI)
Definition: cyclicAMIPolyPatch.H:67