Loading...
Searching...
No Matches
surfaceInflate.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) 2015 OpenFOAM Foundation
9 Copyright (C) 2020-2023 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
27Application
28 surfaceInflate
29
30Group
31 grpSurfaceUtilities
32
33Description
34 Inflates surface. WIP. Checks for overlaps and locally lowers
35 inflation distance
36
37Usage
38 - surfaceInflate [OPTION]
39
40 \param -checkSelfIntersection \n
41 Includes checks for self-intersection.
42
43 \param -nSmooth
44 Specify number of smoothing iterations
45
46 \param -featureAngle
47 Specify a feature angle
48
49
50 E.g. inflate surface by 20mm with 1.5 safety factor:
51 surfaceInflate DTC-scaled.obj 0.02 1.5 -featureAngle 45 -nSmooth 2
52
53\*---------------------------------------------------------------------------*/
54
55#include "argList.H"
56#include "Time.H"
57#include "triSurfaceFields.H"
58#include "triSurfaceMesh.H"
59#include "triSurfaceGeoMesh.H"
60#include "bitSet.H"
61#include "OBJstream.H"
62#include "surfaceFeatures.H"
63
64using namespace Foam;
65
66// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67
68scalar calcVertexNormalWeight
69(
70 const triFace& f,
71 const label pI,
72 const vector& fN,
73 const pointField& points
74)
75{
76 label index = f.find(pI);
77
78 if (index == -1)
79 {
81 << "Point not in face" << abort(FatalError);
82 return 0;
83 }
84
85 const vector e1 = points[f[index]] - points[f[f.fcIndex(index)]];
86 const vector e2 = points[f[index]] - points[f[f.rcIndex(index)]];
87
88 return mag(fN)/(magSqr(e1)*magSqr(e2) + VSMALL);
89}
90
91
92tmp<vectorField> calcVertexNormals(const triSurface& surf)
93{
94 // Weighted average of normals of faces attached to the vertex
95 // Weight = fA / (mag(e0)^2 * mag(e1)^2);
96 auto tpointNormals = tmp<vectorField>::New(surf.nPoints(), Zero);
97 auto& pointNormals = tpointNormals.ref();
98
99 const pointField& points = surf.points();
100 const labelListList& pointFaces = surf.pointFaces();
101 const labelList& meshPoints = surf.meshPoints();
102
103 forAll(pointFaces, pI)
104 {
105 const labelList& pFaces = pointFaces[pI];
106
107 forAll(pFaces, fI)
108 {
109 const label faceI = pFaces[fI];
110 const triFace& f = surf[faceI];
111
112 vector areaNorm = f.areaNormal(points);
113
114 scalar weight = calcVertexNormalWeight
115 (
116 f,
117 meshPoints[pI],
118 areaNorm,
119 points
120 );
121
122 pointNormals[pI] += weight * areaNorm;
123 }
124
125 pointNormals[pI].normalise();
126 }
127
128 return tpointNormals;
129}
130
131
132tmp<vectorField> calcPointNormals
133(
134 const triSurface& s,
135 const bitSet& isFeaturePoint,
136 const List<surfaceFeatures::edgeStatus>& edgeStat
137)
138{
139 //const pointField pointNormals(s.pointNormals());
140 tmp<vectorField> tpointNormals(calcVertexNormals(s));
141 vectorField& pointNormals = tpointNormals.ref();
142
143
144 // feature edges: create edge normals from edgeFaces only.
145 {
146 const labelListList& edgeFaces = s.edgeFaces();
147
148 labelList nNormals(s.nPoints(), Zero);
149 forAll(edgeStat, edgeI)
150 {
151 if (edgeStat[edgeI] != surfaceFeatures::NONE)
152 {
153 const edge& e = s.edges()[edgeI];
154 forAll(e, i)
155 {
156 if (!isFeaturePoint.test(e[i]))
157 {
158 pointNormals[e[i]] = Zero;
159 }
160 }
161 }
162 }
163
164 forAll(edgeStat, edgeI)
165 {
166 if (edgeStat[edgeI] != surfaceFeatures::NONE)
167 {
168 const labelList& eFaces = edgeFaces[edgeI];
169
170 // Get average edge normal
171 vector n = Zero;
172 for (const label facei : eFaces)
173 {
174 n += s.faceNormals()[facei];
175 }
176 n /= eFaces.size();
177
178
179 // Sum to points
180 const edge& e = s.edges()[edgeI];
181 forAll(e, i)
182 {
183 if (!isFeaturePoint.test(e[i]))
184 {
185 pointNormals[e[i]] += n;
186 nNormals[e[i]]++;
187 }
188 }
189 }
190 }
191
192 forAll(nNormals, pointI)
193 {
194 if (nNormals[pointI] > 0)
195 {
196 pointNormals[pointI].normalise();
197 }
198 }
199 }
200
201
202 forAll(pointNormals, pointI)
203 {
204 if (mag(mag(pointNormals[pointI])-1) > SMALL)
205 {
207 << "unitialised"
208 << exit(FatalError);
209 }
210 }
211
212 return tpointNormals;
213}
214
215
216void detectSelfIntersections
217(
218 const triSurfaceMesh& s,
219 bitSet& isEdgeIntersecting
220)
221{
222 const edgeList& edges = s.edges();
224 const labelList& meshPoints = s.meshPoints();
225 const tmp<pointField> tpoints(s.points());
226 const pointField& points = tpoints();
227
228 isEdgeIntersecting.resize_nocopy(edges.size());
229 isEdgeIntersecting = false;
230
231 forAll(edges, edgeI)
232 {
233 const edge& e = edges[edgeI];
234
235 pointIndexHit hitInfo
236 (
237 tree.findLine
238 (
239 points[meshPoints[e[0]]],
240 points[meshPoints[e[1]]],
241 treeDataTriSurface::findSelfIntersectOp(tree, edgeI)
242 )
243 );
244
245 if (hitInfo.hit())
246 {
247 isEdgeIntersecting.set(edgeI);
248 }
249 }
250}
251
252
253label detectIntersectionPoints
254(
255 const scalar tolerance,
256 const triSurfaceMesh& s,
257
258 const scalar extendFactor,
259 const pointField& initialPoints,
260 const vectorField& displacement,
261
262 const bool checkSelfIntersect,
263 const bitSet& initialIsEdgeIntersecting,
264
265 bitSet& isPointOnHitEdge,
266 scalarField& scale
267)
268{
269 const pointField initialLocalPoints(initialPoints, s.meshPoints());
270 const List<labelledTri>& localFaces = s.localFaces();
271
272
273 label nHits = 0;
274 isPointOnHitEdge.setSize(s.nPoints());
275 isPointOnHitEdge = false;
276
277
278 // 1. Extrusion offset vectors intersecting new surface location
279 {
280 scalar tol = Foam::max(tolerance, 10*s.tolerance());
281
282 // Collect all the edge vectors. Slightly shorten the edges to prevent
283 // finding lots of intersections. The fast triangle intersection routine
284 // has problems with rays passing through a point of the
285 // triangle.
286 pointField start(initialLocalPoints+tol*displacement);
287 pointField end(initialLocalPoints+extendFactor*displacement);
288
289 List<pointIndexHit> hits;
290 s.findLineAny(start, end, hits);
291
292 forAll(hits, pointI)
293 {
294 if
295 (
296 hits[pointI].hit()
297 && !localFaces[hits[pointI].index()].found(pointI)
298 )
299 {
300 scale[pointI] = Foam::max(0.0, scale[pointI]-0.2);
301
302 isPointOnHitEdge.set(pointI);
303 nHits++;
304 }
305 }
306 }
307
308
309 // 2. (new) surface self intersections
310 if (checkSelfIntersect)
311 {
312 bitSet isEdgeIntersecting;
313 detectSelfIntersections(s, isEdgeIntersecting);
314
315 const edgeList& edges = s.edges();
316 const tmp<pointField> tpoints(s.points());
317 const pointField& points = tpoints();
318
319 forAll(edges, edgeI)
320 {
321 const edge& e = edges[edgeI];
322
323 if (isEdgeIntersecting[edgeI] && !initialIsEdgeIntersecting[edgeI])
324 {
325 if (isPointOnHitEdge.set(e[0]))
326 {
327 label start = s.meshPoints()[e[0]];
328 const point& pt = points[start];
329
330 Pout<< "Additional self intersection at "
331 << pt
332 << endl;
333
334 scale[e[0]] = Foam::max(0.0, scale[e[0]]-0.2);
335 nHits++;
336 }
337 if (isPointOnHitEdge.set(e[1]))
338 {
339 label end = s.meshPoints()[e[1]];
340 const point& pt = points[end];
341
342 Pout<< "Additional self intersection at "
343 << pt
344 << endl;
345
346 scale[e[1]] = Foam::max(0.0, scale[e[1]]-0.2);
347 nHits++;
348 }
349 }
350 }
351 }
352
353 return nHits;
354}
355
356
358(
359 const triSurface& s,
360 const scalarField& fld,
361 const scalarField& edgeWeights
362)
363{
364 auto tres = tmp<scalarField>::New(s.nPoints(), Zero);
365 auto& res = tres.ref();
366
367 scalarField sumWeight(s.nPoints(), Zero);
368
369 const edgeList& edges = s.edges();
370
371 forAll(edges, edgeI)
372 {
373 const edge& e = edges[edgeI];
374 const scalar w = edgeWeights[edgeI];
375
376 res[e[0]] += w*fld[e[1]];
377 sumWeight[e[0]] += w;
378
379 res[e[1]] += w*fld[e[0]];
380 sumWeight[e[1]] += w;
381 }
382
383 // Average
384 // ~~~~~~~
385
386 forAll(res, pointI)
387 {
388 if (mag(sumWeight[pointI]) < VSMALL)
389 {
390 // Unconnected point. Take over original value
391 res[pointI] = fld[pointI];
392 }
393 else
394 {
395 res[pointI] /= sumWeight[pointI];
396 }
397 }
398
399 return tres;
400}
401
402
403void minSmooth
404(
405 const triSurface& s,
406 const bitSet& isAffectedPoint,
407 const scalarField& fld,
408 scalarField& newFld
409)
410{
411
412 const edgeList& edges = s.edges();
413 const pointField& points = s.points();
414 const labelList& mp = s.meshPoints();
415
416 scalarField edgeWeights(edges.size());
417 forAll(edges, edgeI)
418 {
419 const edge& e = edges[edgeI];
420 scalar w = mag(points[mp[e[0]]]-points[mp[e[1]]]);
421
422 edgeWeights[edgeI] = 1.0/(Foam::max(w, SMALL));
423 }
424
425 tmp<scalarField> tavgFld = avg(s, fld, edgeWeights);
426
427 const scalarField& avgFld = tavgFld();
428
429 forAll(fld, pointI)
430 {
431 if (isAffectedPoint.test(pointI))
432 {
433 newFld[pointI] = Foam::min
434 (
435 fld[pointI],
436 0.5*fld[pointI] + 0.5*avgFld[pointI]
437 //avgFld[pointI]
438 );
439 }
440 }
441}
442
443
444void lloydsSmoothing
445(
446 const label nSmooth,
447 triSurface& s,
448 const bitSet& isFeaturePoint,
449 const List<surfaceFeatures::edgeStatus>& edgeStat,
450 const bitSet& isAffectedPoint
451)
452{
453 const labelList& meshPoints = s.meshPoints();
454 const edgeList& edges = s.edges();
455
456
457 bitSet isSmoothPoint(isAffectedPoint);
458 // Extend isSmoothPoint
459 {
460 bitSet newIsSmoothPoint(isSmoothPoint);
461 forAll(edges, edgeI)
462 {
463 const edge& e = edges[edgeI];
464 if (isSmoothPoint.test(e[0]))
465 {
466 newIsSmoothPoint.set(e[1]);
467 }
468 if (isSmoothPoint.test(e[1]))
469 {
470 newIsSmoothPoint.set(e[0]);
471 }
472 }
473 Info<< "From points-to-smooth " << isSmoothPoint.count()
474 << " to " << newIsSmoothPoint.count()
475 << endl;
476 isSmoothPoint.transfer(newIsSmoothPoint);
477 }
478
479 // Do some smoothing (Lloyds algorithm) around problematic points
480 for (label i = 0; i < nSmooth; i++)
481 {
482 const labelListList& pointFaces = s.pointFaces();
483 const pointField& faceCentres = s.faceCentres();
484
485 pointField newPoints(s.points());
486 forAll(isSmoothPoint, pointI)
487 {
488 if (isSmoothPoint[pointI] && !isFeaturePoint[pointI])
489 {
490 const labelList& pFaces = pointFaces[pointI];
491
492 point avg(Zero);
493 forAll(pFaces, pFaceI)
494 {
495 avg += faceCentres[pFaces[pFaceI]];
496 }
497 avg /= pFaces.size();
498 newPoints[meshPoints[pointI]] = avg;
499 }
500 }
501
502
503 // Move points on feature edges only according to feature edges.
504
505 const pointField& points = s.points();
506
507 vectorField pointSum(s.nPoints(), Zero);
508 labelList nPointSum(s.nPoints(), Zero);
509
510 forAll(edges, edgeI)
511 {
512 if (edgeStat[edgeI] != surfaceFeatures::NONE)
513 {
514 const edge& e = edges[edgeI];
515 const point& start = points[meshPoints[e[0]]];
516 const point& end = points[meshPoints[e[1]]];
517 point eMid = 0.5*(start+end);
518 pointSum[e[0]] += eMid;
519 nPointSum[e[0]]++;
520 pointSum[e[1]] += eMid;
521 nPointSum[e[1]]++;
522 }
523 }
524
525 forAll(pointSum, pointI)
526 {
527 if
528 (
529 isSmoothPoint[pointI]
530 && isFeaturePoint[pointI]
531 && nPointSum[pointI] >= 2
532 )
533 {
534 newPoints[meshPoints[pointI]] =
535 pointSum[pointI]/nPointSum[pointI];
536 }
537 }
538
539
540 s.movePoints(newPoints);
541
542
543 // Extend isSmoothPoint
544 {
545 bitSet newIsSmoothPoint(isSmoothPoint);
546 forAll(edges, edgeI)
547 {
548 const edge& e = edges[edgeI];
549 if (isSmoothPoint[e[0]])
550 {
551 newIsSmoothPoint.set(e[1]);
552 }
553 if (isSmoothPoint[e[1]])
554 {
555 newIsSmoothPoint.set(e[0]);
556 }
557 }
558 Info<< "From points-to-smooth " << isSmoothPoint.count()
559 << " to " << newIsSmoothPoint.count()
560 << endl;
561 isSmoothPoint.transfer(newIsSmoothPoint);
562 }
563 }
564}
565
566
567
568// Main program:
569
570int main(int argc, char *argv[])
571{
573 (
574 "Inflates surface according to point normals."
575 );
576
579 (
580 "Creates inflated version of surface using point normals. "
581 "Takes surface, distance to inflate and additional safety factor"
582 );
584 (
585 "checkSelfIntersection",
586 "Also check for self-intersection"
587 );
589 (
590 "nSmooth",
591 "integer",
592 "Number of smoothing iterations (default 20)"
593 );
595 (
596 "featureAngle",
597 "scalar",
598 "Feature angle"
599 );
601 (
602 "debug",
603 "Switch on additional debug information"
604 );
605
606 argList::addArgument("input", "The input surface file");
607 argList::addArgument("distance", "The inflate distance");
608 argList::addArgument("factor", "The extend safety factor [1,10]");
609
610 argList::noFunctionObjects(); // Never use function objects
611
612 #include "setRootCase.H"
613 #include "createTime.H"
614
615 const auto inputName = args.get<word>(1);
616 const auto distance = args.get<scalar>(2);
617 const auto extendFactor = args.get<scalar>(3);
618 const bool checkSelfIntersect = args.found("checkSelfIntersection");
619 const auto nSmooth = args.getOrDefault<label>("nSmooth", 10);
620 const auto featureAngle = args.getOrDefault<scalar>("featureAngle", 180);
621 const bool debug = args.found("debug");
622
623
624 Info<< "Inflating surface " << inputName
625 << " according to point normals" << nl
626 << " distance : " << distance << nl
627 << " safety factor : " << extendFactor << nl
628 << " self intersection test : " << checkSelfIntersect << nl
629 << " smoothing iterations : " << nSmooth << nl
630 << " feature angle : " << featureAngle << nl
631 << endl;
632
633
634 if (extendFactor < 1 || extendFactor > 10)
635 {
637 << "Illegal safety factor " << extendFactor
638 << ". It is usually 1..2"
639 << exit(FatalError);
640 }
641
642
643
644 // Load triSurfaceMesh
646 (
648 (
649 inputName, // name
650 runTime.constant(), // instance
651 "triSurface", // local
652 runTime, // registry
655 )
656 );
657
658 // Mark features
659 const surfaceFeatures features(s, 180.0-featureAngle);
660
661 Info<< "Detected features:" << nl
662 << " feature points : " << features.featurePoints().size()
663 << " out of " << s.nPoints() << nl
664 << " feature edges : " << features.featureEdges().size()
665 << " out of " << s.nEdges() << nl
666 << endl;
667
668 bitSet isFeaturePoint(s.nPoints(), features.featurePoints());
669
670 const List<surfaceFeatures::edgeStatus> edgeStat(features.toStatus());
671
672
673 const pointField initialPoints(s.points());
674
675
676 // Construct scale
677 Info<< "Constructing field scale\n" << endl;
679 (
681 (
682 "scale", // name
683 runTime.timeName(), // instance
684 s, // registry
687 ),
688 s,
689 scalar(1),
691 );
692
693
694 // Construct unit normals
695
696 Info<< "Calculating vertex normals\n" << endl;
697 const pointField pointNormals
698 (
699 calcPointNormals
700 (
701 s,
702 isFeaturePoint,
703 edgeStat
704 )
705 );
706
707
708 // Construct pointDisplacement
709 Info<< "Constructing field pointDisplacement\n" << endl;
710 triSurfacePointVectorField pointDisplacement
711 (
713 (
714 "pointDisplacement", // name
715 runTime.timeName(), // instance
716 s, // registry
719 ),
720 s,
721 dimLength,
722 // - calculate with primitive fields (#2758)
723 (distance*scale.field()*pointNormals)
724 );
725
726
727 const labelList& meshPoints = s.meshPoints();
728
729
730 // Any point on any intersected edge in any of the iterations
731 bitSet isScaledPoint(s.nPoints());
732
733
734 // Detect any self intersections on initial mesh
735 bitSet initialIsEdgeIntersecting;
736 if (checkSelfIntersect)
737 {
738 detectSelfIntersections(s, initialIsEdgeIntersecting);
739 }
740 else
741 {
742 // Mark all edges as already self intersecting so avoid detecting any
743 // new ones
744 initialIsEdgeIntersecting.setSize(s.nEdges(), true);
745 }
746
747
748 // Inflate
749 while (runTime.loop())
750 {
751 Info<< "Time = " << runTime.timeName() << nl << endl;
752
753 // Move to new location
754 pointField newPoints(initialPoints);
755 forAll(meshPoints, i)
756 {
757 newPoints[meshPoints[i]] += pointDisplacement[i];
758 }
759 s.movePoints(newPoints);
760 Info<< "Bounding box : " << s.bounds() << endl;
761
762
763 // Work out scale from pointDisplacement
764 forAll(scale, pointI)
765 {
766 if (s.pointFaces()[pointI].size())
767 {
768 scale[pointI] = mag(pointDisplacement[pointI])/distance;
769 }
770 else
771 {
772 scale[pointI] = 1.0;
773 }
774 }
775
776
777 Info<< "Scale : " << gAverage(scale) << endl;
778 Info<< "Displacement : " << gAverage(pointDisplacement) << endl;
779
780
781
782 // Detect any intersections and scale back
783 bitSet isAffectedPoint;
784 label nIntersections = detectIntersectionPoints
785 (
786 1e-9, // intersection tolerance
787 s,
788 extendFactor,
789 initialPoints,
790 pointDisplacement,
791
792 checkSelfIntersect,
793 initialIsEdgeIntersecting,
794
795 isAffectedPoint,
796 scale
797 );
798 Info<< "Detected " << nIntersections << " intersections" << nl
799 << endl;
800
801 if (nIntersections == 0)
802 {
803 runTime.write();
804 break;
805 }
806
807
808 // Accumulate all affected points
809 forAll(isAffectedPoint, pointI)
810 {
811 if (isAffectedPoint.test(pointI))
812 {
813 isScaledPoint.set(pointI);
814 }
815 }
816
817 // Smear out lowering of scale so any edges not found are
818 // still included
819 for (label i = 0; i < nSmooth; i++)
820 {
821 triSurfacePointScalarField oldScale(scale);
822 oldScale.rename("oldScale");
823 minSmooth
824 (
825 s,
826 bitSet(s.nPoints(), true),
827 oldScale,
828 scale
829 );
830 }
831
832
833 // From scale update the pointDisplacement
834 // - calculate with primitive fields (#2758)
835 pointDisplacement.normalise();
836 pointDisplacement.field() *= distance*scale.field();
837
838
839 // Do some smoothing (Lloyds algorithm)
840 lloydsSmoothing(nSmooth, s, isFeaturePoint, edgeStat, isAffectedPoint);
841
842 // Update pointDisplacement
843 const tmp<pointField> tpoints(s.points());
844 const pointField& pts = tpoints();
845
846 forAll(meshPoints, i)
847 {
848 label meshPointI = meshPoints[i];
849 pointDisplacement[i] = pts[meshPointI]-initialPoints[meshPointI];
850 }
851
852
853 // Write
854 runTime.write();
855
856 runTime.printExecutionTime(Info);
857 }
858
859
860 Info<< "Detected overall intersecting points " << isScaledPoint.count()
861 << " out of " << s.nPoints() << nl << endl;
862
863
864 if (debug)
865 {
866 OBJstream str(runTime.path()/"isScaledPoint.obj");
867
868 for (const label pointi : isScaledPoint)
869 {
870 str.write(initialPoints[meshPoints[pointi]]);
871 }
872 }
873
874
875 Info<< "End\n" << endl;
876
877 return 0;
878}
879
880
881// ************************************************************************* //
bool found
label n
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
void normalise()
Inplace normalise this field. Generally a no-op except for vector fields.
Definition Field.C:600
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
void resize_nocopy(const label numElem)
Currently identical to resize. Subject to future change (Oct-2021).
void setSize(const label n, unsigned int val=0u)
Alias for resize().
Definition PackedList.H:819
label nPoints() const
Number of points supporting patch faces.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & points() const noexcept
Return reference to global points.
const labelListList & pointFaces() const
Return point-face addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
Definition argList.C:562
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition argList.C:366
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition argList.C:389
static void noParallel()
Remove the parallel options.
Definition argList.C:599
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition argList.C:400
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
Definition bitSet.H:334
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
Non-pointer based hierarchical recursive searching.
Holds feature edges/points of surface.
@ NONE
Not a classified feature edge.
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
IOoject and searching on triSurface.
Triangulated surface description with patch information.
Definition triSurface.H:74
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
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
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))
const char * end
Definition SVGTools.H:223
const dimensionedScalar mp
Proton mass.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
scalar distance(const vector &p1, const vector &p2)
Definition curveTools.C:12
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
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere).
DimensionedField< vector, triSurfacePointGeoMesh > triSurfacePointVectorField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
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:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
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...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
DimensionedField< scalar, triSurfacePointGeoMesh > triSurfacePointScalarField
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} pointMap[start]=pointMap[end]=Foam::min(start, end);} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};constexpr label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< DynamicList< face > > pFaces[nBCs]
face triFace(3)
labelList f(nPoints)
const pointField & pts
Tree tree(triangles.begin(), triangles.end())
Foam::argList args(argc, argv)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Fields for triSurface.