cyclicPolyPatch.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "cyclicPolyPatch.H"
28 #include "polyBoundaryMesh.H"
29 #include "polyMesh.H"
30 #include "demandDrivenData.H"
31 #include "OFstream.H"
32 #include "matchPoints.H"
33 #include "EdgeMap.H"
34 #include "Time.H"
35 #include "transformField.H"
36 #include "SubField.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
43 
46 }
47 
48 
49 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50 
52 {
54 }
55 
56 
58 (
59  const primitivePatch& referPatch,
60  pointField& nbrCtrs,
61  vectorField& nbrAreas,
62  pointField& nbrCc
63 )
64 {}
65 
66 
68 {
69  static_cast<cyclicTransform&>(*this) =
71  (
72  name(),
73  faceCentres(),
74  faceAreas(),
75  *this,
76  nbrPatchName(),
77  nbrPatch().faceCentres(),
78  nbrPatch().faceAreas(),
79  nbrPatch(),
80  matchTolerance()
81  );
82 }
83 
84 
86 (
87  PstreamBuffers& pBufs,
88  const pointField& p
89 )
90 {
92 }
93 
94 
96 (
97  PstreamBuffers& pBufs,
98  const pointField& p
99 )
100 {
101  polyPatch::movePoints(pBufs, p);
102 }
103 
104 
106 {
108 }
109 
110 
112 {
113  polyPatch::topoChange(pBufs);
114  deleteDemandDrivenData(coupledPointsPtr_);
115  deleteDemandDrivenData(coupledEdgesPtr_);
116 }
117 
118 
120 {
121  polyPatch::rename(newNames);
122  nbrPatch().nbrPatchName_ = newNames[index()];
123 }
124 
125 
126 void Foam::cyclicPolyPatch::reorder(const labelUList& newToOldIndex)
127 {
128  polyPatch::reorder(newToOldIndex);
129  if (nbrPatchIndex_ != -1)
130  {
131  nbrPatchIndex_ = findIndex(newToOldIndex, nbrPatchIndex_);
132  }
133 }
134 
135 
136 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
137 
139 (
140  const word& name,
141  const label size,
142  const label start,
143  const label index,
144  const polyBoundaryMesh& bm,
145  const word& patchType
146 )
147 :
148  coupledPolyPatch(name, size, start, index, bm, patchType),
149  cyclicTransform(false),
150  nbrPatchName_(word::null),
151  nbrPatchIndex_(-1),
152  coupledPointsPtr_(nullptr),
153  coupledEdgesPtr_(nullptr),
154  ownToNbrOrderDataPtr_(nullptr),
155  ownToNbrCyclicOrderDataPtr_(nullptr),
156  ownToNbrDebugOrderDataPtr_(nullptr)
157 {}
158 
159 
161 (
162  const word& name,
163  const label size,
164  const label start,
165  const label index,
166  const polyBoundaryMesh& bm,
167  const word& patchType,
168  const word& nbrPatchName,
170 )
171 :
172  coupledPolyPatch(name, size, start, index, bm, patchType),
174  nbrPatchName_(nbrPatchName),
175  nbrPatchIndex_(-1),
176  coupledPointsPtr_(nullptr),
177  coupledEdgesPtr_(nullptr),
178  ownToNbrOrderDataPtr_(nullptr),
179  ownToNbrCyclicOrderDataPtr_(nullptr),
180  ownToNbrDebugOrderDataPtr_(nullptr)
181 {}
182 
183 
185 (
186  const word& name,
187  const dictionary& dict,
188  const label index,
189  const polyBoundaryMesh& bm,
190  const word& patchType,
191  const bool cyclicTransformDefaultIsNone
192 )
193 :
194  coupledPolyPatch(name, dict, index, bm, patchType),
195  cyclicTransform(dict, cyclicTransformDefaultIsNone),
196  nbrPatchName_(dict.lookupOrDefault("neighbourPatch", word::null)),
197  coupleGroup_(dict),
198  nbrPatchIndex_(-1),
199  coupledPointsPtr_(nullptr),
200  coupledEdgesPtr_(nullptr),
201  ownToNbrOrderDataPtr_(nullptr),
202  ownToNbrCyclicOrderDataPtr_(nullptr),
203  ownToNbrDebugOrderDataPtr_(nullptr)
204 {
205  if (nbrPatchName_ == word::null && !coupleGroup_.valid())
206  {
208  (
209  dict
210  ) << "No \"neighbourPatch\" provided." << endl
211  << exit(FatalIOError);
212  }
213 
214  if (nbrPatchName_ == name)
215  {
217  << "Neighbour patch name " << nbrPatchName_
218  << " cannot be the same as this patch " << name
219  << exit(FatalIOError);
220  }
221 
222  // Neighbour patch might not be valid yet so no transformation
223  // calculation possible.
224 }
225 
226 
228 (
229  const cyclicPolyPatch& pp,
230  const polyBoundaryMesh& bm
231 )
232 :
233  coupledPolyPatch(pp, bm),
234  cyclicTransform(pp),
235  nbrPatchName_(pp.nbrPatchName_),
236  coupleGroup_(pp.coupleGroup_),
237  nbrPatchIndex_(-1),
238  coupledPointsPtr_(nullptr),
239  coupledEdgesPtr_(nullptr),
240  ownToNbrOrderDataPtr_(nullptr),
241  ownToNbrCyclicOrderDataPtr_(nullptr),
242  ownToNbrDebugOrderDataPtr_(nullptr)
243 {
244  // Neighbour patch might not be valid yet so no transformation
245  // calculation possible.
246 }
247 
248 
250 (
251  const cyclicPolyPatch& pp,
252  const polyBoundaryMesh& bm,
253  const label index,
254  const label newSize,
255  const label newStart,
256  const word& neiName
257 )
258 :
259  coupledPolyPatch(pp, bm, index, newSize, newStart),
260  cyclicTransform(pp),
261  nbrPatchName_(neiName),
262  coupleGroup_(pp.coupleGroup_),
263  nbrPatchIndex_(-1),
264  coupledPointsPtr_(nullptr),
265  coupledEdgesPtr_(nullptr),
266  ownToNbrOrderDataPtr_(nullptr),
267  ownToNbrCyclicOrderDataPtr_(nullptr),
268  ownToNbrDebugOrderDataPtr_(nullptr)
269 {
270  if (neiName == name())
271  {
273  << "Neighbour patch name " << neiName
274  << " cannot be the same as this patch " << name()
275  << exit(FatalError);
276  }
277 
278  // Neighbour patch might not be valid yet so no transformation
279  // calculation possible.
280 }
281 
282 
283 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
284 
286 {
287  deleteDemandDrivenData(coupledPointsPtr_);
288  deleteDemandDrivenData(coupledEdgesPtr_);
289 }
290 
291 
292 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
293 
295 {
296  if (nbrPatchName_.empty())
297  {
298  // Try and use patchGroup to find neighbourPatch and neighbourRegion
299  label patchID = coupleGroup_.findOtherPatchID(*this);
300 
301  nbrPatchName_ = boundaryMesh()[patchID].name();
302  }
303  return nbrPatchName_;
304 }
305 
306 
308 {
309  if (nbrPatchIndex_ == -1)
310  {
311  nbrPatchIndex_ = this->boundaryMesh().findIndex(nbrPatchName());
312 
313  if (nbrPatchIndex_ == -1)
314  {
316  << "Illegal neighbourPatch name " << nbrPatchName()
317  << endl << "Valid patch names are "
318  << this->boundaryMesh().names()
319  << exit(FatalError);
320  }
321 
322  // Check that it is a cyclic
323  const cyclicPolyPatch& nbrPatch = refCast<const cyclicPolyPatch>
324  (
325  this->boundaryMesh()[nbrPatchIndex_]
326  );
327 
328  if (nbrPatch.nbrPatchName() != name())
329  {
331  << "Patch " << name()
332  << " specifies neighbour patch " << nbrPatchName()
333  << endl << " but that in return specifies "
334  << nbrPatch.nbrPatchName()
335  << endl;
336  }
337  }
338  return nbrPatchIndex_;
339 }
340 
341 
343 {
344  if (!coupledPointsPtr_)
345  {
346  const faceList& nbrLocalFaces = nbrPatch().localFaces();
347  const labelList& nbrMeshPoints = nbrPatch().meshPoints();
348 
349  // Now all we know is that relative face index in *this is same
350  // as coupled face in nbrPatch and also that the 0th vertex
351  // corresponds.
352 
353  // From local point to nbrPatch or -1.
354  labelList coupledPoint(nPoints(), -1);
355 
356  forAll(*this, patchFacei)
357  {
358  const face& fA = localFaces()[patchFacei];
359  const face& fB = nbrLocalFaces[patchFacei];
360 
361  forAll(fA, indexA)
362  {
363  label patchPointA = fA[indexA];
364 
365  if (coupledPoint[patchPointA] == -1)
366  {
367  label indexB = (fB.size() - indexA) % fB.size();
368 
369  // Filter out points on wedge axis
370  if (meshPoints()[patchPointA] != nbrMeshPoints[fB[indexB]])
371  {
372  coupledPoint[patchPointA] = fB[indexB];
373  }
374  }
375  }
376  }
377 
378  coupledPointsPtr_ = new edgeList(nPoints());
379  edgeList& connected = *coupledPointsPtr_;
380 
381  // Extract coupled points.
382  label connectedI = 0;
383 
384  forAll(coupledPoint, i)
385  {
386  if (coupledPoint[i] != -1)
387  {
388  connected[connectedI++] = edge(i, coupledPoint[i]);
389  }
390  }
391 
392  connected.setSize(connectedI);
393 
394  if (debug)
395  {
396  OFstream str
397  (
398  boundaryMesh().mesh().time().path()
399  /name() + "_coupledPoints.obj"
400  );
401  label vertI = 0;
402 
403  Pout<< "Writing file " << str.name() << " with coordinates of "
404  << "coupled points" << endl;
405 
406  forAll(connected, i)
407  {
408  const point& a = points()[meshPoints()[connected[i][0]]];
409  const point& b = points()[nbrMeshPoints[connected[i][1]]];
410 
411  str<< "v " << a.x() << ' ' << a.y() << ' ' << a.z() << nl;
412  str<< "v " << b.x() << ' ' << b.y() << ' ' << b.z() << nl;
413  vertI += 2;
414 
415  str<< "l " << vertI-1 << ' ' << vertI << nl;
416  }
417  }
418  }
419  return *coupledPointsPtr_;
420 }
421 
422 
424 {
425  if (!coupledEdgesPtr_)
426  {
427  const edgeList& pointCouples = coupledPoints();
428 
429  // Build map from points on *this (A) to points on neighbourpatch (B)
430  Map<label> aToB(2*pointCouples.size());
431 
432  forAll(pointCouples, i)
433  {
434  const edge& e = pointCouples[i];
435 
436  aToB.insert(e[0], e[1]);
437  }
438 
439  // Map from edge on A to points (in B indices)
440  EdgeMap<label> edgeMap(nEdges());
441 
442  forAll(*this, patchFacei)
443  {
444  const labelList& fEdges = faceEdges()[patchFacei];
445 
446  forAll(fEdges, i)
447  {
448  label edgeI = fEdges[i];
449 
450  const edge& e = edges()[edgeI];
451 
452  // Convert edge end points to corresponding points on B side.
453  Map<label>::const_iterator fnd0 = aToB.find(e[0]);
454  if (fnd0 != aToB.end())
455  {
456  Map<label>::const_iterator fnd1 = aToB.find(e[1]);
457  if (fnd1 != aToB.end())
458  {
459  edgeMap.insert(edge(fnd0(), fnd1()), edgeI);
460  }
461  }
462  }
463  }
464 
465  // Use the edgeMap to get the edges on the B side.
466 
467  const cyclicPolyPatch& nbrPatch = this->nbrPatch();
468  const labelList& nbrMp = nbrPatch.meshPoints();
469  const labelList& mp = meshPoints();
470 
471 
472 
473  coupledEdgesPtr_ = new edgeList(edgeMap.size());
474  edgeList& coupledEdges = *coupledEdgesPtr_;
475  label coupleI = 0;
476 
477  forAll(nbrPatch, patchFacei)
478  {
479  const labelList& fEdges = nbrPatch.faceEdges()[patchFacei];
480 
481  forAll(fEdges, i)
482  {
483  label edgeI = fEdges[i];
484 
485  const edge& e = nbrPatch.edges()[edgeI];
486 
487  // Look up A edge from HashTable.
488  EdgeMap<label>::iterator iter = edgeMap.find(e);
489 
490  if (iter != edgeMap.end())
491  {
492  label edgeA = iter();
493  const edge& eA = edges()[edgeA];
494 
495  // Store correspondence. Filter out edges on wedge axis.
496  if
497  (
498  edge(mp[eA[0]], mp[eA[1]])
499  != edge(nbrMp[e[0]], nbrMp[e[1]])
500  )
501  {
502  coupledEdges[coupleI++] = edge(edgeA, edgeI);
503  }
504 
505  // Remove so we build unique list
506  edgeMap.erase(iter);
507  }
508  }
509  }
510  coupledEdges.setSize(coupleI);
511 
512 
513  // Some checks
514 
515  forAll(coupledEdges, i)
516  {
517  const edge& e = coupledEdges[i];
518 
519  if (e[0] < 0 || e[1] < 0)
520  {
522  << "Problem : at position " << i
523  << " illegal couple:" << e
524  << abort(FatalError);
525  }
526  }
527 
528  if (debug)
529  {
530  OFstream str
531  (
532  boundaryMesh().mesh().time().path()
533  /name() + "_coupledEdges.obj"
534  );
535  label vertI = 0;
536 
537  Pout<< "Writing file " << str.name() << " with centres of "
538  << "coupled edges" << endl;
539 
540  forAll(coupledEdges, i)
541  {
542  const edge& e = coupledEdges[i];
543 
544  const point& a = edges()[e[0]].centre(localPoints());
545  const point& b = nbrPatch.edges()[e[1]].centre
546  (
547  nbrPatch.localPoints()
548  );
549 
550  str<< "v " << a.x() << ' ' << a.y() << ' ' << a.z() << nl;
551  str<< "v " << b.x() << ' ' << b.y() << ' ' << b.z() << nl;
552  vertI += 2;
553 
554  str<< "l " << vertI-1 << ' ' << vertI << nl;
555  }
556  }
557  }
558  return *coupledEdgesPtr_;
559 }
560 
561 
563 (
565  const primitivePatch& pp
566 ) const
567 {
568  if (pp.empty())
569  {
570  return;
571  }
572 
573  if (owner())
574  {
575  ownToNbrOrderDataPtr_ = new ownToNbrOrderData();
576  if (coupledPolyPatch::debug)
577  {
578  ownToNbrDebugOrderDataPtr_ = new ownToNbrDebugOrderData();
579  }
580 
582  (
583  ownToNbrOrderDataPtr_(),
584  ownToNbrDebugOrderDataPtr_,
585  pp
586  );
587 
588  const scalarField magAreas(mag(pp.faceAreas()));
589 
590  ownToNbrCyclicOrderDataPtr_ = new ownToNbrCyclicOrderData();
591  ownToNbrCyclicOrderDataPtr_->ctr =
592  sum(pp.faceCentres()*magAreas)/sum(magAreas);
593  ownToNbrCyclicOrderDataPtr_->area = sum(pp.faceAreas());
594  }
595 }
596 
597 
599 (
600  PstreamBuffers& pBufs,
601  const primitivePatch& pp,
603  labelList& rotation
604 ) const
605 {
606  if (pp.empty())
607  {
608  return false;
609  }
610 
611  ownToNbrOrderData ownToNbr;
612  autoPtr<ownToNbrDebugOrderData> ownToNbrDebugPtr(nullptr);
613 
614  if (!owner())
615  {
616  ownToNbr = nbrPatch().ownToNbrOrderDataPtr_();
617  ownToNbrDebugPtr = nbrPatch().ownToNbrDebugOrderDataPtr_;
618 
619  cyclicTransform ct
620  (
621  name(),
622  pp.faceCentres(),
623  pp.faceAreas(),
624  *this,
625  nbrPatchName(),
626  pointField(1, nbrPatch().ownToNbrCyclicOrderDataPtr_->ctr),
627  vectorField(1, nbrPatch().ownToNbrCyclicOrderDataPtr_->area),
628  nbrPatch(),
629  matchTolerance()
630  );
631 
632  ownToNbr.transform(ct.transform());
633  if (ownToNbrDebugPtr.valid())
634  {
635  ownToNbrDebugPtr->transform(ct.transform());
636  }
637  }
638 
639  return
641  (
642  ownToNbr,
643  ownToNbrDebugPtr,
644  pp,
645  faceMap,
646  rotation
647  );
648 }
649 
650 
652 {
654 
655  if (!nbrPatchName_.empty())
656  {
657  writeEntry(os, "neighbourPatch", nbrPatchName_);
658  }
659 
660  coupleGroup_.write(os);
661 
663 }
664 
665 
666 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
bool erase(const iterator &)
Erase a hashedEntry specified by given iterator.
Definition: HashTable.C:445
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:167
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
Output to file stream.
Definition: OFstream.H:86
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:120
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
const Field< PointType > & faceAreas() const
Return face areas for patch.
const labelListList & faceEdges() const
Return face-edge addressing.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
const Field< PointType > & faceCentres() const
Return face centres for patch.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & x() const
Definition: VectorI.H:75
const Vector< Cmpt > & centre(const Foam::List< Vector< Cmpt >> &) const
Return *this (used for point which is a typedef to Vector<scalar>.
Definition: VectorI.H:116
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
bool valid() const
Is a valid patchGroup.
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
virtual bool order(const ownToNbrOrderData &ownToNbr, const autoPtr< ownToNbrDebugOrderData > &ownToNbrDebugPtr, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for the given primitivePatch.
virtual void initOrder(ownToNbrOrderData &ownToNbr, autoPtr< ownToNbrDebugOrderData > &ownToNbrDebugPtr, const primitivePatch &) const
Initialise ordering for the given primitivePatch. Fills the.
Cyclic plane patch.
cyclicPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType)
Construct from components.
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
virtual void rename(const wordList &newNames)
Reset the patch name.
const edgeList & coupledEdges() const
Return connected edges (from patch local to neighbour patch local).
const word & nbrPatchName() const
Neighbour patch name.
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
const edgeList & coupledPoints() const
Return connected points (from patch local to neighbour patch local)
virtual void reorder(const labelUList &newToOldIndex)
Reset the patch index.
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialise ordering for primitivePatch. Does not.
virtual void initCalcGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
virtual label nbrPatchIndex() const
Neighbour patchID.
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
virtual void movePoints(PstreamBuffers &, const pointField &)
Correct patches after moving points.
virtual void initTopoChange(PstreamBuffers &)
Initialise the update of the patch topology.
virtual ~cyclicPolyPatch()
Destructor.
virtual void topoChange(PstreamBuffers &)
Update of the patch topology.
Cyclic plane transformation.
void write(Ostream &os) const
Write the data to a dictionary.
const transformer & transform() const
Return transformation between the coupled patches.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
const word & name() const
Return name.
Foam::polyBoundaryMesh.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
virtual void rename(const wordList &newNames)
Reset the patch name.
Definition: polyPatch.C:82
virtual void reorder(const labelUList &newToOldIndex)
Reset the patch index.
Definition: polyPatch.C:88
virtual void movePoints(const pointField &p)
Correct patches after moving points.
Definition: polyPatch.C:57
virtual void initCalcGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: polyPatch.H:97
virtual void initTopoChange(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: polyPatch.H:115
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
Definition: polyPatch.H:108
virtual void topoChange(PstreamBuffers &)
Update of the patch topology.
Definition: polyPatch.C:69
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
Template functions to aid in the implementation of demand driven data.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const pointField & points
label nPoints
volScalarField & b
Definition: createFields.H:25
Determine correspondence between points. See below.
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar mp
Proton mass.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const doubleScalar e
Definition: doubleScalar.H:106
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:59
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
void deleteDemandDrivenData(DataType *&dataPtr)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:504
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
IOerror FatalIOError
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
List< edge > edgeList
Definition: edgeList.H:38
static const char nl
Definition: Ostream.H:266
dictionary dict
volScalarField & p
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
Data to pass from owner.initOrder to nbr.order if debugging.
Data to pass from owner.initOrder to nbr.order.
void transform(const transformer &tr)
Spatial transformation functions for primitive fields.