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-2022 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 #include "unitConversion.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
44 
47 }
48 
49 
50 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
51 
53 {
55 }
56 
57 
59 (
60  const primitivePatch& referPatch,
61  pointField& nbrCtrs,
62  vectorField& nbrAreas,
63  pointField& nbrCc
64 )
65 {}
66 
67 
69 {
70  static_cast<cyclicTransform&>(*this) =
72  (
73  name(),
74  faceCentres(),
75  faceAreas(),
76  *this,
77  nbrPatchName(),
78  nbrPatch().faceCentres(),
79  nbrPatch().faceAreas(),
80  nbrPatch(),
81  matchTolerance()
82  );
83 }
84 
85 
87 (
88  PstreamBuffers& pBufs,
89  const pointField& p
90 )
91 {
93 }
94 
95 
97 (
98  PstreamBuffers& pBufs,
99  const pointField& p
100 )
101 {
102  polyPatch::movePoints(pBufs, p);
103 }
104 
105 
107 {
109 }
110 
111 
113 {
114  polyPatch::topoChange(pBufs);
115  deleteDemandDrivenData(coupledPointsPtr_);
116  deleteDemandDrivenData(coupledEdgesPtr_);
117 }
118 
119 
121 {
122  polyPatch::rename(newNames);
123  nbrPatch().nbrPatchName_ = newNames[index()];
124 }
125 
126 
127 void Foam::cyclicPolyPatch::reorder(const labelUList& newToOldIndex)
128 {
129  polyPatch::reorder(newToOldIndex);
130  if (nbrPatchID_ != -1)
131  {
132  nbrPatchID_ = findIndex(newToOldIndex, nbrPatchID_);
133  }
134 }
135 
136 
137 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
138 
140 (
141  const word& name,
142  const label size,
143  const label start,
144  const label index,
145  const polyBoundaryMesh& bm,
146  const word& patchType
147 )
148 :
149  coupledPolyPatch(name, size, start, index, bm, patchType),
150  cyclicTransform(false),
151  nbrPatchName_(word::null),
152  nbrPatchID_(-1),
153  coupledPointsPtr_(nullptr),
154  coupledEdgesPtr_(nullptr),
155  ownToNbrOrderDataPtr_(nullptr),
156  ownToNbrCyclicOrderDataPtr_(nullptr),
157  ownToNbrDebugOrderDataPtr_(nullptr)
158 {}
159 
160 
162 (
163  const word& name,
164  const label size,
165  const label start,
166  const label index,
167  const polyBoundaryMesh& bm,
168  const word& patchType,
169  const word& nbrPatchName,
171 )
172 :
173  coupledPolyPatch(name, size, start, index, bm, patchType),
175  nbrPatchName_(nbrPatchName),
176  nbrPatchID_(-1),
177  coupledPointsPtr_(nullptr),
178  coupledEdgesPtr_(nullptr),
179  ownToNbrOrderDataPtr_(nullptr),
180  ownToNbrCyclicOrderDataPtr_(nullptr),
181  ownToNbrDebugOrderDataPtr_(nullptr)
182 {}
183 
184 
186 (
187  const word& name,
188  const dictionary& dict,
189  const label index,
190  const polyBoundaryMesh& bm,
191  const word& patchType,
192  const bool cyclicTransformDefaultIsNone
193 )
194 :
195  coupledPolyPatch(name, dict, index, bm, patchType),
196  cyclicTransform(dict, cyclicTransformDefaultIsNone),
197  nbrPatchName_(dict.lookupOrDefault("neighbourPatch", word::null)),
198  coupleGroup_(dict),
199  nbrPatchID_(-1),
200  coupledPointsPtr_(nullptr),
201  coupledEdgesPtr_(nullptr),
202  ownToNbrOrderDataPtr_(nullptr),
203  ownToNbrCyclicOrderDataPtr_(nullptr),
204  ownToNbrDebugOrderDataPtr_(nullptr)
205 {
206  if (nbrPatchName_ == word::null && !coupleGroup_.valid())
207  {
209  (
210  dict
211  ) << "No \"neighbourPatch\" provided." << endl
212  << exit(FatalIOError);
213  }
214 
215  if (nbrPatchName_ == name)
216  {
218  << "Neighbour patch name " << nbrPatchName_
219  << " cannot be the same as this patch " << name
220  << exit(FatalIOError);
221  }
222 
223  // Neighbour patch might not be valid yet so no transformation
224  // calculation possible.
225 }
226 
227 
229 (
230  const cyclicPolyPatch& pp,
231  const polyBoundaryMesh& bm
232 )
233 :
234  coupledPolyPatch(pp, bm),
235  cyclicTransform(pp),
236  nbrPatchName_(pp.nbrPatchName_),
237  coupleGroup_(pp.coupleGroup_),
238  nbrPatchID_(-1),
239  coupledPointsPtr_(nullptr),
240  coupledEdgesPtr_(nullptr),
241  ownToNbrOrderDataPtr_(nullptr),
242  ownToNbrCyclicOrderDataPtr_(nullptr),
243  ownToNbrDebugOrderDataPtr_(nullptr)
244 {
245  // Neighbour patch might not be valid yet so no transformation
246  // calculation possible.
247 }
248 
249 
251 (
252  const cyclicPolyPatch& pp,
253  const polyBoundaryMesh& bm,
254  const label index,
255  const label newSize,
256  const label newStart,
257  const word& neiName
258 )
259 :
260  coupledPolyPatch(pp, bm, index, newSize, newStart),
261  cyclicTransform(pp),
262  nbrPatchName_(neiName),
263  coupleGroup_(pp.coupleGroup_),
264  nbrPatchID_(-1),
265  coupledPointsPtr_(nullptr),
266  coupledEdgesPtr_(nullptr),
267  ownToNbrOrderDataPtr_(nullptr),
268  ownToNbrCyclicOrderDataPtr_(nullptr),
269  ownToNbrDebugOrderDataPtr_(nullptr)
270 {
271  if (neiName == name())
272  {
274  << "Neighbour patch name " << neiName
275  << " cannot be the same as this patch " << name()
276  << exit(FatalError);
277  }
278 
279  // Neighbour patch might not be valid yet so no transformation
280  // calculation possible.
281 }
282 
283 
284 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
285 
287 {
288  deleteDemandDrivenData(coupledPointsPtr_);
289  deleteDemandDrivenData(coupledEdgesPtr_);
290 }
291 
292 
293 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
294 
296 {
297  if (nbrPatchName_.empty())
298  {
299  // Try and use patchGroup to find neighbourPatch and neighbourRegion
300  label patchID = coupleGroup_.findOtherPatchID(*this);
301 
302  nbrPatchName_ = boundaryMesh()[patchID].name();
303  }
304  return nbrPatchName_;
305 }
306 
307 
309 {
310  if (nbrPatchID_ == -1)
311  {
312  nbrPatchID_ = this->boundaryMesh().findPatchID(nbrPatchName());
313 
314  if (nbrPatchID_ == -1)
315  {
317  << "Illegal neighbourPatch name " << nbrPatchName()
318  << endl << "Valid patch names are "
319  << this->boundaryMesh().names()
320  << exit(FatalError);
321  }
322 
323  // Check that it is a cyclic
324  const cyclicPolyPatch& nbrPatch = refCast<const cyclicPolyPatch>
325  (
326  this->boundaryMesh()[nbrPatchID_]
327  );
328 
329  if (nbrPatch.nbrPatchName() != name())
330  {
332  << "Patch " << name()
333  << " specifies neighbour patch " << nbrPatchName()
334  << endl << " but that in return specifies "
335  << nbrPatch.nbrPatchName()
336  << endl;
337  }
338  }
339  return nbrPatchID_;
340 }
341 
342 
344 {
345  if (!coupledPointsPtr_)
346  {
347  const faceList& nbrLocalFaces = nbrPatch().localFaces();
348  const labelList& nbrMeshPoints = nbrPatch().meshPoints();
349 
350  // Now all we know is that relative face index in *this is same
351  // as coupled face in nbrPatch and also that the 0th vertex
352  // corresponds.
353 
354  // From local point to nbrPatch or -1.
355  labelList coupledPoint(nPoints(), -1);
356 
357  forAll(*this, patchFacei)
358  {
359  const face& fA = localFaces()[patchFacei];
360  const face& fB = nbrLocalFaces[patchFacei];
361 
362  forAll(fA, indexA)
363  {
364  label patchPointA = fA[indexA];
365 
366  if (coupledPoint[patchPointA] == -1)
367  {
368  label indexB = (fB.size() - indexA) % fB.size();
369 
370  // Filter out points on wedge axis
371  if (meshPoints()[patchPointA] != nbrMeshPoints[fB[indexB]])
372  {
373  coupledPoint[patchPointA] = fB[indexB];
374  }
375  }
376  }
377  }
378 
379  coupledPointsPtr_ = new edgeList(nPoints());
380  edgeList& connected = *coupledPointsPtr_;
381 
382  // Extract coupled points.
383  label connectedI = 0;
384 
385  forAll(coupledPoint, i)
386  {
387  if (coupledPoint[i] != -1)
388  {
389  connected[connectedI++] = edge(i, coupledPoint[i]);
390  }
391  }
392 
393  connected.setSize(connectedI);
394 
395  if (debug)
396  {
397  OFstream str
398  (
399  boundaryMesh().mesh().time().path()
400  /name() + "_coupledPoints.obj"
401  );
402  label vertI = 0;
403 
404  Pout<< "Writing file " << str.name() << " with coordinates of "
405  << "coupled points" << endl;
406 
407  forAll(connected, i)
408  {
409  const point& a = points()[meshPoints()[connected[i][0]]];
410  const point& b = points()[nbrMeshPoints[connected[i][1]]];
411 
412  str<< "v " << a.x() << ' ' << a.y() << ' ' << a.z() << nl;
413  str<< "v " << b.x() << ' ' << b.y() << ' ' << b.z() << nl;
414  vertI += 2;
415 
416  str<< "l " << vertI-1 << ' ' << vertI << nl;
417  }
418  }
419  }
420  return *coupledPointsPtr_;
421 }
422 
423 
425 {
426  if (!coupledEdgesPtr_)
427  {
428  const edgeList& pointCouples = coupledPoints();
429 
430  // Build map from points on *this (A) to points on neighbourpatch (B)
431  Map<label> aToB(2*pointCouples.size());
432 
433  forAll(pointCouples, i)
434  {
435  const edge& e = pointCouples[i];
436 
437  aToB.insert(e[0], e[1]);
438  }
439 
440  // Map from edge on A to points (in B indices)
441  EdgeMap<label> edgeMap(nEdges());
442 
443  forAll(*this, patchFacei)
444  {
445  const labelList& fEdges = faceEdges()[patchFacei];
446 
447  forAll(fEdges, i)
448  {
449  label edgeI = fEdges[i];
450 
451  const edge& e = edges()[edgeI];
452 
453  // Convert edge end points to corresponding points on B side.
454  Map<label>::const_iterator fnd0 = aToB.find(e[0]);
455  if (fnd0 != aToB.end())
456  {
457  Map<label>::const_iterator fnd1 = aToB.find(e[1]);
458  if (fnd1 != aToB.end())
459  {
460  edgeMap.insert(edge(fnd0(), fnd1()), edgeI);
461  }
462  }
463  }
464  }
465 
466  // Use the edgeMap to get the edges on the B side.
467 
468  const cyclicPolyPatch& nbrPatch = this->nbrPatch();
469  const labelList& nbrMp = nbrPatch.meshPoints();
470  const labelList& mp = meshPoints();
471 
472 
473 
474  coupledEdgesPtr_ = new edgeList(edgeMap.size());
475  edgeList& coupledEdges = *coupledEdgesPtr_;
476  label coupleI = 0;
477 
478  forAll(nbrPatch, patchFacei)
479  {
480  const labelList& fEdges = nbrPatch.faceEdges()[patchFacei];
481 
482  forAll(fEdges, i)
483  {
484  label edgeI = fEdges[i];
485 
486  const edge& e = nbrPatch.edges()[edgeI];
487 
488  // Look up A edge from HashTable.
489  EdgeMap<label>::iterator iter = edgeMap.find(e);
490 
491  if (iter != edgeMap.end())
492  {
493  label edgeA = iter();
494  const edge& eA = edges()[edgeA];
495 
496  // Store correspondence. Filter out edges on wedge axis.
497  if
498  (
499  edge(mp[eA[0]], mp[eA[1]])
500  != edge(nbrMp[e[0]], nbrMp[e[1]])
501  )
502  {
503  coupledEdges[coupleI++] = edge(edgeA, edgeI);
504  }
505 
506  // Remove so we build unique list
507  edgeMap.erase(iter);
508  }
509  }
510  }
511  coupledEdges.setSize(coupleI);
512 
513 
514  // Some checks
515 
516  forAll(coupledEdges, i)
517  {
518  const edge& e = coupledEdges[i];
519 
520  if (e[0] < 0 || e[1] < 0)
521  {
523  << "Problem : at position " << i
524  << " illegal couple:" << e
525  << abort(FatalError);
526  }
527  }
528 
529  if (debug)
530  {
531  OFstream str
532  (
533  boundaryMesh().mesh().time().path()
534  /name() + "_coupledEdges.obj"
535  );
536  label vertI = 0;
537 
538  Pout<< "Writing file " << str.name() << " with centres of "
539  << "coupled edges" << endl;
540 
541  forAll(coupledEdges, i)
542  {
543  const edge& e = coupledEdges[i];
544 
545  const point& a = edges()[e[0]].centre(localPoints());
546  const point& b = nbrPatch.edges()[e[1]].centre
547  (
548  nbrPatch.localPoints()
549  );
550 
551  str<< "v " << a.x() << ' ' << a.y() << ' ' << a.z() << nl;
552  str<< "v " << b.x() << ' ' << b.y() << ' ' << b.z() << nl;
553  vertI += 2;
554 
555  str<< "l " << vertI-1 << ' ' << vertI << nl;
556  }
557  }
558  }
559  return *coupledEdgesPtr_;
560 }
561 
562 
564 (
566  const primitivePatch& pp
567 ) const
568 {
569  if (pp.empty())
570  {
571  return;
572  }
573 
574  if (owner())
575  {
576  ownToNbrOrderDataPtr_ = new ownToNbrOrderData();
577  if (coupledPolyPatch::debug)
578  {
579  ownToNbrDebugOrderDataPtr_ = new ownToNbrDebugOrderData();
580  }
581 
583  (
584  ownToNbrOrderDataPtr_(),
585  ownToNbrDebugOrderDataPtr_,
586  pp
587  );
588 
589  const scalarField magAreas(mag(pp.faceAreas()));
590 
591  ownToNbrCyclicOrderDataPtr_ = new ownToNbrCyclicOrderData();
592  ownToNbrCyclicOrderDataPtr_->ctr =
593  sum(pp.faceCentres()*magAreas)/sum(magAreas);
594  ownToNbrCyclicOrderDataPtr_->area = sum(pp.faceAreas());
595  }
596 }
597 
598 
600 (
601  PstreamBuffers& pBufs,
602  const primitivePatch& pp,
604  labelList& rotation
605 ) const
606 {
607  if (pp.empty())
608  {
609  return false;
610  }
611 
612  ownToNbrOrderData ownToNbr;
613  autoPtr<ownToNbrDebugOrderData> ownToNbrDebugPtr(nullptr);
614 
615  if (!owner())
616  {
617  ownToNbr = nbrPatch().ownToNbrOrderDataPtr_();
618  ownToNbrDebugPtr = nbrPatch().ownToNbrDebugOrderDataPtr_;
619 
620  cyclicTransform ct
621  (
622  name(),
623  pp.faceCentres(),
624  pp.faceAreas(),
625  *this,
626  nbrPatchName(),
627  pointField(1, nbrPatch().ownToNbrCyclicOrderDataPtr_->ctr),
628  vectorField(1, nbrPatch().ownToNbrCyclicOrderDataPtr_->area),
629  nbrPatch(),
630  matchTolerance()
631  );
632 
633  ownToNbr.transform(ct.transform());
634  if (ownToNbrDebugPtr.valid())
635  {
636  ownToNbrDebugPtr->transform(ct.transform());
637  }
638  }
639 
640  return
642  (
643  ownToNbr,
644  ownToNbrDebugPtr,
645  pp,
646  faceMap,
647  rotation
648  );
649 }
650 
651 
653 {
655 
656  if (!nbrPatchName_.empty())
657  {
658  writeEntry(os, "neighbourPatch", nbrPatchName_);
659  }
660 
661  coupleGroup_.write(os);
662 
664 }
665 
666 
667 // ************************************************************************* //
#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:371
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:142
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
A list of faces which address into the list of points.
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).
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 nbrPatchID() 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:160
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:318
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const pointField & points
label nPoints
volScalarField & b
Definition: createFields.H:27
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:105
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
void deleteDemandDrivenData(DataPtr &dataPtr)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:483
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
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
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.
Unit conversion functions.