Loading...
Searching...
No Matches
processorFaPatch.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) 2016-2017 Wikki Ltd
9 Copyright (C) 2019-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "processorFaPatch.H"
30#include "processorPolyPatch.H" // For newName()
32#include "transformField.H"
33#include "faMesh.H"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38namespace Foam
44// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
45
47{
48 neighbPointsPtr_.clear();
49 nonGlobalPatchPointsPtr_.clear();
50}
51
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
56(
57 const word& name,
58 const labelUList& edgeLabels,
59 const label index,
60 const faBoundaryMesh& bm,
61 const label nbrPolyPatchi,
62 const label myProcNo,
63 const label neighbProcNo,
64 const word& patchType
65)
66:
67 coupledFaPatch(name, edgeLabels, index, bm, nbrPolyPatchi, patchType),
68 myProcNo_(myProcNo),
69 neighbProcNo_(neighbProcNo),
70 neighbEdgeCentres_(),
71 neighbEdgeLengths_(),
72 neighbEdgeFaceCentres_(),
73 neighbPointsPtr_(nullptr),
74 nonGlobalPatchPointsPtr_(nullptr)
75{}
76
77
79(
80 const labelUList& edgeLabels,
81 const label index,
82 const faBoundaryMesh& bm,
83 const label nbrPolyPatchi,
84 const label myProcNo,
85 const label neighbProcNo,
86 const word& patchType
87)
88:
90 (
91 processorPolyPatch::newName(myProcNo, neighbProcNo),
92 edgeLabels,
93 index,
94 bm,
95 nbrPolyPatchi,
98 patchType
99 )
100{}
101
102
104(
105 const word& name,
106 const dictionary& dict,
107 const label index,
108 const faBoundaryMesh& bm,
109 const word& patchType
110)
111:
112 coupledFaPatch(name, dict, index, bm, patchType),
113 myProcNo_(dict.get<label>("myProcNo")),
114 neighbProcNo_(dict.get<label>("neighbProcNo")),
115 neighbEdgeCentres_(),
116 neighbEdgeLengths_(),
117 neighbEdgeFaceCentres_(),
118 neighbPointsPtr_(nullptr),
119 nonGlobalPatchPointsPtr_(nullptr)
120{}
121
122
123// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
125Foam::label Foam::processorFaPatch::comm() const
126{
127 return boundaryMesh().mesh().comm();
128}
129
130
132{
133 // If it is not running parallel or there are no global points
134 // create a 1->1 map
135
136 // Can not use faGlobalMeshData at this point yet
137 // - use polyMesh globalData instead
138
139 const auto& pMeshGlobalData = boundaryMesh().mesh().mesh().globalData();
140
141 if (!Pstream::parRun() || !pMeshGlobalData.nGlobalPoints())
142 {
143 // 1 -> 1 mapping
144 nonGlobalPatchPointsPtr_.reset(new labelList(identity(nPoints())));
145 }
146 else
147 {
148 nonGlobalPatchPointsPtr_.reset(new labelList(nPoints()));
149 labelList& ngpp = *nonGlobalPatchPointsPtr_;
150
151 // Get reference to shared points
152 const labelList& sharedPoints = pMeshGlobalData.sharedPointLabels();
153
154 const labelList& faMeshPatchPoints = pointLabels();
155
156 const labelList& meshPoints =
157 boundaryMesh().mesh().patch().meshPoints();
158
159 label nNonShared = 0;
160
161 forAll(faMeshPatchPoints, pointi)
162 {
163 const label mpi = meshPoints[faMeshPatchPoints[pointi]];
164 if (!sharedPoints.found(mpi))
165 {
166 ngpp[nNonShared] = pointi;
167 ++nNonShared;
168 }
170
171 ngpp.resize(nNonShared);
172 }
173}
174
175
177{
178 if (Pstream::parRun())
179 {
180 if (neighbProcNo() >= pBufs.nProcs())
181 {
183 << "On patch " << name()
184 << " trying to access out of range neighbour processor "
185 << neighbProcNo() << ". This can happen if" << nl
186 << " trying to run on an incorrect number of processors"
187 << exit(FatalError);
188 }
189
190 UOPstream toNeighbProc(neighbProcNo(), pBufs);
191
192 toNeighbProc
194 << edgeLengths()
195 << edgeFaceCentres();
196 }
197}
198
199
200void Foam::processorFaPatch::calcGeometry(PstreamBuffers& pBufs)
201{
202 if (Pstream::parRun())
203 {
204 {
205 UIPstream fromNeighbProc(neighbProcNo(), pBufs);
206
207 fromNeighbProc
208 >> neighbEdgeCentres_
209 >> neighbEdgeLengths_
210 >> neighbEdgeFaceCentres_;
211 }
212
213 #ifdef FULLDEBUG
214 const scalarField& magEl = magEdgeLengths();
215
216 forAll(magEl, edgei)
217 {
218 scalar nmagEl = mag(neighbEdgeLengths_[edgei]);
219 scalar avEl = (magEl[edgei] + nmagEl)/2.0;
220
221 if (mag(magEl[edgei] - nmagEl)/avEl > 1e-6)
222 {
224 << "edge " << edgei
225 << " length does not match neighbour by "
226 << 100*mag(magEl[edgei] - nmagEl)/avEl
227 << "% -- possible edge ordering problem" << nl;
228 }
229 }
230 #endif
231
232 calcTransformTensors
233 (
234 edgeCentres(),
235 neighbEdgeCentres_,
238 );
239 }
240}
241
242
244(
245 PstreamBuffers& pBufs,
246 const pointField& p
248{
249 faPatch::movePoints(pBufs, p);
250 initGeometry(pBufs);
251}
252
253
255(
256 PstreamBuffers& pBufs,
258)
259{
261}
262
263
265{
266 // For completeness
268
269 neighbPointsPtr_.clear();
270
271 if (Pstream::parRun())
272 {
273 if (neighbProcNo() >= pBufs.nProcs())
274 {
276 << "On patch " << name()
277 << " trying to access out of range neighbour processor "
278 << neighbProcNo() << ". This can happen if" << nl
279 << " trying to run on an incorrect number of processors"
280 << exit(FatalError);
281 }
282
283 // Express all points as patch edge and index in edge.
284 labelList patchEdge(nPoints());
285 labelList indexInEdge(nPoints());
286
287 const edgeList::subList patchEdges =
288 patchSlice(boundaryMesh().mesh().edges());
289
290 const labelListList& ptEdges = pointEdges();
291
292 for (label patchPointI = 0; patchPointI < nPoints(); ++patchPointI)
293 {
294 label edgeI = ptEdges[patchPointI][0];
295
296 patchEdge[patchPointI] = edgeI;
297
298 const edge& e = patchEdges[edgeI];
299
300 indexInEdge[patchPointI] = e.find(pointLabels()[patchPointI]);
301 }
302
303 UOPstream toNeighbProc(neighbProcNo(), pBufs);
304
305 toNeighbProc
306 << patchEdge
307 << indexInEdge;
308 }
309}
310
311
312void Foam::processorFaPatch::updateMesh(PstreamBuffers& pBufs)
313{
314 // For completeness
315 faPatch::updateMesh(pBufs);
316
317 neighbPointsPtr_.clear();
318
319 if (Pstream::parRun())
320 {
321 labelList nbrPatchEdge;
322 labelList nbrIndexInEdge;
323
324 {
325 // Note cannot predict exact size since edgeList not (yet) sent as
326 // binary entity but as List of edges.
327
328 UIPstream fromNeighbProc(neighbProcNo(), pBufs);
329
330 fromNeighbProc
331 >> nbrPatchEdge
332 >> nbrIndexInEdge;
333 }
334
335 if (nbrPatchEdge.size() == nPoints())
336 {
337 // Convert neighbour edges and indices into face back into
338 // my edges and points.
339 neighbPointsPtr_.reset(new labelList(nPoints()));
340 labelList& neighbPoints = *neighbPointsPtr_;
341
342 const edgeList::subList patchEdges =
343 patchSlice(boundaryMesh().mesh().edges());
344
345 forAll(nbrPatchEdge, nbrPointI)
346 {
347 // Find edge and index in edge on this side.
348 const edge& e = patchEdges[nbrPatchEdge[nbrPointI]];
349
350 const label index = 1 - nbrIndexInEdge[nbrPointI];
351
352 const label patchPointI = pointLabels().find(e[index]);
353
354 neighbPoints[patchPointI] = nbrPointI;
355 }
356 }
357 else
358 {
359 // Differing number of points. Probably patch includes
360 // part of a cyclic.
361 }
362 }
363}
364
365
366Foam::tmp<Foam::vectorField> Foam::processorFaPatch::neighbEdgeNormals() const
368 auto tresult = tmp<vectorField>::New(neighbEdgeLengths_);
369 tresult.ref().normalise();
370 return tresult;
371}
372
373
375{
376 if (!neighbPointsPtr_)
377 {
378 // Was probably created from cyclic patch and hence the
379 // number of edges or points might differ on both
380 // sides of the processor patch since one side might have
381 // it merged with another bit of geometry
382
384 << "No extended addressing calculated for patch " << name()
385 << nl
386 << "This can happen if the number of points on both"
387 << " sides of the two coupled patches differ." << nl
388 << "This happens if the processorPatch was constructed from"
389 << " part of a cyclic patch."
391 }
392
393 return *neighbPointsPtr_;
394}
395
396
398{
399 if (Pstream::parRun())
400 {
401 const vectorField& skew = skewCorrectionVectors();
402 const scalarField& PN = lPN();
403
404 const scalarField lEN
405 (
406 mag(neighbEdgeFaceCentres() - edgeCentres() + skew)
407 );
408
409 forAll(w, i)
410 {
411 if (mag(PN[i]) > SMALL)
412 {
413 w[i] = lEN[i]/PN[i];
414 }
415 }
416 }
417 else
418 {
419 w = scalar(1);
420 }
421}
422
423
425{
426 const vectorField& skew = skewCorrectionVectors();
427
428 if (Pstream::parRun())
429 {
430 lPN =
431 mag(edgeCentres() - skew - edgeFaceCentres()) // lPE
432 + mag(neighbEdgeFaceCentres() - edgeCentres() + skew); // lEN
433
434 for (scalar& lpn: lPN)
435 {
436 if (mag(lpn) < SMALL)
437 {
438 lpn = SMALL;
439 }
440 }
441 }
442 else
443 {
445 mag(edgeCentres() - skew - edgeFaceCentres()) // lPE
446 + mag(edgeFaceCentres() - edgeCentres() + skew); // lEN
447 }
448}
449
450
452{
453 if (Pstream::parRun())
454 {
455 const edgeList& edges = boundaryMesh().mesh().edges();
456 const pointField& points = boundaryMesh().mesh().points();
457 const vectorField& lengths = edgeLengths();
458 const scalarField& PN = lPN();
459
460 tmp<vectorField> tdelta = processorFaPatch::delta();
461 vectorField& unitDelta = tdelta.ref();
462
463 forAll(dc, i)
464 {
465 vector edgeNormal =
466 normalised(lengths[i] ^ edges[i + start()].vec(points));
467
468 unitDelta[i].removeCollinear(edgeNormal);
469 unitDelta[i].normalise();
470
471 edgeNormal = normalised(lengths[i]);
472
473 const scalar alpha = PN[i]*(edgeNormal & unitDelta[i]);
474 if (mag(alpha) > SMALL)
475 {
476 dc[i] = scalar(1)/max(alpha, 0.05*PN[i]);
477 }
478 }
479 }
480 else
482 dc = scalar(1)
483 /((edgeNormals() & coupledFaPatch::delta()) + VSMALL);
484 }
485}
486
487
489{
490 if (Pstream::parRun())
491 {
492 const edgeList& edges = boundaryMesh().mesh().edges();
493 const pointField& points = boundaryMesh().mesh().points();
494 const vectorField& lengths = edgeLengths();
495
496 tmp<vectorField> tdelta = processorFaPatch::delta();
497 vectorField& unitDelta = tdelta.ref();
498
499 forAll(cv, i)
500 {
501 vector edgeNormal =
502 normalised(lengths[i] ^ edges[i + start()].vec(points));
503
504 unitDelta[i].removeCollinear(edgeNormal);
505 unitDelta.normalise();
506
507 edgeNormal = normalised(lengths[i]);
508
509 const scalar alpha = unitDelta[i] & edgeNormal;
510 scalar dc = SMALL;
511 if (mag(alpha) > SMALL)
512 {
513 dc = scalar(1)/alpha;
514 }
515
516 cv[i] = edgeNormal - dc*unitDelta[i];
517 }
518 }
519 else
520 {
521 cv = Zero;
522 }
523}
524
525
527{
528 if (Pstream::parRun())
529 {
530 // Do the transformation if necessary
531 if (parallel())
532 {
533 return
535 - (
536 neighbEdgeCentres()
537 - neighbEdgeFaceCentres()
538 );
539 }
540 else
541 {
542 return
544 - transform
545 (
546 forwardT(),
547 (
548 neighbEdgeCentres()
549 - neighbEdgeFaceCentres()
550 )
551 );
552 }
553 }
554 else
555 {
556 return coupledFaPatch::delta();
557 }
558}
559
560
562{
563 if (!nonGlobalPatchPointsPtr_)
564 {
566 }
567
568 return *nonGlobalPatchPointsPtr_;
569}
570
571
573(
574 const labelUList& internalData
575) const
576{
577 return patchInternalField(internalData);
578}
579
580
582(
583 const labelUList& internalData,
584 const labelUList& edgeFaces
585) const
587 auto tpfld = tmp<labelField>::New(this->size());
588 patchInternalField(internalData, edgeFaces, tpfld.ref());
589 return tpfld;
590}
591
592
594(
595 const Pstream::commsTypes commsType,
596 const labelUList& interfaceData
597) const
598{
599 send(commsType, interfaceData);
600}
601
602
604(
605 const Pstream::commsTypes commsType,
607) const
608{
609 return receive<label>(commsType, this->size());
610}
611
612
614(
615 const Pstream::commsTypes commsType,
624(
625 const Pstream::commsTypes commsType,
627) const
628{
629 return receive<label>(commsType, this->size());
630}
631
632
634(
635 const Pstream::commsTypes commsType,
636 const labelUList&,
638) const
639{
640 return receive<label>(commsType, this->size());
641}
642
643
645{
647 os.writeEntry("myProcNo", myProcNo_);
648 os.writeEntry("neighbProcNo", neighbProcNo_);
649}
650
651
652// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Info<< " "<< writer.output().name()<< nl;}{ const Field< vector > edgeCentres(faMeshTools::flattenEdgeField(aMesh.edgeCentres(), true))
void normalise()
Inplace normalise this field. Generally a no-op except for vector fields.
Definition Field.C:600
SubList< edge > subList
Definition List.H:129
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & points() const noexcept
Return reference to global points.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
int nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
Input inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UIPstream.H:313
bool found(const T &val, label pos=0) const
Same as contains().
Definition UList.H:983
bool get(const label i) const
Definition UList.H:868
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label find(const T &val) const
Find index of the first occurrence of the value.
Definition UList.C:160
Output inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UOPstream.H:408
commsTypes
Communications types.
Definition UPstream.H:81
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
const bMesh & mesh() const
coupledFaPatch is an abstract base class for patches that couple regions of the computational domain ...
bool parallel() const
Are the cyclic planes parallel.
coupledFaPatch(const word &name, const labelUList &edgeLabels, const label index, const faBoundaryMesh &bm, const label nbrPolyPatchIndex, const word &patchType)
Construct from components.
void calcTransformTensors(const vector &Cf, const vector &Cr, const vector &nf, const vector &nr) const
Calculate the uniform transformation tensors.
virtual tmp< vectorField > delta() const =0
Return delta (P to N) vectors across coupled patch.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
Finite area boundary mesh, which is a faPatch list with registered IO, a reference to the associated ...
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition faPatch.H:76
virtual label size() const
Patch size is the number of edge labels, but can be overloaded.
Definition faPatch.H:392
label nPoints() const
Number of patch points.
Definition faPatch.H:343
const labelListList & pointEdges() const
Return patch point-edge addressing.
Definition faPatch.C:278
virtual void write(Ostream &) const
Write.
Definition faPatch.C:565
const labelList & edgeLabels() const noexcept
Return the list of edges.
Definition faPatch.H:335
const scalarField & magEdgeLengths() const
Return edge length magnitudes, like the faMesh::magLe() method.
Definition faPatch.C:450
const scalarField & lPN() const
Return patch geodesic distance between P and N.
Definition faPatch.C:502
const vectorField & edgeLengths() const
Return edge length vectors, like the faMesh::Le() method.
Definition faPatch.C:444
tmp< vectorField > edgeFaceCentres() const
Return neighbour face centres.
Definition faPatch.C:466
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition faPatch.H:168
const labelUList & edgeFaces() const
Return edge-face addressing.
Definition faPatch.C:424
virtual void movePoints(PstreamBuffers &, const pointField &)
Correct patch after moving points.
Definition faPatch.C:547
const List< T >::subList patchSlice(const UList< T > &values) const
This patch slice from the complete list of values, which has size mesh::nEdges(), using the virtual p...
Definition faPatch.H:411
tmp< vectorField > edgeNormals() const
Return edge unit normals, like the faMesh::unitLe() method.
Definition faPatch.C:456
void patchInternalField(const UList< Type > &internalData, const labelUList &addressing, UList< Type > &pfld) const
Extract internal field next to patch using specified addressing.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition faPatch.H:174
label start() const
Patch start in edge list.
Definition faPatch.C:191
label index() const noexcept
The index of this patch in the boundaryMesh.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
virtual void initTransfer(const Pstream::commsTypes commsType, const labelUList &interfaceData) const
Initialise interface data transfer.
virtual ~processorFaPatch()
Destructor.
processorFaPatch(const word &name, const labelUList &edgeLabels, const label index, const faBoundaryMesh &bm, const label nbrPolyPatchi, const label myProcNo, const label neighbProcNo, const word &patchType=typeName)
Construct from components with specified name.
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Return neighbour field.
void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
void makeWeights(scalarField &) const
Make patch weighting factors.
void makeLPN(scalarField &) const
Make patch geodesic distance between P and N.
void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
const labelList & neighbPoints() const
Return neighbour point labels. This is for my local point the.
void makeNonGlobalPatchPoints() const
Find non-globa patch points.
void makeDeltaCoeffs(scalarField &) const
Make patch face - neighbour cell distances.
int neighbProcNo() const noexcept
Return neighbour processor number.
virtual tmp< labelField > transfer(const Pstream::commsTypes commsType, const labelUList &interfaceData) const
Transfer and return neighbour field.
void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
virtual void write(Ostream &os) const
Write the patch data as a dictionary.
tmp< vectorField > neighbEdgeNormals() const
Return processor-neighbour patch edge unit normals.
virtual label comm() const
Return communicator used for communication.
virtual tmp< vectorField > delta() const
Return delta (P to N) vectors across coupled patch.
const vectorField & neighbEdgeCentres() const noexcept
Return processor-neighbour patch edge centres.
virtual void initInternalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &internalData) const
Initialise neighbour field transfer.
int myProcNo() const noexcept
Return processor number.
virtual tmp< labelField > interfaceInternalField(const labelUList &internalData) const
Return the values of the given internal data adjacent to the interface as a field.
void movePoints(PstreamBuffers &, const pointField &)
Correct patches after moving points.
const vectorField & neighbEdgeFaceCentres() const noexcept
Return processor-neighbour patch neighbour face centres.
const labelList & nonGlobalPatchPoints() const
Return the set of labels of the processor patch points which are.
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
void makeCorrectionVectors(vectorField &) const
Make non-orthogonality correction vectors.
virtual const tensorField & forwardT() const
Return face transformation tensor.
void receive(const UPstream::commsTypes commsType, UList< Type > &f) const
Raw receive function.
void send(const UPstream::commsTypes commsType, const UList< Type > &f) const
Raw send function.
Neighbour processor patch.
Skew-correction vectors for the skewness-corrected interpolation scheme.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
auto & name
const pointField & points
label nPoints
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
dimensionedTensor skew(const dimensionedTensor &dt)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList pointLabels(nPoints, -1)
volScalarField & alpha
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Spatial transformation functions for primitive fields.