Loading...
Searching...
No Matches
faceAreaWeightAMI2D.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) 2020,2022,2024 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "faceAreaWeightAMI2D.H"
29#include "profiling.H"
31#include "triangle2D.H"
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
40 (
44 );
45}
46
47// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
48
50(
51 const label srcFacei,
52 const labelList& tgtFaceCandidates,
53 const boundBox& srcFaceBb
54) const
55{
56 Info<< "NO MATCH for source face " << srcFacei << endl;
57 {
58 const auto& src = this->srcPatch();
59 const auto& tgt = this->tgtPatch(); // might be the extended patch!
60
61 OFstream os("no_match_" + Foam::name(srcFacei) + ".obj");
62
63 const pointField& srcPoints = src.points();
64 const pointField& tgtPoints = tgt.points();
65
66 label np = 0;
67
68 // Write source face
69 const face& srcF = src[srcFacei];
70 string faceStr = "f";
71 for (const label pointi : srcF)
72 {
73 const point& p = srcPoints[pointi];
74 os << "v " << p.x() << " " << p.y() << " " << p.z() << nl;
75 ++np;
76 faceStr += " " + Foam::name(np);
77 }
78 os << faceStr.c_str() << " " << (np - srcF.size() + 1) << nl;
79
80 // Write target faces as lines
81 for (const label tgtFacei : tgtFaceCandidates)
82 {
83 const face& tgtF = tgt[tgtFacei];
84 forAll(tgtF, pointi)
85 {
86 const point& p = tgtPoints[tgtF[pointi]];
87 os << "v " << p.x() << " " << p.y() << " " << p.z() << nl;
88 ++np;
89 if (pointi)
90 {
91 os << "l " << np-1 << " " << np << nl;
92 }
93 }
94 os << "l " << (np - tgtF.size() + 1) << " " << np << nl;
95 }
96 }
97
98 {
99 OFstream os("no_match_" + Foam::name(srcFacei) + "_bb.obj");
100
101 const pointField points(srcFaceBb.points());
102 for (const point& p : points)
103 {
104 os << "v " << p.x() << " " << p.y() << " " << p.z() << endl;
105 }
106 os << "l 1 2" << nl;
107 os << "l 2 4" << nl;
108 os << "l 4 3" << nl;
109 os << "l 3 1" << nl;
110 os << "l 5 6" << nl;
111 os << "l 6 8" << nl;
112 os << "l 8 7" << nl;
113 os << "l 7 5" << nl;
114 os << "l 5 1" << nl;
115 os << "l 6 2" << nl;
116 os << "l 8 4" << nl;
117 os << "l 7 3" << nl;
118 }
119}
120
121
123(
124 const label srcFacei,
125 const label tgtFacei,
126 DynamicList<label>& srcAddr,
127 DynamicList<scalar>& srcWght,
128 DynamicList<vector>& srcCtr,
129 DynamicList<label>& tgtAddr,
130 DynamicList<scalar>& tgtWght
131) const
132{
133 addProfiling(ami, "faceAreaWeightAMI2D::calcInterArea");
134
135 // Quick reject if either face has zero area/too far away/wrong orientation
136 if (!isCandidate(srcFacei, tgtFacei))
137 {
138 return;
139 }
140
141 const auto& srcPatch = this->srcPatch();
142 const auto& tgtPatch = this->tgtPatch();
143
144 const pointField& srcPoints = srcPatch.points();
145 const pointField& tgtPoints = tgtPatch.points();
146
147 const auto& srcTris = srcTris_[srcFacei];
148 const auto& tgtTris = tgtTris_[tgtFacei];
149
150 const auto& srcFaceNormals = srcPatch.faceNormals();
151
152 //triangle2D::debug = 1;
153
154 scalar area = 0;
155 vector centroid = Zero;
156
157 for (const auto& tris : srcTris)
158 {
159 const vector& origin = srcPoints[tris[0]];
160 const vector p10(srcPoints[tris[1]] - origin);
161 const vector p20(srcPoints[tris[2]] - origin);
162 const vector axis1(p10/(mag(p10) + ROOTVSMALL));
163 const vector axis2(srcFaceNormals[srcFacei]^axis1);
164
166 (
167 vector2D(0, 0),
168 vector2D(axis1&p10, axis2&p10),
169 vector2D(axis1&p20, axis2&p20)
170 );
171
172 for (const auto& trit : tgtTris)
173 {
174 // Triangle t has opposite orientation wrt triangle s
175 triangle2D t
176 (
177 tgtPoints[trit[0]] - origin,
178 tgtPoints[trit[2]] - origin,
179 tgtPoints[trit[1]] - origin,
180 axis1,
181 axis2
182 );
183
184 scalar da = 0;
185 vector2D c(Zero);
186
187 if (t.snapClosePoints(s) == 3)
188 {
189 c = s.centre();
190 da = mag(s.area());
191 }
192 else
193 {
194 s.interArea(t, c, da);
195 }
196
197 area += da;
198 centroid += da*(origin + c.x()*axis1 + c.y()*axis2);
199 }
200 }
201
202 //triangle2D::debug = 0;
203
204 if (area/srcMagSf_[srcFacei] > triangle2D::relTol)
205 {
206 centroid /= area + ROOTVSMALL;
207
208 srcAddr.append(tgtFacei);
209 srcWght.append(area);
210 srcCtr.append(centroid);
212 tgtAddr.append(srcFacei);
213 tgtWght.append(area);
214 }
215}
216
217
219(
220 const AABBTree<face>& tree,
221 const List<boundBox>& tgtFaceBbs,
222 const boundBox& srcFaceBb
223) const
224{
225 labelHashSet faceIds(16);
226
227 const auto& treeBb = tree.boundBoxes();
228 const auto& treeAddr = tree.addressing();
229
230 forAll(treeBb, boxi)
231 {
232 const auto& tbb = treeBb[boxi];
233
234 if (srcFaceBb.overlaps(tbb))
235 {
236 const auto& boxAddr = treeAddr[boxi];
237
238 for (const auto& tgtFacei : boxAddr)
239 {
240 if (srcFaceBb.overlaps(tgtFaceBbs[tgtFacei]))
241 {
242 faceIds.insert(tgtFacei);
243 }
244 }
245 }
246 }
248 return faceIds.toc();
249}
250
251
252// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
253
255(
256 const dictionary& dict,
257 const bool reverseTarget
259:
261 Cbb_(dict.getCheckOrDefault<scalar>("Cbb", 0.1, scalarMinMax::ge(SMALL)))
262{}
263
264
266(
267 const bool requireMatch,
268 const bool reverseTarget,
269 const scalar lowWeightCorrection,
271 const bool restartUncoveredSourceFace
272)
273:
275 (
276 requireMatch,
277 reverseTarget,
279 triMode
280 ),
281 Cbb_(0.1)
282{}
283
284
286(
287 const faceAreaWeightAMI2D& ami
288)
289:
291 Cbb_(ami.Cbb_)
292{}
293
294
295// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
296
298(
299 const primitivePatch& srcPatch,
300 const primitivePatch& tgtPatch,
301 const autoPtr<searchableSurface>& surfPtr
302)
303{
304 if (upToDate_)
305 {
306 return false;
307 }
308
309 addProfiling(ami, "faceAreaWeightAMI2D::calculate");
310
311 advancingFrontAMI::calculate(srcPatch, tgtPatch, surfPtr);
312
313 const auto& src = this->srcPatch();
314 const auto& tgt = this->tgtPatch(); // might be the extended patch!
315
316 bool validSize = true;
317
318 // Check that patch sizes are valid
319 if (!src.size())
320 {
321 validSize = false;
322 }
323 else if (!tgt.size())
324 {
326 << src.size() << " source faces but no target faces" << endl;
327
328 validSize = false;
329 }
330
331 srcCentroids_.setSize(srcAddress_.size());
332
333 // Temporary storage for addressing and weights
334 List<DynamicList<label>> tgtAddr(tgt.size());
335 List<DynamicList<scalar>> tgtWght(tgt.size());
336
337 // Find an approximate face length scale
338 const scalar lf(Cbb_*Foam::sqrt(gAverage(srcMagSf_)));
339
340 if (validSize)
341 {
342 // Create the tgt tree
343 const bool equalBinSize = true;
344 const label maxLevel = 10;
345 const label minBinSize = 4;
347 (
348 tgt,
349 tgt.points(),
350 equalBinSize,
351 maxLevel,
352 minBinSize
353 );
354
355 const auto& tgtPoints = tgt.points();
356 List<boundBox> tgtFaceBbs(tgt.size());
357 forAll(tgt, facei)
358 {
359 tgtFaceBbs[facei] = boundBox(tgtPoints, tgt[facei], false);
360 }
361
362 DynamicList<label> nonOverlapFaces;
363
364 const auto& srcPoints = src.points();
365
366 forAll(src, srcFacei)
367 {
368 const face& srcFace = src[srcFacei];
369
370 treeBoundBox srcFaceBb(srcPoints, srcFace);
371 srcFaceBb.grow(lf);
372
373 const labelList tgtFaces
374 (
375 overlappingTgtFaces(tree, tgtFaceBbs, srcFaceBb)
376 );
377
378 DynamicList<label> srcAddr(tgtFaces.size());
379 DynamicList<scalar> srcWght(tgtFaces.size());
380 DynamicList<point> srcCtr(tgtFaces.size());
381
382 for (const label tgtFacei : tgtFaces)
383 {
384 storeInterArea
385 (
386 srcFacei,
387 tgtFacei,
388 srcAddr,
389 srcWght,
390 srcCtr,
391 tgtAddr[tgtFacei],
392 tgtWght[tgtFacei]
393 );
394 }
395
396 if (mustMatchFaces() && srcAddr.empty())
397 {
398 if (debug) writeNoMatch(srcFacei, tgtFaces, srcFaceBb);
399
400 // FatalErrorInFunction
401 // << "Unable to find match for source face " << srcFace
402 // << exit(FatalError);
403 }
404
405 srcAddress_[srcFacei].transfer(srcAddr);
406 srcWeights_[srcFacei].transfer(srcWght);
407 srcCentroids_[srcFacei].transfer(srcCtr);
408 }
409
410 srcNonOverlap_.transfer(nonOverlapFaces);
411
412 if (debug && !srcNonOverlap_.empty())
413 {
414 Pout<< " AMI: " << srcNonOverlap_.size()
415 << " non-overlap faces identified"
416 << endl;
417 }
418 }
419
420 // Transfer data to persistent storage
421
422 forAll(tgtAddr, i)
423 {
424 tgtAddress_[i].transfer(tgtAddr[i]);
425 tgtWeights_[i].transfer(tgtWght[i]);
426 }
427
428 if (distributed() && comm() != -1)
429 {
430 const label myRank = UPstream::myProcNo(comm());
431 // Allocate unique tag for all comms
432 const int oldTag = UPstream::incrMsgType();
433
434 const primitivePatch& srcPatch0 = this->srcPatch0();
435 const primitivePatch& tgtPatch0 = this->tgtPatch0();
436
437 // Create global indexing for each original patch
438 const globalIndex globalSrcFaces(srcPatch0.size(), comm());
439 const globalIndex globalTgtFaces(tgtPatch0.size(), comm());
440
441 for (labelList& addressing : srcAddress_)
442 {
443 for (label& addr : addressing)
444 {
445 addr = extendedTgtFaceIDs_[addr];
446 }
447 }
448
449 for (labelList& addressing : tgtAddress_)
450 {
451 globalSrcFaces.inplaceToGlobal(myRank, addressing);
452 }
453
454 // Send data back to originating procs. Note that contributions
455 // from different processors get added (ListOps::appendEqOp)
456
458 (
461 tgtPatch0.size(),
462 extendedTgtMapPtr_->constructMap(),
463 false, // has flip
464 extendedTgtMapPtr_->subMap(),
465 false, // has flip
466 tgtAddress_,
467 labelList(),
469 flipOp(), // flip operation
470 UPstream::msgType()+77431,
471 comm()
472 );
473
475 (
478 tgtPatch0.size(),
479 extendedTgtMapPtr_->constructMap(),
480 false,
481 extendedTgtMapPtr_->subMap(),
482 false,
483 tgtWeights_,
484 scalarList(),
486 flipOp(), // flip operation
487 UPstream::msgType()+77432,
488 comm()
489 );
490
491 // Note: using patch face areas calculated by the AMI method
492 extendedTgtMapPtr_->reverseDistribute(tgtPatch0.size(), tgtMagSf_);
493
494 // Cache maps and reset addresses
495 List<Map<label>> cMapSrc;
496 srcMapPtr_.reset
497 (
498 new mapDistribute
499 (
500 globalSrcFaces,
501 tgtAddress_,
502 cMapSrc,
503 UPstream::msgType()+77433,
504 comm()
505 )
506 );
507
508 List<Map<label>> cMapTgt;
509 tgtMapPtr_.reset
510 (
511 new mapDistribute
512 (
513 globalTgtFaces,
514 srcAddress_,
515 cMapTgt,
516 UPstream::msgType()+77434,
517 comm()
518 )
519 );
520
521 // Reset tag
522 UPstream::msgType(oldTag);
523 }
524
525 // Convert the weights from areas to normalised values
526 normaliseWeights(requireMatch_, true);
527
528 nonConformalCorrection();
530 upToDate_ = true;
531
532 return upToDate_;
533}
534
535
537{
539 os.writeEntryIfDifferent<scalar>("Cbb", 0.1, Cbb_);
540}
541
542
543// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Templated tree of axis-aligned bounding boxes (AABB).
Definition AABBTree.H:116
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
bool requireMatch_
Flag to indicate that the two patches must be matched/an overlap exists between them.
bool mustMatchFaces() const
Return true if requireMatch and but not lowWeightCorrection.
bool reverseTarget() const noexcept
Access to the reverseTarget flag.
labelListList srcAddress_
Addresses of target faces per source face.
label comm() const noexcept
Communicator (local or otherwise) for parallel operations.
bool distributed() const noexcept
Distributed across processors (singlePatchProc == -1).
scalarList tgtMagSf_
Target face areas.
autoPtr< mapDistribute > srcMapPtr_
Source map pointer - parallel running only.
bool requireMatch() const noexcept
Return the requireMatch flag.
const primitivePatch & tgtPatch0() const
Return the original tgt patch with optionally updated points.
bool upToDate_
Up-to-date flag.
autoPtr< mapDistribute > tgtMapPtr_
Target map pointer - parallel running only.
labelListList tgtAddress_
Addresses of source faces per target face.
pointListList srcCentroids_
Centroid of target faces per source face.
static void normaliseWeights(const scalarList &patchAreas, const word &patchName, const labelListList &addr, scalarListList &wght, scalarField &wghtSum, const bool conformal, const bool output, const scalar lowWeightTol, const label comm)
Normalise the (area) weights - suppresses numerical error in weights calculation.
scalar lowWeightCorrection() const
Threshold weight below which interpolation is deactivated.
scalarListList tgtWeights_
Weights of source faces per target face.
scalarListList srcWeights_
Weights of target faces per source face.
scalarList srcMagSf_
Source face areas.
const primitivePatch & srcPatch0() const
Return the original src patch with optionally updated points.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition HashTable.C:141
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
static const List< T > & null() noexcept
Return a null List (reference to a nullObject). Behaves like an empty List.
Definition List.H:138
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Definition UPstream.H:84
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition UPstream.H:1948
Base class for Arbitrary Mesh Interface (AMI) methods.
const primitivePatch & tgtPatch() const
Return const access to the target patch.
List< DynamicList< face > > srcTris_
Storage for src-side triangle decomposition.
virtual bool calculate(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const autoPtr< searchableSurface > &surfPtr=nullptr)
Update addressing, weights and (optional) centroids.
bool isCandidate(const label srcFacei, const label tgtFacei) const
Is source/target a valid pair (i.e. not too far/different.
virtual void nonConformalCorrection()
Correction for non-conformal interpolations, e.g. for ACMI.
virtual void write(Ostream &os) const
Write.
advancingFrontAMI(const dictionary &dict, const bool reverseTarget)
Construct from components.
List< DynamicList< face > > tgtTris_
Storage for tgt-side triangle decomposition.
const primitivePatch & srcPatch() const
Return const access to the source patch.
autoPtr< mapDistribute > extendedTgtMapPtr_
Extended patch map.
labelList extendedTgtFaceIDs_
Extended patch face IDs.
labelList srcNonOverlap_
Labels of faces that are not overlapped by any target faces (should be empty for correct functioning ...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
void grow(const scalar delta)
Expand box by adjusting min/max by specified amount in each dimension.
Definition boundBoxI.H:367
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition boundBoxI.H:439
tmp< pointField > points() const
Corner points in an order corresponding to a 'hex' cell.
Definition boundBox.H:381
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Face area weighted Arbitrary Mesh Interface (AMI) method that performs the intersection of src and tg...
void storeInterArea(const label srcFacei, const label tgtFacei, DynamicList< label > &srcAddr, DynamicList< scalar > &srcWght, DynamicList< vector > &srcCtr, DynamicList< label > &tgtAddr, DynamicList< scalar > &tgtWght) const
Calculate and store the area of intersection between source and target faces.
virtual bool calculate(const primitivePatch &srcPatch, const primitivePatch &tgtPatch, const autoPtr< searchableSurface > &surfPtr=nullptr)
Update addressing, weights and (optional) centroids.
scalar Cbb_
Face bounding box factor.
virtual void write(Ostream &os) const
Write.
labelList overlappingTgtFaces(const AABBTree< face > &tree, const List< boundBox > &tgtFaceBbs, const boundBox &srcFaceBb) const
Return the set of tgt face IDs that overlap the src face bb.
void writeNoMatch(const label srcFacei, const labelList &tgtFaceCandidates, const boundBox &srcFaceBb) const
Helper function to write non-matched source faces to the set of candidate faces.
faceAreaWeightAMI2D(const dictionary &dict, const bool reverseTarget=false)
Construct from dictionary.
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
void inplaceToGlobal(const label proci, labelUList &labels) const
From local to global index on proci (inplace).
static void distribute(const UPstream::commsTypes commsType, const UList< labelPair > &schedule, const label constructSize, const labelListList &subMap, const bool subHasFlip, const labelListList &constructMap, const bool constructHasFlip, List< T > &field, const T &nullValue, const CombineOp &cop, const NegateOp &negOp, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Distribute combine data with specified combine operation and negate operator (for flips).
Class containing processor-to-processor mapping information.
Standard boundBox with extra functionality for use in octree.
2-D triangle and queries
Definition triangle2D.H:51
label snapClosePoints(const triangle2D &triB)
Snap [this] triangle's points to those of triB if they are within absTol.
Definition triangle2D.C:86
static scalar relTol
Relative tolerance.
Definition triangle2D.H:72
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for handling debugging switches.
Definition debug.C:45
const wordList area
Standard area field types (scalar, vector, tensor, etc).
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
List< label > labelList
A List of labels.
Definition List.H:62
Vector2D< scalar > vector2D
A 2D vector of scalars obtained from the generic Vector2D.
Definition vector2D.H:56
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
messageStream Info
Information stream (stdout output on master, null elsewhere).
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...
dictionary dict
Tree tree(triangles.begin(), triangles.end())
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
List helper to append y elements onto the end of x.
Definition ListOps.H:712
Functor to negate primitives. Dummy for most other types.
Definition flipOp.H:67