isoSurfacePointTemplates.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) 2018-2022 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 "isoSurfacePoint.H"
30#include "polyMesh.H"
31#include "syncTools.H"
32#include "surfaceFields.H"
33#include "OFstream.H"
34#include "meshTools.H"
35
36// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37
38template<class Type>
40Foam::isoSurfacePoint::adaptPatchFields
41(
43) const
44{
46 (
48 (
49 fld.name(),
50 fld.instance(),
51 fld.db(),
54 false
55 ),
56 fld, // internal field
57 true // preserveCouples
58 );
59 auto& sliceFld = tslice.ref();
60
61 const fvMesh& mesh = fld.mesh();
62
64
65 auto& sliceFldBf = sliceFld.boundaryFieldRef();
66
67 forAll(patches, patchi)
68 {
69 const polyPatch& pp = patches[patchi];
70
71 if
72 (
73 isA<emptyPolyPatch>(pp)
74 && pp.size() != sliceFldBf[patchi].size()
75 )
76 {
77 // Clear old value. Cannot resize it since is a slice.
78 sliceFldBf.set(patchi, nullptr);
79
80 // Set new value we can change
81 sliceFldBf.set
82 (
83 patchi,
85 (
86 mesh.boundary()[patchi],
87 sliceFld
88 )
89 );
90
91 // Note: cannot use patchInternalField since uses emptyFvPatch::size
92 // Do our own internalField instead.
93 const labelUList& faceCells =
94 mesh.boundary()[patchi].patch().faceCells();
95
96 Field<Type>& pfld = sliceFldBf[patchi];
97 pfld.setSize(faceCells.size());
99 {
100 pfld[i] = sliceFld[faceCells[i]];
101 }
102 }
103 else if (isA<cyclicPolyPatch>(pp))
104 {
105 // Already has interpolate as value
106 }
107 else if (isA<processorPolyPatch>(pp))
108 {
109 fvPatchField<Type>& pfld = const_cast<fvPatchField<Type>&>
110 (
111 sliceFldBf[patchi]
112 );
113
114 const scalarField& w = mesh.weights().boundaryField()[patchi];
115
117 w*pfld.patchInternalField()
118 + (1.0-w)*pfld.patchNeighbourField();
119
120 bitSet isCollocated
121 (
122 collocatedFaces(refCast<const processorPolyPatch>(pp))
123 );
124
125 forAll(isCollocated, i)
126 {
127 if (!isCollocated[i])
128 {
129 pfld[i] = f()[i];
130 }
131 }
132 }
133 }
134 return tslice;
135}
136
137
138template<class Type>
139Type Foam::isoSurfacePoint::generatePoint
140(
141 const scalar s0,
142 const Type& p0,
143 const bool hasSnap0,
144 const Type& snapP0,
145
146 const scalar s1,
147 const Type& p1,
148 const bool hasSnap1,
149 const Type& snapP1
150) const
151{
152 const scalar d = s1-s0;
153
154 if (mag(d) > VSMALL)
155 {
156 const scalar s = (iso_-s0)/d;
157
158 if (hasSnap1 && s >= 0.5 && s <= 1)
159 {
160 return snapP1;
161 }
162 else if (hasSnap0 && s >= 0.0 && s <= 0.5)
163 {
164 return snapP0;
165 }
166 else
167 {
168 return s*p1 + (1.0-s)*p0;
169 }
170 }
171 else
172 {
173 constexpr scalar s = 0.4999;
174
175 return s*p1 + (1.0-s)*p0;
176 }
177}
178
179
180template<class Type>
181void Foam::isoSurfacePoint::generateTriPoints
182(
183 const scalar s0,
184 const Type& p0,
185 const bool hasSnap0,
186 const Type& snapP0,
187
188 const scalar s1,
189 const Type& p1,
190 const bool hasSnap1,
191 const Type& snapP1,
192
193 const scalar s2,
194 const Type& p2,
195 const bool hasSnap2,
196 const Type& snapP2,
197
198 const scalar s3,
199 const Type& p3,
200 const bool hasSnap3,
201 const Type& snapP3,
202
203 DynamicList<Type>& pts
204) const
205{
206 // Note: cannot use simpler isoSurfaceCell::generateTriPoints since
207 // the need here to sometimes pass in remote 'snappedPoints'
208
209 int triIndex = 0;
210 if (s0 < iso_)
211 {
212 triIndex |= 1;
213 }
214 if (s1 < iso_)
215 {
216 triIndex |= 2;
217 }
218 if (s2 < iso_)
219 {
220 triIndex |= 4;
221 }
222 if (s3 < iso_)
223 {
224 triIndex |= 8;
225 }
226
227 // Form the vertices of the triangles for each case
228 switch (triIndex)
229 {
230 case 0x00:
231 case 0x0F:
232 break;
233
234 case 0x01:
235 case 0x0E:
236 {
237 pts.append
238 (
239 generatePoint(s0,p0,hasSnap0,snapP0,s1,p1,hasSnap1,snapP1)
240 );
241 pts.append
242 (
243 generatePoint(s0,p0,hasSnap0,snapP0,s2,p2,hasSnap2,snapP2)
244 );
245 pts.append
246 (
247 generatePoint(s0,p0,hasSnap0,snapP0,s3,p3,hasSnap3,snapP3)
248 );
249 if (triIndex == 0x0E)
250 {
251 // Flip normals
252 const label sz = pts.size();
253 std::swap(pts[sz-2], pts[sz-1]);
254 }
255 }
256 break;
257
258 case 0x02:
259 case 0x0D:
260 {
261 pts.append
262 (
263 generatePoint(s1,p1,hasSnap1,snapP1,s0,p0,hasSnap0,snapP0)
264 );
265 pts.append
266 (
267 generatePoint(s1,p1,hasSnap1,snapP1,s3,p3,hasSnap3,snapP3)
268 );
269 pts.append
270 (
271 generatePoint(s1,p1,hasSnap1,snapP1,s2,p2,hasSnap2,snapP2)
272 );
273 if (triIndex == 0x0D)
274 {
275 // Flip normals
276 const label sz = pts.size();
277 std::swap(pts[sz-2], pts[sz-1]);
278 }
279 }
280 break;
281
282 case 0x03:
283 case 0x0C:
284 {
285 Type p0p2 =
286 generatePoint(s0,p0,hasSnap0,snapP0,s2,p2,hasSnap2,snapP2);
287 Type p1p3 =
288 generatePoint(s1,p1,hasSnap1,snapP1,s3,p3,hasSnap3,snapP3);
289
290 pts.append
291 (
292 generatePoint(s0,p0,hasSnap0,snapP0,s3,p3,hasSnap3,snapP3)
293 );
294 pts.append(p1p3);
295 pts.append(p0p2);
296
297 pts.append(p1p3);
298 pts.append
299 (
300 generatePoint(s1,p1,hasSnap1,snapP1,s2,p2,hasSnap2,snapP2)
301 );
302 pts.append(p0p2);
303
304 if (triIndex == 0x0C)
305 {
306 // Flip normals
307 const label sz = pts.size();
308 std::swap(pts[sz-5], pts[sz-4]);
309 std::swap(pts[sz-2], pts[sz-1]);
310 }
311 }
312 break;
313
314 case 0x04:
315 case 0x0B:
316 {
317 pts.append
318 (
319 generatePoint(s2,p2,hasSnap2,snapP2,s0,p0,hasSnap0,snapP0)
320 );
321 pts.append
322 (
323 generatePoint(s2,p2,hasSnap2,snapP2,s1,p1,hasSnap1,snapP1)
324 );
325 pts.append
326 (
327 generatePoint(s2,p2,hasSnap2,snapP2,s3,p3,hasSnap3,snapP3)
328 );
329
330 if (triIndex == 0x0B)
331 {
332 // Flip normals
333 const label sz = pts.size();
334 std::swap(pts[sz-2], pts[sz-1]);
335 }
336 }
337 break;
338
339 case 0x05:
340 case 0x0A:
341 {
342 Type p0p1 =
343 generatePoint(s0,p0,hasSnap0,snapP0,s1,p1,hasSnap1,snapP1);
344 Type p2p3 =
345 generatePoint(s2,p2,hasSnap2,snapP2,s3,p3,hasSnap3,snapP3);
346
347 pts.append(p0p1);
348 pts.append(p2p3);
349 pts.append
350 (
351 generatePoint(s0,p0,hasSnap0,snapP0,s3,p3,hasSnap3,snapP3)
352 );
353
354 pts.append(p0p1);
355 pts.append
356 (
357 generatePoint(s1,p1,hasSnap1,snapP1,s2,p2,hasSnap2,snapP2)
358 );
359 pts.append(p2p3);
360
361 if (triIndex == 0x0A)
362 {
363 // Flip normals
364 const label sz = pts.size();
365 std::swap(pts[sz-5], pts[sz-4]);
366 std::swap(pts[sz-2], pts[sz-1]);
367 }
368 }
369 break;
370
371 case 0x06:
372 case 0x09:
373 {
374 Type p0p1 =
375 generatePoint(s0,p0,hasSnap0,snapP0,s1,p1,hasSnap1,snapP1);
376 Type p2p3 =
377 generatePoint(s2,p2,hasSnap2,snapP2,s3,p3,hasSnap3,snapP3);
378
379 pts.append(p0p1);
380 pts.append
381 (
382 generatePoint(s1,p1,hasSnap1,snapP1,s3,p3,hasSnap3,snapP3)
383 );
384 pts.append(p2p3);
385
386 pts.append(p0p1);
387 pts.append(p2p3);
388 pts.append
389 (
390 generatePoint(s0,p0,hasSnap0,snapP0,s2,p2,hasSnap2,snapP2)
391 );
392
393 if (triIndex == 0x09)
394 {
395 // Flip normals
396 const label sz = pts.size();
397 std::swap(pts[sz-5], pts[sz-4]);
398 std::swap(pts[sz-2], pts[sz-1]);
399 }
400 }
401 break;
402
403 case 0x08:
404 case 0x07:
405 {
406 pts.append
407 (
408 generatePoint(s3,p3,hasSnap3,snapP3,s0,p0,hasSnap0,snapP0)
409 );
410 pts.append
411 (
412 generatePoint(s3,p3,hasSnap3,snapP3,s2,p2,hasSnap2,snapP2)
413 );
414 pts.append
415 (
416 generatePoint(s3,p3,hasSnap3,snapP3,s1,p1,hasSnap1,snapP1)
417 );
418
419 if (triIndex == 0x07)
420 {
421 // Flip normals
422 const label sz = pts.size();
423 std::swap(pts[sz-2], pts[sz-1]);
424 }
425 }
426 break;
427 }
428}
429
430
431template<class Type>
432Foam::label Foam::isoSurfacePoint::generateFaceTriPoints
433(
434 const volScalarField& cVals,
435 const scalarField& pVals,
436
437 const VolumeField<Type>& cCoords,
438 const Field<Type>& pCoords,
439
440 const DynamicList<Type>& snappedPoints,
441 const labelList& snappedCc,
442 const labelList& snappedPoint,
443 const label facei,
444
445 const scalar neiVal,
446 const Type& neiPt,
447 const bool hasNeiSnap,
448 const Type& neiSnapPt,
449
451 DynamicList<label>& triMeshCells
452) const
453{
454 const label own = mesh_.faceOwner()[facei];
455
456 label oldNPoints = triPoints.size();
457
458 const face& f = mesh_.faces()[facei];
459
460 forAll(f, fp)
461 {
462 label pointi = f[fp];
463 label nextPointi = f[f.fcIndex(fp)];
464
465 generateTriPoints
466 (
467 pVals[pointi],
468 pCoords[pointi],
469 snappedPoint[pointi] != -1,
470 (
471 snappedPoint[pointi] != -1
472 ? snappedPoints[snappedPoint[pointi]]
473 : Type(Zero)
474 ),
475
476 pVals[nextPointi],
477 pCoords[nextPointi],
478 snappedPoint[nextPointi] != -1,
479 (
480 snappedPoint[nextPointi] != -1
481 ? snappedPoints[snappedPoint[nextPointi]]
482 : Type(Zero)
483 ),
484
485 cVals[own],
486 cCoords[own],
487 snappedCc[own] != -1,
488 (
489 snappedCc[own] != -1
490 ? snappedPoints[snappedCc[own]]
491 : Type(Zero)
492 ),
493
494 neiVal,
495 neiPt,
496 hasNeiSnap,
497 neiSnapPt,
498
500 );
501 }
502
503 // Every three triPoints is a triangle
504 label nTris = (triPoints.size()-oldNPoints)/3;
505 for (label i = 0; i < nTris; i++)
506 {
507 triMeshCells.append(own);
508 }
509
510 return nTris;
511}
512
513
514template<class Type>
515void Foam::isoSurfacePoint::generateTriPoints
516(
517 const volScalarField& cVals,
518 const scalarField& pVals,
519
520 const VolumeField<Type>& cCoords,
521 const Field<Type>& pCoords,
522
523 const DynamicList<Type>& snappedPoints,
524 const labelList& snappedCc,
525 const labelList& snappedPoint,
526
528 DynamicList<label>& triMeshCells
529) const
530{
531 const polyBoundaryMesh& patches = mesh_.boundaryMesh();
532 const labelList& own = mesh_.faceOwner();
533 const labelList& nei = mesh_.faceNeighbour();
534
535 if
536 (
537 (cVals.size() != mesh_.nCells())
538 || (pVals.size() != mesh_.nPoints())
539 || (cCoords.size() != mesh_.nCells())
540 || (pCoords.size() != mesh_.nPoints())
541 || (snappedCc.size() != mesh_.nCells())
542 || (snappedPoint.size() != mesh_.nPoints())
543 )
544 {
546 << "Incorrect size." << endl
547 << "mesh: nCells:" << mesh_.nCells()
548 << " points:" << mesh_.nPoints() << endl
549 << "cVals:" << cVals.size() << endl
550 << "cCoords:" << cCoords.size() << endl
551 << "snappedCc:" << snappedCc.size() << endl
552 << "pVals:" << pVals.size() << endl
553 << "pCoords:" << pCoords.size() << endl
554 << "snappedPoint:" << snappedPoint.size() << endl
555 << abort(FatalError);
556 }
557
558
559 // Generate triangle points
560
561 triPoints.clear();
562 triMeshCells.clear();
563
564 for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
565 {
566 if ((faceCutType_[facei] & cutType::ANYCUT) != 0)
567 {
568 generateFaceTriPoints
569 (
570 cVals,
571 pVals,
572
573 cCoords,
574 pCoords,
575
576 snappedPoints,
577 snappedCc,
578 snappedPoint,
579 facei,
580
581 cVals[nei[facei]],
582 cCoords[nei[facei]],
583 snappedCc[nei[facei]] != -1,
584 (
585 snappedCc[nei[facei]] != -1
586 ? snappedPoints[snappedCc[nei[facei]]]
587 : Type(Zero)
588 ),
589
590 triPoints,
591 triMeshCells
592 );
593 }
594 }
595
596
597 // Determine neighbouring snap status
598 boolList neiSnapped(mesh_.nBoundaryFaces(), false);
599 List<Type> neiSnappedPoint(neiSnapped.size(), Type(Zero));
600 for (const polyPatch& pp : patches)
601 {
602 if (pp.coupled())
603 {
604 label facei = pp.start();
605 forAll(pp, i)
606 {
607 label bFacei = facei-mesh_.nInternalFaces();
608 label snappedIndex = snappedCc[own[facei]];
609
610 if (snappedIndex != -1)
611 {
612 neiSnapped[bFacei] = true;
613 neiSnappedPoint[bFacei] = snappedPoints[snappedIndex];
614 }
615 facei++;
616 }
617 }
618 }
619 syncTools::swapBoundaryFaceList(mesh_, neiSnapped);
620 syncTools::swapBoundaryFaceList(mesh_, neiSnappedPoint);
621
622
623 forAll(patches, patchi)
624 {
625 const polyPatch& pp = patches[patchi];
626
627 if (isA<processorPolyPatch>(pp))
628 {
629 const processorPolyPatch& cpp =
630 refCast<const processorPolyPatch>(pp);
631
632 bitSet isCollocated(collocatedFaces(cpp));
633
634 forAll(isCollocated, i)
635 {
636 const label facei = pp.start()+i;
637
638 if ((faceCutType_[facei] & cutType::ANYCUT) != 0)
639 {
640 if (isCollocated[i])
641 {
642 generateFaceTriPoints
643 (
644 cVals,
645 pVals,
646
647 cCoords,
648 pCoords,
649
650 snappedPoints,
651 snappedCc,
652 snappedPoint,
653 facei,
654
655 cVals.boundaryField()[patchi][i],
656 cCoords.boundaryField()[patchi][i],
657 neiSnapped[facei-mesh_.nInternalFaces()],
658 neiSnappedPoint[facei-mesh_.nInternalFaces()],
659
660 triPoints,
661 triMeshCells
662 );
663 }
664 else
665 {
666 generateFaceTriPoints
667 (
668 cVals,
669 pVals,
670
671 cCoords,
672 pCoords,
673
674 snappedPoints,
675 snappedCc,
676 snappedPoint,
677 facei,
678
679 cVals.boundaryField()[patchi][i],
680 cCoords.boundaryField()[patchi][i],
681 false,
682 Type(Zero),
683
684 triPoints,
685 triMeshCells
686 );
687 }
688 }
689 }
690 }
691 else
692 {
693 label facei = pp.start();
694
695 forAll(pp, i)
696 {
697 if ((faceCutType_[facei] & cutType::ANYCUT) != 0)
698 {
699 generateFaceTriPoints
700 (
701 cVals,
702 pVals,
703
704 cCoords,
705 pCoords,
706
707 snappedPoints,
708 snappedCc,
709 snappedPoint,
710 facei,
711
712 cVals.boundaryField()[patchi][i],
713 cCoords.boundaryField()[patchi][i],
714 false, // fc not snapped
715 Type(Zero),
716
717 triPoints,
718 triMeshCells
719 );
720 }
721 facei++;
722 }
723 }
724 }
725
726 triPoints.shrink();
727 triMeshCells.shrink();
728}
729
730
731template<class Type>
734(
735 const label nPoints,
736 const labelList& triPointMergeMap,
737 const labelList& interpolatedPoints,
738 const List<FixedList<label, 3>>& interpolatedOldPoints,
740 const DynamicList<Type>& unmergedValues
741)
742{
743 // One value per point
744 auto tvalues = tmp<Field<Type>>::New(nPoints, Type(Zero));
745 auto& values = tvalues.ref();
746
747
748 // Pass1: unweighted average of merged point values
749 {
750 labelList nValues(values.size(), Zero);
751
752 forAll(unmergedValues, i)
753 {
754 label mergedPointi = triPointMergeMap[i];
755
756 if (mergedPointi >= 0)
757 {
758 values[mergedPointi] += unmergedValues[i];
759 nValues[mergedPointi]++;
760 }
761 }
762
763 forAll(values, i)
764 {
765 if (nValues[i] > 0)
766 {
767 values[i] /= scalar(nValues[i]);
768 }
769 }
770 }
771
772
773 // Pass2: weighted average for remaining values (from clipped triangles)
774
775 forAll(interpolatedPoints, i)
776 {
777 label pointi = interpolatedPoints[i];
778 const FixedList<label, 3>& oldPoints = interpolatedOldPoints[i];
780
781 // Note: zeroing should not be necessary if interpolation only done
782 // for newly introduced points (i.e. not in triPointMergeMap)
783 values[pointi] = Type(Zero);
784 forAll(oldPoints, j)
785 {
786 values[pointi] = w[j]*unmergedValues[oldPoints[j]];
787 }
788 }
789
790 return tvalues;
791}
792
793
794template<class Type>
796Foam::isoSurfacePoint::interpolateTemplate
797(
798 const VolumeField<Type>& cCoords,
799 const Field<Type>& pCoords
800) const
801{
802 // Recalculate boundary values
803 tmp<VolumeSliceField<Type>> c2(adaptPatchFields(cCoords));
804
805
806 DynamicList<Type> triPoints(3*nCutCells_);
807 DynamicList<label> triMeshCells(nCutCells_);
808
809 // Dummy snap data
810 DynamicList<Type> snappedPoints;
811 labelList snappedCc(mesh_.nCells(), -1);
812 labelList snappedPoint(mesh_.nPoints(), -1);
813
814 generateTriPoints
815 (
816 cValsPtr_(),
817 pVals_,
818
819 c2(),
820 pCoords,
821
822 snappedPoints,
823 snappedCc,
824 snappedPoint,
825
826 triPoints,
827 triMeshCells
828 );
829
830 return interpolate
831 (
832 this->points().size(),
833 triPointMergeMap_,
834 interpolatedPoints_,
835 interpolatedOldPoints_,
836 interpolationWeights_,
838 );
839}
840
841
842// ************************************************************************* //
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:391
DynamicList< T, SizeMin > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:434
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
Generic templated field type.
Definition: Field.H:82
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: FixedList.H:81
static constexpr label size() noexcept
Return the number of elements in the FixedList.
Definition: FixedList.H:416
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
void setSize(const label n)
Alias for resize()
Definition: List.H:218
const T * set(const label i) const
Definition: PtrList.H:138
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label fcIndex(const label i) const noexcept
Definition: UListI.H:60
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:453
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:237
Abstract base class for interpolating in 1D.
static autoPtr< isoSurfaceBase > New(const isoSurfaceParams &params, const volScalarField &cellValues, const scalarField &pointValues, const scalar iso, const bitSet &ignoreCells=bitSet())
Create for specified algorithm type.
const polyMesh & mesh() const noexcept
The mesh for which the iso-surface is associated.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
bool interpolate() const noexcept
Same as isPointData()
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:445
A class for managing temporary objects.
Definition: tmp.H:65
Triangle storage. Null constructable (unfortunately triangle<point, point> is not)
Definition: triPoints.H:55
const volScalarField & p0
Definition: EEqn.H:36
const polyBoundaryMesh & patches
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
label nPoints
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))
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
List< bool > boolList
A List of bools.
Definition: List.H:64
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.