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-2018 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 const Foam::scalar Foam::cyclicACMIPolyPatch::tolerance_ = 1e-10;
45 
46 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47 
49 (
51 ) const
52 {
53  if (owner())
54  {
55  const polyPatch& nonOverlapPatch = this->nonOverlapPatch();
56 
57  if (debug)
58  {
59  Pout<< "cyclicACMIPolyPatch::resetAMI : recalculating weights"
60  << " for " << name() << " and " << nonOverlapPatch.name()
61  << endl;
62  }
63 
64  if (boundaryMesh().mesh().hasCellCentres())
65  {
66  if (debug)
67  {
68  Pout<< "cyclicACMIPolyPatch::resetAMI : clearing cellCentres"
69  << " for " << name() << " and " << nonOverlapPatch.name()
70  << endl;
71  }
72 
73  //WarningInFunction
74  // << "The mesh already has cellCentres calculated when"
75  // << " resetting ACMI " << name() << "." << endl
76  // << "This is a problem since ACMI adapts the face areas"
77  // << " (to close cells) so this has" << endl
78  // << "to be done before cell centre calculation." << endl
79  // << "This can happen if e.g. the cyclicACMI is after"
80  // << " any processor patches in the boundary." << endl;
81  const_cast<polyMesh&>
82  (
83  boundaryMesh().mesh()
84  ).primitiveMesh::clearGeom();
85  }
86 
87 
88  // Trigger re-building of faceAreas
89  (void)boundaryMesh().mesh().faceAreas();
90 
91 
92  // Calculate the AMI using partial face-area-weighted. This leaves
93  // the weights as fractions of local areas (sum(weights) = 1 means
94  // face is fully covered)
95  cyclicAMIPolyPatch::resetAMI
96  (
97  AMIPatchToPatchInterpolation::imPartialFaceAreaWeight
98  );
99 
101  const_cast<AMIPatchToPatchInterpolation&>(this->AMI());
102 
103  // Output some stats. AMIInterpolation will have already output the
104  // average weights ("sum(weights) min:1 max:1 average:1")
105  {
106  const scalarField& wghtsSum = AMI.srcWeightsSum();
107 
108  label nUncovered = 0;
109  label nCovered = 0;
110  forAll(wghtsSum, facei)
111  {
112  scalar sum = wghtsSum[facei];
113  if (sum < tolerance_)
114  {
115  nUncovered++;
116  }
117  else if (sum > scalar(1)-tolerance_)
118  {
119  nCovered++;
120  }
121  }
122  reduce(nUncovered, sumOp<label>());
123  reduce(nCovered, sumOp<label>());
124  label nTotal = returnReduce(wghtsSum.size(), sumOp<label>());
125 
126  Info<< "ACMI: Patch source uncovered/blended/covered = "
127  << nUncovered << ", " << nTotal-nUncovered-nCovered
128  << ", " << nCovered << endl;
129  }
130  {
131  const scalarField& wghtsSum = AMI.tgtWeightsSum();
132 
133  label nUncovered = 0;
134  label nCovered = 0;
135  forAll(wghtsSum, facei)
136  {
137  scalar sum = wghtsSum[facei];
138  if (sum < tolerance_)
139  {
140  nUncovered++;
141  }
142  else if (sum > scalar(1)-tolerance_)
143  {
144  nCovered++;
145  }
146  }
147  reduce(nUncovered, sumOp<label>());
148  reduce(nCovered, sumOp<label>());
149  label nTotal = returnReduce(wghtsSum.size(), sumOp<label>());
150 
151  Info<< "ACMI: Patch target uncovered/blended/covered = "
152  << nUncovered << ", " << nTotal-nUncovered-nCovered
153  << ", " << nCovered << endl;
154  }
155 
156  srcMask_ =
157  min(scalar(1) - tolerance_, max(tolerance_, AMI.srcWeightsSum()));
158 
159  tgtMask_ =
160  min(scalar(1) - tolerance_, max(tolerance_, AMI.tgtWeightsSum()));
161 
162 
163  // Adapt owner side areas. Note that in uncoupled situations (e.g.
164  // decomposePar) srcMask, tgtMask can be zero size.
165  if (srcMask_.size())
166  {
167  vectorField::subField Sf = faceAreas();
168  vectorField::subField noSf = nonOverlapPatch.faceAreas();
169 
170  forAll(Sf, facei)
171  {
172  Sf[facei] *= srcMask_[facei];
173  noSf[facei] *= 1.0 - srcMask_[facei];
174  }
175  }
176  // Adapt slave side areas
177  if (tgtMask_.size())
178  {
179  const cyclicACMIPolyPatch& cp =
180  refCast<const cyclicACMIPolyPatch>(this->neighbPatch());
181  const polyPatch& pp = cp.nonOverlapPatch();
182 
183  vectorField::subField Sf = cp.faceAreas();
184  vectorField::subField noSf = pp.faceAreas();
185 
186  forAll(Sf, facei)
187  {
188  Sf[facei] *= tgtMask_[facei];
189  noSf[facei] *= 1.0 - tgtMask_[facei];
190  }
191  }
192 
193  // Re-normalise the weights since the effect of overlap is already
194  // accounted for in the area.
195  {
196  scalarListList& srcWeights = AMI.srcWeights();
197  scalarField& srcWeightsSum = AMI.srcWeightsSum();
198  forAll(srcWeights, i)
199  {
200  scalarList& wghts = srcWeights[i];
201  if (wghts.size())
202  {
203  scalar& sum = srcWeightsSum[i];
204 
205  forAll(wghts, j)
206  {
207  wghts[j] /= sum;
208  }
209  sum = 1.0;
210  }
211  }
212  }
213  {
214  scalarListList& tgtWeights = AMI.tgtWeights();
215  scalarField& tgtWeightsSum = AMI.tgtWeightsSum();
216  forAll(tgtWeights, i)
217  {
218  scalarList& wghts = tgtWeights[i];
219  if (wghts.size())
220  {
221  scalar& sum = tgtWeightsSum[i];
222  forAll(wghts, j)
223  {
224  wghts[j] /= sum;
225  }
226  sum = 1.0;
227  }
228  }
229  }
230 
231  // Set the updated flag
232  updated_ = true;
233  }
234 }
235 
236 
238 {
239  if (debug)
240  {
241  Pout<< "cyclicACMIPolyPatch::initGeometry : " << name() << endl;
242  }
243 
244  // Note: calculates transformation and triggers face centre calculation
246 
247  // Initialise the AMI early to make sure we adapt the face areas before the
248  // cell centre calculation gets triggered.
249  resetAMI();
250 }
251 
252 
254 {
255  if (debug)
256  {
257  Pout<< "cyclicACMIPolyPatch::calcGeometry : " << name() << endl;
258  }
260 }
261 
262 
264 (
265  PstreamBuffers& pBufs,
266  const pointField& p
267 )
268 {
269  if (debug)
270  {
271  Pout<< "cyclicACMIPolyPatch::initMovePoints : " << name() << endl;
272  }
273 
274  // Note: calculates transformation and triggers face centre calculation
276 
277  // Initialise the AMI early. See initGeometry.
278  resetAMI();
279 }
280 
281 
283 (
284  PstreamBuffers& pBufs,
285  const pointField& p
286 )
287 {
288  if (debug)
289  {
290  Pout<< "cyclicACMIPolyPatch::movePoints : " << name() << endl;
291  }
293 }
294 
295 
297 {
298  if (debug)
299  {
300  Pout<< "cyclicACMIPolyPatch::initUpdateMesh : " << name() << endl;
301  }
303 }
304 
305 
307 {
308  if (debug)
309  {
310  Pout<< "cyclicACMIPolyPatch::updateMesh : " << name() << endl;
311  }
313 }
314 
315 
317 {
318  if (debug)
319  {
320  Pout<< "cyclicACMIPolyPatch::clearGeom : " << name() << endl;
321  }
323 }
324 
325 
327 {
328  return srcMask_;
329 }
330 
331 
333 {
334  return tgtMask_;
335 }
336 
337 
338 // * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
339 
341 (
342  const word& name,
343  const label size,
344  const label start,
345  const label index,
346  const polyBoundaryMesh& bm,
347  const word& patchType,
349 )
350 :
351  cyclicAMIPolyPatch(name, size, start, index, bm, patchType, transform),
352  nonOverlapPatchName_(word::null),
353  nonOverlapPatchID_(-1),
354  srcMask_(),
355  tgtMask_(),
356  updated_(false)
357 {
358  AMIRequireMatch_ = false;
359 
360  // Non-overlapping patch might not be valid yet so cannot determine
361  // associated patchID
362 }
363 
364 
366 (
367  const word& name,
368  const dictionary& dict,
369  const label index,
370  const polyBoundaryMesh& bm,
371  const word& patchType
372 )
373 :
374  cyclicAMIPolyPatch(name, dict, index, bm, patchType),
375  nonOverlapPatchName_(dict.lookup("nonOverlapPatch")),
376  nonOverlapPatchID_(-1),
377  srcMask_(),
378  tgtMask_(),
379  updated_(false)
380 {
381  AMIRequireMatch_ = false;
382 
383  if (nonOverlapPatchName_ == name)
384  {
386  << "Non-overlapping patch name " << nonOverlapPatchName_
387  << " cannot be the same as this patch " << name
388  << exit(FatalIOError);
389  }
390 
391  // Non-overlapping patch might not be valid yet so cannot determine
392  // associated patchID
393 }
394 
395 
397 (
398  const cyclicACMIPolyPatch& pp,
399  const polyBoundaryMesh& bm
400 )
401 :
402  cyclicAMIPolyPatch(pp, bm),
403  nonOverlapPatchName_(pp.nonOverlapPatchName_),
404  nonOverlapPatchID_(-1),
405  srcMask_(),
406  tgtMask_(),
407  updated_(false)
408 {
409  AMIRequireMatch_ = false;
410 
411  // Non-overlapping patch might not be valid yet so cannot determine
412  // associated patchID
413 }
414 
415 
417 (
418  const cyclicACMIPolyPatch& pp,
419  const polyBoundaryMesh& bm,
420  const label index,
421  const label newSize,
422  const label newStart,
423  const word& nbrPatchName,
424  const word& nonOverlapPatchName
425 )
426 :
427  cyclicAMIPolyPatch(pp, bm, index, newSize, newStart, nbrPatchName),
428  nonOverlapPatchName_(nonOverlapPatchName),
429  nonOverlapPatchID_(-1),
430  srcMask_(),
431  tgtMask_(),
432  updated_(false)
433 {
434  AMIRequireMatch_ = false;
435 
436  if (nonOverlapPatchName_ == name())
437  {
439  << "Non-overlapping patch name " << nonOverlapPatchName_
440  << " cannot be the same as this patch " << name()
441  << exit(FatalError);
442  }
443 
444  // Non-overlapping patch might not be valid yet so cannot determine
445  // associated patchID
446 }
447 
448 
450 (
451  const cyclicACMIPolyPatch& pp,
452  const polyBoundaryMesh& bm,
453  const label index,
454  const labelUList& mapAddressing,
455  const label newStart
456 )
457 :
458  cyclicAMIPolyPatch(pp, bm, index, mapAddressing, newStart),
459  nonOverlapPatchName_(pp.nonOverlapPatchName_),
460  nonOverlapPatchID_(-1),
461  srcMask_(),
462  tgtMask_(),
463  updated_(false)
464 {
465  AMIRequireMatch_ = false;
466 }
467 
468 
469 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
470 
472 {}
473 
474 
475 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
476 
478 {
479  const polyPatch& pp = this->boundaryMesh()[neighbPatchID()];
480  return refCast<const cyclicACMIPolyPatch>(pp);
481 }
482 
483 
485 {
486  if (nonOverlapPatchID_ == -1)
487  {
488  nonOverlapPatchID_ =
489  this->boundaryMesh().findPatchID(nonOverlapPatchName_);
490 
491  if (nonOverlapPatchID_ == -1)
492  {
494  << "Illegal non-overlapping patch name " << nonOverlapPatchName_
495  << nl << "Valid patch names are "
496  << this->boundaryMesh().names()
497  << exit(FatalError);
498  }
499 
500  if (nonOverlapPatchID_ < index())
501  {
503  << "Boundary ordering error: " << type()
504  << " patch must be defined prior to its non-overlapping patch"
505  << nl
506  << type() << " patch: " << name() << ", ID:" << index() << nl
507  << "Non-overlap patch: " << nonOverlapPatchName_
508  << ", ID:" << nonOverlapPatchID_ << nl
509  << exit(FatalError);
510  }
511 
512  const polyPatch& noPp = this->boundaryMesh()[nonOverlapPatchID_];
513 
514  bool ok = true;
515 
516  if (size() == noPp.size())
517  {
518  const scalarField magSf(mag(faceAreas()));
519  const scalarField noMagSf(mag(noPp.faceAreas()));
520 
521  forAll(magSf, facei)
522  {
523  scalar ratio = mag(magSf[facei]/(noMagSf[facei] + ROOTVSMALL));
524 
525  if (ratio - 1 > tolerance_)
526  {
527  ok = false;
528  break;
529  }
530  }
531  }
532  else
533  {
534  ok = false;
535  }
536 
537  if (!ok)
538  {
540  << "Inconsistent ACMI patches " << name() << " and "
541  << noPp.name() << ". Patches should have identical topology"
542  << exit(FatalError);
543  }
544  }
545 
546  return nonOverlapPatchID_;
547 }
548 
549 
551 (
552  PstreamBuffers& pBufs,
553  const primitivePatch& pp
554 ) const
555 {
557 }
558 
559 
561 (
562  PstreamBuffers& pBufs,
563  const primitivePatch& pp,
565  labelList& rotation
566 ) const
567 {
568  return cyclicAMIPolyPatch::order(pBufs, pp, faceMap, rotation);
569 }
570 
571 
573 {
575 
576  os.writeEntry("nonOverlapPatch", nonOverlapPatchName_);
577 }
578 
579 
580 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::AMIInterpolation::interpolationMethod
interpolationMethod
Enumeration specifying interpolation method.
Definition: AMIInterpolation.H:90
Foam::cyclicAMIPolyPatch::order
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
Definition: cyclicAMIPolyPatch.C:980
Foam::boundaryMesh::mesh
const bMesh & mesh() const
Definition: boundaryMesh.H:205
Foam::cyclicACMIPolyPatch::~cyclicACMIPolyPatch
virtual ~cyclicACMIPolyPatch()
Destructor.
Definition: cyclicACMIPolyPatch.C:471
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::cyclicACMIPolyPatch::calcGeometry
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
Definition: cyclicACMIPolyPatch.C:253
Foam::cyclicACMIPolyPatch::movePoints
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
Definition: cyclicACMIPolyPatch.C:283
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeFast.C:94
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
cp
const volScalarField & cp
Definition: setRegionSolidFields.H:8
Foam::cyclicAMIPolyPatch::initOrder
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
Definition: cyclicAMIPolyPatch.C:972
Foam::AMIInterpolation::tgtWeightsSum
const scalarField & tgtWeightsSum() const
Definition: AMIInterpolationI.H:176
Foam::cyclicAMIPolyPatch::initGeometry
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: cyclicAMIPolyPatch.C:403
Foam::PstreamBuffers
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Definition: PstreamBuffers.H:88
Foam::AMIInterpolation::srcWeightsSum
const scalarField & srcWeightsSum() const
Definition: AMIInterpolationI.H:104
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:519
Foam::AMIInterpolation::tgtWeights
const scalarListList & tgtWeights() const
Return const access to target patch weights.
Definition: AMIInterpolationI.H:160
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)
Construct from (base coupled patch) components.
Definition: cyclicACMIPolyPatch.C:341
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:290
Foam::cyclicAMIPolyPatch::initMovePoints
virtual void initMovePoints(PstreamBuffers &pBufs, const pointField &)
Initialise the patches for moving points.
Definition: cyclicAMIPolyPatch.C:425
Foam::cyclicAMIPolyPatch::updateMesh
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: cyclicAMIPolyPatch.C:464
Foam::PrimitivePatch::faceAreas
const Field< PointType > & faceAreas() const
Return face area vectors for patch.
Definition: PrimitivePatch.C:537
Foam::SubField
SubField is a Field obtained as a section of another Field.
Definition: Field.H:64
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::cyclicACMIPolyPatch::updateMesh
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: cyclicACMIPolyPatch.C:306
Foam::cyclicACMIPolyPatch::write
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: cyclicACMIPolyPatch.C:572
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::cyclicAMIPolyPatch::calcGeometry
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
Definition: cyclicAMIPolyPatch.C:418
Foam::Field< scalar >
Foam::cyclicACMIPolyPatch::resetAMI
virtual void resetAMI(const AMIPatchToPatchInterpolation::interpolationMethod &AMIMethod=AMIPatchToPatchInterpolation::imFaceAreaWeight) const
Reset the AMI interpolator.
Definition: cyclicACMIPolyPatch.C:49
Foam::cyclicAMIPolyPatch::write
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: cyclicAMIPolyPatch.C:1044
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:66
Foam::cyclicACMIPolyPatch::srcMask
virtual const scalarField & srcMask() const
Return the mask/weighting for the source patch.
Definition: cyclicACMIPolyPatch.C:326
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::cyclicAMIPolyPatch::movePoints
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
Definition: cyclicAMIPolyPatch.C:444
Foam::cyclicACMIPolyPatch::clearGeom
virtual void clearGeom()
Clear geometry.
Definition: cyclicACMIPolyPatch.C:316
Foam::cyclicACMIPolyPatch::order
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
Definition: cyclicACMIPolyPatch.C:561
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:551
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:419
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::cyclicACMIPolyPatch::nonOverlapPatchID
virtual label nonOverlapPatchID() const
Non-overlapping patch ID.
Definition: cyclicACMIPolyPatch.C:484
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
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
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:264
Time.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::cyclicACMIPolyPatch::neighbPatch
virtual const cyclicACMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
Definition: cyclicACMIPolyPatch.C:477
Foam::cyclicAMIPolyPatch::initUpdateMesh
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: cyclicAMIPolyPatch.C:455
Foam::cyclicAMIPolyPatch::clearGeom
virtual void clearGeom()
Clear geometry.
Definition: cyclicAMIPolyPatch.C:470
Foam::List< scalarList >
Foam::polyPatch::faceAreas
const vectorField::subField faceAreas() const
Return face normals.
Definition: polyPatch.C:317
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::start
label ListType::const_reference const label start
Definition: ListOps.H:408
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:81
Foam::UList< label >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
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:219
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:61
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
Foam::cyclicACMIPolyPatch::tgtMask
virtual const scalarField & tgtMask() const
Return the mask/weighting for the target patch.
Definition: cyclicACMIPolyPatch.C:332
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:53
Foam::coupledPolyPatch::transformType
transformType
Definition: coupledPolyPatch.H:59
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::patchIdentifier::name
const word & name() const
Return the patch name.
Definition: patchIdentifier.H:109
Foam::AMIInterpolation::srcWeights
const scalarListList & srcWeights() const
Return const access to source patch weights.
Definition: AMIInterpolationI.H:88
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:90
Foam::cyclicACMIPolyPatch::initGeometry
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: cyclicACMIPolyPatch.C:237
Foam::cyclicACMIPolyPatch::initUpdateMesh
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: cyclicACMIPolyPatch.C:296
Foam::cyclicAMIPolyPatch
Cyclic patch for Arbitrary Mesh Interface (AMI)
Definition: cyclicAMIPolyPatch.H:53