edgeIntersections.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2015-2020 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 "edgeIntersections.H"
30#include "triSurfaceSearch.H"
31#include "OFstream.H"
32#include "triSurface.H"
33#include "pointIndexHit.H"
34#include "treeDataTriSurface.H"
35#include "indexedOctree.H"
36#include "meshTools.H"
37#include "plane.H"
38#include "Random.H"
39#include "unitConversion.H"
40#include "treeBoundBox.H"
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
44namespace Foam
45{
47}
48
49Foam::scalar Foam::edgeIntersections::alignedCos_ = Foam::cos(degToRad(89.0));
50
51
52// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53
54void Foam::edgeIntersections::checkEdges(const triSurface& surf)
55{
56 const pointField& localPoints = surf.localPoints();
57 const edgeList& edges = surf.edges();
58
59 treeBoundBox bb(localPoints);
60
61 scalar minSize = SMALL * bb.minDim();
62
63 forAll(edges, edgeI)
64 {
65 const edge& e = edges[edgeI];
66
67 scalar eMag = e.mag(localPoints);
68
69 if (eMag < minSize)
70 {
72 << "Edge " << edgeI << " vertices " << e
73 << " coords:" << localPoints[e[0]] << ' '
74 << localPoints[e[1]] << " is very small compared to bounding"
75 << " box dimensions " << bb << endl
76 << "This might lead to problems in intersection"
77 << endl;
78 }
79 }
80}
81
82
83// Update intersections for selected edges.
84void Foam::edgeIntersections::intersectEdges
85(
86 const triSurface& surf1,
87 const pointField& points1, // surf1 meshPoints (not localPoints!)
88 const triSurfaceSearch& querySurf2,
89 const scalarField& surf1PointTol, // surf1 tolerance per point
90 const labelList& edgeLabels
91)
92{
93 const triSurface& surf2 = querySurf2.surface();
94 const vectorField& normals2 = surf2.faceNormals();
95
96 const labelList& meshPoints = surf1.meshPoints();
97
98 if (debug)
99 {
100 Pout<< "Calculating intersection of " << edgeLabels.size() << " edges"
101 << " out of " << surf1.nEdges() << " with " << surf2.size()
102 << " triangles ..." << endl;
103 }
104
105 pointField start(edgeLabels.size());
106 pointField end(edgeLabels.size());
107 vectorField edgeDirs(edgeLabels.size());
108
109 // Go through all edges, calculate intersections
110 forAll(edgeLabels, i)
111 {
112 label edgeI = edgeLabels[i];
113
114 if (debug)// && (i % 1000 == 0))
115 {
116 Pout<< "Intersecting edge " << edgeI << " with surface" << endl;
117 }
118
119 const edge& e = surf1.edges()[edgeI];
120
121 const point& pStart = points1[meshPoints[e.start()]];
122 const point& pEnd = points1[meshPoints[e.end()]];
123
124 const vector n = normalised(pEnd - pStart);
125
126 // Start tracking somewhat before pStart and up to somewhat after p1.
127 // Note that tolerances here are smaller than those used to classify
128 // hit below.
129 // This will cause this hit to be marked as degenerate and resolved
130 // later on.
131 start[i] = pStart - 0.5*surf1PointTol[e[0]]*n;
132 end[i] = pEnd + 0.5*surf1PointTol[e[1]]*n;
133
134 edgeDirs[i] = n;
135 }
136
137 List<List<pointIndexHit>> edgeIntersections;
138 querySurf2.findLineAll
139 (
140 start,
141 end,
142 edgeIntersections
143 );
144
145 label nHits = 0;
146
147 // Classify the hits
148 forAll(edgeLabels, i)
149 {
150 const label edgeI = edgeLabels[i];
151
152 labelList& intersectionTypes = classification_[edgeI];
153 intersectionTypes.setSize(edgeIntersections[i].size(), -1);
154
155 this->operator[](edgeI).transfer(edgeIntersections[i]);
156
157 forAll(intersectionTypes, hitI)
158 {
159 const pointIndexHit& pHit = this->operator[](edgeI)[hitI];
160 label& hitType = intersectionTypes[hitI];
161
162 if (!pHit.hit())
163 {
164 continue;
165 }
166
167 const edge& e = surf1.edges()[edgeI];
168
169 // Classify point on surface1 edge.
170 if (mag(pHit.hitPoint() - start[i]) < surf1PointTol[e[0]])
171 {
172 // Intersection is close to edge start
173 hitType = 0;
174 }
175 else if (mag(pHit.hitPoint() - end[i]) < surf1PointTol[e[1]])
176 {
177 // Intersection is close to edge end
178 hitType = 1;
179 }
180 else if (mag(edgeDirs[i] & normals2[pHit.index()]) < alignedCos_)
181 {
182 // Edge is almost coplanar with the face
183
184 Pout<< "Flat angle edge:" << edgeI
185 << " face:" << pHit.index()
186 << " cos:" << mag(edgeDirs[i] & normals2[pHit.index()])
187 << endl;
188
189 hitType = 2;
190 }
191
192 if (debug)
193 {
194 Info<< " hit " << pHit << " classify = " << hitType << endl;
195 }
196
197 nHits++;
198 }
199 }
200
201 if (debug)
202 {
203 Pout<< "Found " << nHits << " intersections of edges with surface ..."
204 << endl;
205 }
206}
207
208
209// If edgeI intersections are close to endpoint of edge shift endpoints
210// slightly along edge
211// Updates
212// - points1 with new endpoint position
213// - affectedEdges with all edges affected by moving the point
214// Returns true if changed anything.
215bool Foam::edgeIntersections::inlinePerturb
216(
217 const triSurface& surf1,
218 const scalarField& surf1PointTol, // surf1 tolerance per point
219 const label edgeI,
220 Random& rndGen,
221 pointField& points1,
222 boolList& affectedEdges
223) const
224{
225 bool hasPerturbed = false;
226
227 // Check if edge close to endpoint. Note that we only have to check
228 // the intersection closest to the edge endpoints (i.e. first and last in
229 // edgeEdgs)
230
231 const labelList& edgeEnds = classification_[edgeI];
232
233 if (edgeEnds.size())
234 {
235 bool perturbStart = false;
236 bool perturbEnd = false;
237
238 // Check first intersection.
239 if (edgeEnds.first() == 0)
240 {
241 perturbStart = true;
242 }
243
244 if (edgeEnds.last() == 1)
245 {
246 perturbEnd = true;
247 }
248
249
250 if (perturbStart || perturbEnd)
251 {
252 const edge& e = surf1.edges()[edgeI];
253
254 label v0 = surf1.meshPoints()[e[0]];
255 label v1 = surf1.meshPoints()[e[1]];
256
257 const vector n = normalised(points1[v1] - points1[v0]);
258
259 if (perturbStart)
260 {
261 // Perturb with something (hopefully) larger than tolerance.
262 scalar t = 4.0*(rndGen.sample01<scalar>() - 0.5);
263 points1[v0] += t*surf1PointTol[e[0]]*n;
264
265 const labelList& pEdges = surf1.pointEdges()[e[0]];
266
267 forAll(pEdges, i)
268 {
269 affectedEdges[pEdges[i]] = true;
270 }
271 }
272 if (perturbEnd)
273 {
274 // Perturb with something larger than tolerance.
275 scalar t = 4.0*(rndGen.sample01<scalar>() - 0.5);
276 points1[v1] += t*surf1PointTol[e[1]]*n;
277
278 const labelList& pEdges = surf1.pointEdges()[e[1]];
279
280 forAll(pEdges, i)
281 {
282 affectedEdges[pEdges[i]] = true;
283 }
284 }
285 hasPerturbed = true;
286 }
287 }
288
289 return hasPerturbed;
290}
291
292
293// Perturb single edge endpoint when perpendicular to face
294bool Foam::edgeIntersections::rotatePerturb
295(
296 const triSurface& surf1,
297 const scalarField& surf1PointTol, // surf1 tolerance per point
298 const label edgeI,
299
300 Random& rndGen,
301 pointField& points1,
302 boolList& affectedEdges
303) const
304{
305 const labelList& meshPoints = surf1.meshPoints();
306
307 const labelList& edgeEnds = classification_[edgeI];
308
309 bool hasPerturbed = false;
310
311 forAll(edgeEnds, i)
312 {
313 if (edgeEnds[i] == 2)
314 {
315 const edge& e = surf1.edges()[edgeI];
316
317 // Endpoint to modify. Choose either start or end.
318 label pointi = e[rndGen.bit()];
319 //label pointi = e[0];
320
321 // Generate random vector slightly larger than tolerance.
322 vector rndVec = rndGen.sample01<vector>() - vector(0.5, 0.5, 0.5);
323
324 // Make sure rndVec only perp to edge
325 vector n(points1[meshPoints[e[1]]] - points1[meshPoints[e[0]]]);
326 scalar magN = mag(n) + VSMALL;
327 n /= magN;
328
329 rndVec.removeCollinear(n);
330 rndVec.normalise();
331
332 // Scale to be moved by tolerance.
333 rndVec *= 0.01*magN;
334
335 Pout<< "rotating: shifting endpoint " << meshPoints[pointi]
336 << " of edge:" << edgeI << " verts:"
337 << points1[meshPoints[e[0]]] << ' '
338 << points1[meshPoints[e[1]]]
339 << " by " << rndVec
340 << " tol:" << surf1PointTol[pointi] << endl;
341
342 points1[meshPoints[pointi]] += rndVec;
343
344 // Mark edges affected by change to point
345 const labelList& pEdges = surf1.pointEdges()[pointi];
346
347 forAll(pEdges, i)
348 {
349 affectedEdges[pEdges[i]] = true;
350 }
351
352 hasPerturbed = true;
353
354 // Enough done for current edge; no need to test other intersections
355 // of this edge.
356 break;
357 }
358 }
359
360 return hasPerturbed;
361}
362
363
364// Perturb edge when close to face
365bool Foam::edgeIntersections::offsetPerturb
366(
367 const triSurface& surf1,
368 const triSurface& surf2,
369 const label edgeI,
370
371 Random& rndGen,
372 pointField& points1,
373 boolList& affectedEdges
374) const
375{
376 const labelList& meshPoints = surf1.meshPoints();
377
378 const edge& e = surf1.edges()[edgeI];
379
380 const List<pointIndexHit>& hits = operator[](edgeI);
381
382
383 bool hasPerturbed = false;
384
385 // For all hits on edge
386 for (const pointIndexHit& pHit : hits)
387 {
388 // Classify point on face of surface2
389 const label surf2Facei = pHit.index();
390
391 const triSurface::face_type& f2 = surf2.localFaces()[surf2Facei];
392 const pointField& surf2Pts = surf2.localPoints();
393
394 const point ctr = f2.centre(surf2Pts);
395
396 label nearType, nearLabel;
397
398 f2.nearestPointClassify(pHit.hitPoint(), surf2Pts, nearType, nearLabel);
399
400 if (nearType == triPointRef::POINT || nearType == triPointRef::EDGE)
401 {
402 // Shift edge towards tri centre
403 vector offset =
404 0.01*rndGen.sample01<scalar>()*(ctr - pHit.hitPoint());
405
406 // shift e[0]
407 points1[meshPoints[e[0]]] += offset;
408
409 // Mark edges affected by change to e0
410 const labelList& pEdges0 = surf1.pointEdges()[e[0]];
411
412 forAll(pEdges0, i)
413 {
414 affectedEdges[pEdges0[i]] = true;
415 }
416
417 // shift e[1]
418 points1[meshPoints[e[1]]] += offset;
419
420 // Mark edges affected by change to e1
421 const labelList& pEdges1 = surf1.pointEdges()[e[1]];
422
423 forAll(pEdges1, i)
424 {
425 affectedEdges[pEdges1[i]] = true;
426 }
427
428 hasPerturbed = true;
429
430 // No need to test any other hits on this edge
431 break;
432 }
433 }
434
435 return hasPerturbed;
436}
437
438
439// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
440
442:
444 classification_()
445{}
446
447
449(
450 const triSurface& surf1,
451 const triSurfaceSearch& query2,
452 const scalarField& surf1PointTol
453)
454:
455 List<List<pointIndexHit>>(surf1.nEdges()),
456 classification_(surf1.nEdges())
457{
458 checkEdges(surf1);
459
460 // Current set of edges to test
461 labelList edgesToTest(identity(surf1.nEdges()));
462
463 // Determine intersections for edgesToTest
464 intersectEdges
465 (
466 surf1,
467 surf1.points(), // surf1 meshPoints (not localPoints!)
468 query2,
469 surf1PointTol,
470 edgesToTest
471 );
472}
473
474
476(
477 const List<List<pointIndexHit>>& intersections,
478 const labelListList& classification
479)
480:
481 List<List<pointIndexHit>>(intersections),
482 classification_(classification)
483{}
484
485
486// * * * * * * * * * * * * * * * Static Functions * * * * * * * * * * * * * //
487
489{
490 const pointField& localPoints = surf.localPoints();
491 const labelListList& pointEdges = surf.pointEdges();
492 const edgeList& edges = surf.edges();
493
494 scalarField minLen(localPoints.size());
495
496 forAll(minLen, pointi)
497 {
498 const labelList& pEdges = pointEdges[pointi];
499
500 scalar minDist = GREAT;
501
502 forAll(pEdges, i)
503 {
504 minDist = min(minDist, edges[pEdges[i]].mag(localPoints));
505 }
506
507 minLen[pointi] = minDist;
508 }
509 return minLen;
510}
511
512
513// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
514
516(
517 const label nIters,
518 const triSurface& surf1,
519 const triSurfaceSearch& query2,
520 const scalarField& surf1PointTol,
521 pointField& points1
522)
523{
524 const triSurface& surf2 = query2.surface();
525
526 Random rndGen(356574);
527
528 // Current set of edges to (re)test
529 labelList edgesToTest(surf1.nEdges());
530
531 // Start off with all edges
532 forAll(edgesToTest, i)
533 {
534 edgesToTest[i] = i;
535 }
536
537
538 label iter = 0;
539
540 for (; iter < nIters; iter++)
541 {
542 // Go through all edges to (re)test and perturb points if they are
543 // degenerate hits. Mark off edges that need to be recalculated.
544
545 boolList affectedEdges(surf1.nEdges(), false);
546 label nShifted = 0;
547 label nRotated = 0;
548 label nOffset = 0;
549
550 forAll(edgesToTest, i)
551 {
552 label edgeI = edgesToTest[i];
553
554 // If edge not already marked for retesting
555 if (!affectedEdges[edgeI])
556 {
557 // 1. Check edges close to endpoint and perturb if necessary.
558
559 bool shiftedEdgeEndPoints =
560 inlinePerturb
561 (
562 surf1,
563 surf1PointTol,
564 edgeI,
565 rndGen,
566 points1,
567 affectedEdges
568 );
569
570 nShifted += (shiftedEdgeEndPoints ? 1 : 0);
571
572 if (!shiftedEdgeEndPoints)
573 {
574 bool rotatedEdge =
575 rotatePerturb
576 (
577 surf1,
578 surf1PointTol,
579 edgeI,
580 rndGen,
581 points1,
582 affectedEdges
583 );
584
585 nRotated += (rotatedEdge ? 1 : 0);
586
587 if (!rotatedEdge)
588 {
589 // 2. we're sure now that the edge actually pierces the
590 // face. Now check the face for intersections close its
591 // points/edges
592
593 bool offsetEdgePoints =
594 offsetPerturb
595 (
596 surf1,
597 surf2,
598 edgeI,
599 rndGen,
600 points1,
601 affectedEdges
602 );
603
604 nOffset += (offsetEdgePoints ? 1 : 0);
605 }
606 }
607 }
608 }
609
610 if (debug)
611 {
612 Pout<< "Edges to test : " << nl
613 << " total:" << edgesToTest.size() << nl
614 << " resolved by:" << nl
615 << " shifting : " << nShifted << nl
616 << " rotating : " << nRotated << nl
617 << " offsetting : " << nOffset << nl
618 << endl;
619 }
620
621
622 if (nShifted == 0 && nRotated == 0 && nOffset == 0)
623 {
624 // Nothing changed in current iteration. Current hit pattern is
625 // without any degenerates.
626 break;
627 }
628
629 // Repack affected edges
630 labelList newEdgesToTest(surf1.nEdges());
631 label newEdgeI = 0;
632
633 forAll(affectedEdges, edgeI)
634 {
635 if (affectedEdges[edgeI])
636 {
637 newEdgesToTest[newEdgeI++] = edgeI;
638 }
639 }
640 newEdgesToTest.setSize(newEdgeI);
641
642 if (debug)
643 {
644 Pout<< "Edges to test:" << nl
645 << " was : " << edgesToTest.size() << nl
646 << " is : " << newEdgesToTest.size() << nl
647 << endl;
648 }
649
650 // Transfer and test.
651 edgesToTest.transfer(newEdgesToTest);
652
653 if (edgesToTest.empty())
654 {
656 }
657
658 // Re intersect moved edges.
659 intersectEdges
660 (
661 surf1,
662 points1, // surf1 meshPoints (not localPoints!)
663 query2,
664 surf1PointTol,
665 edgesToTest
666 );
667 }
668
669 return iter;
670}
671
672
674(
675 const edgeIntersections& subInfo,
676 const labelList& edgeMap,
677 const labelList& faceMap,
678 const bool merge
679)
680{
681 forAll(subInfo, subI)
682 {
683 const List<pointIndexHit>& subHits = subInfo[subI];
684 const labelList& subClass = subInfo.classification()[subI];
685
686 label edgeI = edgeMap[subI];
687 List<pointIndexHit>& intersections = operator[](edgeI);
688 labelList& intersectionTypes = classification_[edgeI];
689
690 // Count unique hits. Assume edge can hit face only once
691 label sz = 0;
692 if (merge)
693 {
694 sz = intersections.size();
695 }
696
697 label nNew = 0;
698 forAll(subHits, i)
699 {
700 const pointIndexHit& subHit = subHits[i];
701
702 bool foundFace = false;
703 for (label interI = 0; interI < sz; interI++)
704 {
705 if (intersections[interI].index() == faceMap[subHit.index()])
706 {
707 foundFace = true;
708 break;
709 }
710 }
711 if (!foundFace)
712 {
713 nNew++;
714 }
715 }
716
717
718 intersections.setSize(sz+nNew);
719 intersectionTypes.setSize(sz+nNew);
720 nNew = sz;
721
722 forAll(subHits, i)
723 {
724 const pointIndexHit& subHit = subHits[i];
725
726 bool foundFace = false;
727 for (label interI = 0; interI < sz; interI++)
728 {
729 if (intersections[interI].index() == faceMap[subHit.index()])
730 {
731 foundFace = true;
732 break;
733 }
734 }
735
736 if (!foundFace)
737 {
738 intersections[nNew] = pointIndexHit
739 (
740 subHit.hit(),
741 subHit.rawPoint(),
742 faceMap[subHit.index()]
743 );
744 intersectionTypes[nNew] = subClass[i];
745 nNew++;
746 }
747 }
748 }
749}
750
751
752// ************************************************************************* //
label n
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
void transfer(List< T > &list)
Definition: List.C:447
void setSize(const label n)
Alias for resize()
Definition: List.H:218
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:66
const point_type & rawPoint() const noexcept
The point, no checks. Same as point()
label index() const noexcept
Return the hit index.
bool hit() const noexcept
Is there a hit?
label nEdges() const
Number of edges in patch.
const labelListList & pointEdges() const
Return point-edge addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const Field< point_type > & localPoints() const
Return pointField of points in patch.
const Field< point_type > & points() const noexcept
Return reference to global points.
Random number generator.
Definition: Random.H:60
int bit()
Return a random bit.
Definition: RandomI.H:38
Type sample01()
Return a sample whose components lie in the range [0,1].
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
iterator end() noexcept
Return an iterator to end traversing the UList.
Definition: UListI.H:350
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
const Vector< Cmpt > & centre(const Foam::UList< Vector< Cmpt > > &) const
Return *this (used for point which is a typedef to Vector<scalar>.
Definition: VectorI.H:114
virtual bool merge() const
Holder of intersections of edges of a surface with another surface. Optionally shuffles around points...
const labelListList & classification() const
For every intersection the classification status.
edgeIntersections()
Construct null.
static scalarField minEdgeLength(const triSurface &surf)
Calculate min edge length for every surface point.
label removeDegenerates(const label nIters, const triSurface &surf1, const triSurfaceSearch &query2, const scalarField &surf1PointTol, pointField &points1)
Resolve ties. Shuffles points so all edge - face intersections.
static scalar alignedCos_
Cosine between edge and face normal when considered parallel.
Helper class to search on triSurface.
const triSurface & surface() const
Return reference to the surface.
Triangulated surface description with patch information.
Definition: triSurface.H:79
labelledTri face_type
The face type (same as the underlying PrimitivePatch)
Definition: triSurface.H:209
@ POINT
Close to point.
Definition: triangle.H:96
@ EDGE
Close to edge.
Definition: triangle.H:97
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const labelList nEdges(UPstream::listGatherValues< label >(aMesh.nEdges()))
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
List< label > labelList
A List of labels.
Definition: List.H:66
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
vector point
Point is a vector.
Definition: point.H:43
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< bool > boolList
A List of bools.
Definition: List.H:64
error FatalError
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Vector< scalar > vector
Definition: vector.H:61
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:158
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Unit conversion functions.
Random rndGen
Definition: createFields.H:23