CV2D.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) 2013-2016 OpenFOAM Foundation
9 Copyright (C) 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 "CV2D.H"
30#include "Random.H"
31#include "transform.H"
32#include "IFstream.H"
33#include "uint.H"
34
35namespace Foam
36{
38}
39
40// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41
42void Foam::CV2D::insertBoundingBox()
43{
44 Info<< "insertBoundingBox: creating bounding mesh" << endl;
45 scalar bigSpan = 10*meshControls().span();
46 insertPoint(point2D(-bigSpan, -bigSpan), Vb::FAR_POINT);
47 insertPoint(point2D(-bigSpan, bigSpan), Vb::FAR_POINT);
48 insertPoint(point2D(bigSpan, -bigSpan), Vb::FAR_POINT);
49 insertPoint(point2D(bigSpan, bigSpan), Vb::FAR_POINT);
50}
51
52
53void Foam::CV2D::fast_restore_Delaunay(Vertex_handle vh)
54{
55 int i;
56 Face_handle f = vh->face(), next, start(f);
57
58 do
59 {
60 i=f->index(vh);
61 if (!is_infinite(f))
62 {
63 if (!internal_flip(f, cw(i))) external_flip(f, i);
64 if (f->neighbor(i) == start) start = f;
65 }
66 f = f->neighbor(cw(i));
67 } while (f != start);
68}
69
70
71void Foam::CV2D::external_flip(Face_handle& f, int i)
72{
73 Face_handle n = f->neighbor(i);
74
75 if
76 (
77 CGAL::ON_POSITIVE_SIDE
78 != side_of_oriented_circle(n, f->vertex(i)->point())
79 ) return;
80
81 flip(f, i);
82 i = n->index(f->vertex(i));
83 external_flip(n, i);
84}
85
86
87bool Foam::CV2D::internal_flip(Face_handle& f, int i)
88{
89 Face_handle n = f->neighbor(i);
90
91 if
92 (
93 CGAL::ON_POSITIVE_SIDE
94 != side_of_oriented_circle(n, f->vertex(i)->point())
95 )
96 {
97 return false;
98 }
99
100 flip(f, i);
101
102 return true;
103}
104
105
106// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
107
109(
110 const Time& runTime,
111 const dictionary& cvMeshDict
112)
113:
114 Delaunay(),
115 runTime_(runTime),
116 rndGen_(64293*Pstream::myProcNo()),
117 allGeometry_
118 (
119 IOobject
120 (
121 "cvSearchableSurfaces",
122 runTime_.constant(),
123 "triSurface",
124 runTime_,
125 IOobject::MUST_READ,
126 IOobject::NO_WRITE
127 ),
128 cvMeshDict.subDict("geometry"),
129 cvMeshDict.getOrDefault("singleRegionName", true)
130 ),
131 qSurf_
132 (
133 runTime_,
134 rndGen_,
135 allGeometry_,
136 cvMeshDict.subDict("surfaceConformation")
137 ),
138 controls_(cvMeshDict, qSurf_.globalBounds()),
139 cellSizeControl_
140 (
141 runTime,
142 cvMeshDict.subDict("motionControl").subDict("shapeControlFunctions"),
143 qSurf_,
144 controls_.minCellSize()
145 ),
146 relaxationModel_
147 (
148 relaxationModel::New
149 (
150 cvMeshDict.subDict("motionControl"),
151 runTime
152 )
153 ),
154 z_
155 (
156 cvMeshDict.subDict("surfaceConformation")
157 .get<Foam::point>("locationInMesh").z()
158 ),
159 startOfInternalPoints_(0),
160 startOfSurfacePointPairs_(0),
161 startOfBoundaryConformPointPairs_(0),
162 featurePoints_()
163{
164 Info<< meshControls() << endl;
165
166 insertBoundingBox();
167 insertFeaturePoints();
168}
169
170
171// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
172
174{}
175
176
177// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
178
180(
181 const point2DField& points,
182 const scalar nearness
183)
184{
185 Info<< "insertInitialPoints(const point2DField& points): ";
186
187 startOfInternalPoints_ = number_of_vertices();
188 label nVert = startOfInternalPoints_;
189
190 // Add the points and index them
191 for (const point2D& p : points)
192 {
193 if (qSurf_.wellInside(toPoint3D(p), nearness))
194 {
195 insert(toPoint(p))->index() = nVert++;
196 }
197 else
198 {
199 Warning
200 << "Rejecting point " << p << " outside surface" << endl;
201 }
202 }
203
204 Info<< nVert << " vertices inserted" << endl;
205
206 if (meshControls().objOutput())
207 {
208 // Checking validity of triangulation
209 assert(is_valid());
210
211 writeTriangles("initial_triangles.obj", true);
212 writeFaces("initial_faces.obj", true);
213 }
214}
215
216
217void Foam::CV2D::insertPoints(const fileName& pointFileName)
218{
219 IFstream pointsFile(pointFileName);
220
221 if (pointsFile.good())
222 {
223 insertPoints
224 (
225 point2DField(pointsFile),
226 0.5*meshControls().minCellSize2()
227 );
228 }
229 else
230 {
232 << "Could not open pointsFile " << pointFileName
233 << exit(FatalError);
234 }
235}
236
237
239{
240 Info<< "insertInitialGrid: ";
241
242 startOfInternalPoints_ = number_of_vertices();
243 label nVert = startOfInternalPoints_;
244
245 scalar x0 = qSurf_.globalBounds().min().x();
246 scalar xR = qSurf_.globalBounds().max().x() - x0;
247 int ni = int(xR/meshControls().minCellSize()) + 1;
248 scalar deltax = xR/ni;
249
250 scalar y0 = qSurf_.globalBounds().min().y();
251 scalar yR = qSurf_.globalBounds().max().y() - y0;
252 int nj = int(yR/meshControls().minCellSize()) + 1;
253 scalar deltay = yR/nj;
254
255 Random rndGen(1321);
256 scalar pert = meshControls().randomPerturbation()*min(deltax, deltay);
257
258 for (int i=0; i<ni; i++)
259 {
260 for (int j=0; j<nj; j++)
261 {
262 Foam::point p(x0 + i*deltax, y0 + j*deltay, 0);
263
264 if (meshControls().randomiseInitialGrid())
265 {
266 p.x() += pert*(rndGen.sample01<scalar>() - 0.5);
267 p.y() += pert*(rndGen.sample01<scalar>() - 0.5);
268 }
269
270 if (qSurf_.wellInside(p, 0.5*meshControls().minCellSize2()))
271 {
272 insert(Point(p.x(), p.y()))->index() = nVert++;
273 }
274 }
275 }
276
277 Info<< nVert << " vertices inserted" << endl;
278
279 if (meshControls().objOutput())
280 {
281 // Checking validity of triangulation
282 assert(is_valid());
283
284 writeTriangles("initial_triangles.obj", true);
285 writeFaces("initial_faces.obj", true);
286 }
287}
288
289
291{
292 startOfSurfacePointPairs_ = number_of_vertices();
293
294 if (meshControls().insertSurfaceNearestPointPairs())
295 {
296 insertSurfaceNearestPointPairs();
297 }
298
299 write("nearest");
300
301 // Insertion of point-pairs for near-points may cause protrusions
302 // so insertBoundaryConformPointPairs must be executed last
303 if (meshControls().insertSurfaceNearPointPairs())
304 {
305 insertSurfaceNearPointPairs();
306 }
307
308 startOfBoundaryConformPointPairs_ = number_of_vertices();
309}
310
311
313{
314 if (!meshControls().insertSurfaceNearestPointPairs())
315 {
316 markNearBoundaryPoints();
317 }
318
319 // Mark all the faces as SAVE_CHANGED
320 for
321 (
322 Triangulation::Finite_faces_iterator fit = finite_faces_begin();
323 fit != finite_faces_end();
324 fit++
325 )
326 {
327 fit->faceIndex() = Fb::SAVE_CHANGED;
328 }
329
330 for (label iter=1; iter<=meshControls().maxBoundaryConformingIter(); iter++)
331 {
332 label nIntersections = insertBoundaryConformPointPairs
333 (
334 "surfaceIntersections_" + Foam::name(iter) + ".obj"
335 );
336
337 if (nIntersections == 0)
338 {
339 break;
340 }
341 else
342 {
343 Info<< "BC iteration " << iter << ": "
344 << nIntersections << " point-pairs inserted" << endl;
345 }
346
347 // Any faces changed by insertBoundaryConformPointPairs will now
348 // be marked CHANGED, mark those as SAVE_CHANGED and those that
349 // remained SAVE_CHANGED as UNCHANGED
350 for
351 (
352 Triangulation::Finite_faces_iterator fit = finite_faces_begin();
353 fit != finite_faces_end();
354 fit++
355 )
356 {
357 if (fit->faceIndex() == Fb::SAVE_CHANGED)
358 {
359 fit->faceIndex() = Fb::UNCHANGED;
360 }
361 else if (fit->faceIndex() == Fb::CHANGED)
362 {
363 fit->faceIndex() = Fb::SAVE_CHANGED;
364 }
365 }
366 }
367
368 Info<< nl;
369
370 write("boundary");
371}
372
373
375{
376 for
377 (
378 Triangulation::Finite_vertices_iterator vit = finite_vertices_begin();
379 vit != finite_vertices_end();
380 ++vit
381 )
382 {
383 if (vit->index() >= startOfSurfacePointPairs_)
384 {
385 remove(vit);
386 }
387 }
388}
389
390
392{
393 const scalar relaxation = relaxationModel_->relaxation();
394
395 Info<< "Relaxation = " << relaxation << endl;
396
397 Field<point2D> dualVertices(number_of_faces());
398
399 label dualVerti = 0;
400
401 // Find the dual point of each tetrahedron and assign it an index.
402 for
403 (
404 Triangulation::Finite_faces_iterator fit = finite_faces_begin();
405 fit != finite_faces_end();
406 ++fit
407 )
408 {
409 fit->faceIndex() = -1;
410
411 if
412 (
413 fit->vertex(0)->internalOrBoundaryPoint()
414 || fit->vertex(1)->internalOrBoundaryPoint()
415 || fit->vertex(2)->internalOrBoundaryPoint()
416 )
417 {
418 fit->faceIndex() = dualVerti;
419
420 dualVertices[dualVerti] = toPoint2D(circumcenter(fit));
421
422 dualVerti++;
423 }
424 }
425
426 dualVertices.setSize(dualVerti);
427
428 Field<vector2D> displacementAccumulator
429 (
430 startOfSurfacePointPairs_,
431 vector2D::zero
432 );
433
434 // Calculate target size and alignment for vertices
435 scalarField sizes
436 (
437 number_of_vertices(),
438 meshControls().minCellSize()
439 );
440
441 Field<vector2D> alignments
442 (
443 number_of_vertices(),
444 vector2D(1, 0)
445 );
446
447 for
448 (
449 Triangulation::Finite_vertices_iterator vit = finite_vertices_begin();
450 vit != finite_vertices_end();
451 ++vit
452 )
453 {
454 if (vit->internalOrBoundaryPoint())
455 {
456 point2D vert = toPoint2D(vit->point());
457
458 // alignment and size determination
459 pointIndexHit pHit;
460 label hitSurface = -1;
461
462 qSurf_.findSurfaceNearest
463 (
464 toPoint3D(vert),
465 meshControls().span2(),
466 pHit,
467 hitSurface
468 );
469
470 if (pHit.hit())
471 {
472 vectorField norm(1);
473 allGeometry_[hitSurface].getNormal
474 (
475 List<pointIndexHit>(1, pHit),
476 norm
477 );
478
479 alignments[vit->index()] = toPoint2D(norm[0]);
480
481 sizes[vit->index()] =
482 cellSizeControl_.cellSize
483 (
484 toPoint3D(vit->point())
485 );
486 }
487 }
488 }
489
490 // Info<< "Calculated alignments" << endl;
491
492 scalar cosAlignmentAcceptanceAngle = 0.68;
493
494 // Upper and lower edge length ratios for weight
495 scalar u = 1.0;
496 scalar l = 0.7;
497
498 bitSet pointToBeRetained(startOfSurfacePointPairs_, true);
499
500 std::list<Point> pointsToInsert;
501
502 for
503 (
504 Triangulation::Finite_edges_iterator eit = finite_edges_begin();
505 eit != finite_edges_end();
506 eit++
507 )
508 {
509 Vertex_handle vA = eit->first->vertex(cw(eit->second));
510 Vertex_handle vB = eit->first->vertex(ccw(eit->second));
511
512 if (!vA->internalOrBoundaryPoint() || !vB->internalOrBoundaryPoint())
513 {
514 continue;
515 }
516
517 const point2D& dualV1 = dualVertices[eit->first->faceIndex()];
518 const point2D& dualV2 =
519 dualVertices[eit->first->neighbor(eit->second)->faceIndex()];
520
521 scalar dualEdgeLength = mag(dualV1 - dualV2);
522
523 point2D dVA = toPoint2D(vA->point());
524 point2D dVB = toPoint2D(vB->point());
525
526 Field<vector2D> alignmentDirsA(2);
527
528 alignmentDirsA[0] = alignments[vA->index()];
529 alignmentDirsA[1] = vector2D
530 (
531 -alignmentDirsA[0].y(),
532 alignmentDirsA[0].x()
533 );
534
535 Field<vector2D> alignmentDirsB(2);
536
537 alignmentDirsB[0] = alignments[vB->index()];
538 alignmentDirsB[1] = vector2D
539 (
540 -alignmentDirsB[0].y(),
541 alignmentDirsB[0].x()
542 );
543
544 Field<vector2D> alignmentDirs(alignmentDirsA);
545
546 forAll(alignmentDirsA, aA)
547 {
548 const vector2D& a(alignmentDirsA[aA]);
549
550 scalar maxDotProduct = 0.0;
551
552 forAll(alignmentDirsB, aB)
553 {
554 const vector2D& b(alignmentDirsB[aB]);
555
556 scalar dotProduct = a & b;
557
558 if (mag(dotProduct) > maxDotProduct)
559 {
560 maxDotProduct = mag(dotProduct);
561
562 alignmentDirs[aA] = a + sign(dotProduct)*b;
563
564 alignmentDirs[aA].normalise();
565 }
566 }
567 }
568
569 vector2D rAB = dVA - dVB;
570
571 scalar rABMag = mag(rAB);
572
573 forAll(alignmentDirs, aD)
574 {
575 vector2D& alignmentDir = alignmentDirs[aD];
576
577 if ((rAB & alignmentDir) < 0)
578 {
579 // swap the direction of the alignment so that has the
580 // same sense as rAB
581 alignmentDir *= -1;
582 }
583
584 scalar alignmentDotProd = ((rAB/rABMag) & alignmentDir);
585
586 if (alignmentDotProd > cosAlignmentAcceptanceAngle)
587 {
588 scalar targetFaceSize =
589 0.5*(sizes[vA->index()] + sizes[vB->index()]);
590
591 // Test for changing aspect ratio on second alignment (first
592 // alignment is nearest surface normal)
593 // if (aD == 1)
594 // {
595 // targetFaceSize *= 2.0;
596 // }
597
598 alignmentDir *= 0.5*targetFaceSize;
599
600 vector2D delta = alignmentDir - 0.5*rAB;
601
602 if (dualEdgeLength < 0.7*targetFaceSize)
603 {
604 delta *= 0;
605 }
606 else if (dualEdgeLength < targetFaceSize)
607 {
608 delta *=
609 (
610 dualEdgeLength
611 /(targetFaceSize*(u - l))
612 - 1/((u/l) - 1)
613 );
614 }
615
616 if
617 (
618 vA->internalPoint()
619 && vB->internalPoint()
620 && rABMag > 1.75*targetFaceSize
621 && dualEdgeLength > 0.05*targetFaceSize
622 && alignmentDotProd > 0.93
623 )
624 {
625 // Point insertion
626 pointsToInsert.push_back(toPoint(0.5*(dVA + dVB)));
627 }
628 else if
629 (
630 (vA->internalPoint() || vB->internalPoint())
631 && rABMag < 0.65*targetFaceSize
632 )
633 {
634 // Point removal
635
636 // Only insert a point at the midpoint of the short edge
637 // if neither attached point has already been identified
638 // to be removed.
639 if
640 (
641 pointToBeRetained.test(vA->index())
642 && pointToBeRetained.test(vB->index())
643 )
644 {
645 pointsToInsert.push_back(toPoint(0.5*(dVA + dVB)));
646 }
647
648 if (vA->internalPoint())
649 {
650 pointToBeRetained.unset(vA->index());
651 }
652
653 if (vB->internalPoint())
654 {
655 pointToBeRetained.unset(vB->index());
656 }
657 }
658 else
659 {
660 if (vA->internalPoint())
661 {
662 displacementAccumulator[vA->index()] += delta;
663 }
664
665 if (vB->internalPoint())
666 {
667 displacementAccumulator[vB->index()] += -delta;
668 }
669 }
670 }
671 }
672 }
673
674 vector2D totalDisp = sum(displacementAccumulator);
675 scalar totalDist = sum(mag(displacementAccumulator));
676
677 // Relax the calculated displacement
678 displacementAccumulator *= relaxation;
679
680 label numberOfNewPoints = pointsToInsert.size();
681
682 for
683 (
684 Triangulation::Finite_vertices_iterator vit = finite_vertices_begin();
685 vit != finite_vertices_end();
686 ++vit
687 )
688 {
689 if (vit->internalPoint())
690 {
691 if (pointToBeRetained.test(vit->index()))
692 {
693 pointsToInsert.push_front
694 (
695 toPoint
696 (
697 toPoint2D(vit->point())
698 + displacementAccumulator[vit->index()]
699 )
700 );
701 }
702 }
703 }
704
705 // Clear the triangulation and reinsert the bounding box and feature points.
706 // This is faster than removing and moving points.
707 this->clear();
708
709 insertBoundingBox();
710
711 reinsertFeaturePoints();
712
713 startOfInternalPoints_ = number_of_vertices();
714
715 label nVert = startOfInternalPoints_;
716
717 Info<< "Inserting " << numberOfNewPoints << " new points" << endl;
718
719 // Use the range insert as it is faster than individually inserting points.
720 insert(pointsToInsert.begin(), pointsToInsert.end());
721
722 for
723 (
724 Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
725 vit != finite_vertices_end();
726 ++vit
727 )
728 {
729 if
730 (
731 vit->type() == Vb::INTERNAL_POINT
732 && vit->index() == Vb::INTERNAL_POINT
733 )
734 {
735 vit->index() = nVert++;
736 }
737 }
738
739 Info<< " Total displacement = " << totalDisp << nl
740 << " Total distance = " << totalDist << nl
741 << " Points added = " << pointsToInsert.size()
742 << endl;
743
744 write("internal");
745
746 insertSurfacePointPairs();
747
748 boundaryConform();
749
750
751 // Old Method
752 /*
753 for
754 (
755 Triangulation::Finite_vertices_iterator vit = finite_vertices_begin();
756 vit != finite_vertices_end();
757 ++vit
758 )
759 {
760 if (vit->internalPoint())
761 {
762 // Current dual-cell defining vertex ("centre")
763 point2DFromPoint defVert0 = toPoint2D(vit->point());
764
765 Triangulation::Edge_circulator ec = incident_edges(vit);
766 Triangulation::Edge_circulator ecStart = ec;
767
768 // Circulate around the edges to find the first which is not
769 // infinite
770 do
771 {
772 if (!is_infinite(ec)) break;
773 } while (++ec != ecStart);
774
775 // Store the start-end of the first non-infinite edge
776 point2D de0 = toPoint2D(circumcenter(ec->first));
777
778 // Keep track of the maximum edge length^2
779 scalar maxEdgeLen2 = 0.0;
780
781 // Keep track of the index of the longest edge
782 label edgecd0i = -1;
783
784 // Edge counter
785 label edgei = 0;
786
787 do
788 {
789 if (!is_infinite(ec))
790 {
791 // Get the end of the current edge
792 point2D de1 = toPoint2D
793 (
794 circumcenter(ec->first->neighbor(ec->second))
795 );
796
797 // Store the current edge vector
798 edges[edgei] = de1 - de0;
799
800 // Store the edge mid-point in the vertices array
801 vertices[edgei] = 0.5*(de1 + de0);
802
803 // Move the current edge end into the edge start for the
804 // next iteration
805 de0 = de1;
806
807 // Keep track of the longest edge
808
809 scalar edgeLen2 = magSqr(edges[edgei]);
810
811 if (edgeLen2 > maxEdgeLen2)
812 {
813 maxEdgeLen2 = edgeLen2;
814 edgecd0i = edgei;
815 }
816
817 edgei++;
818 }
819 } while (++ec != ecStart);
820
821 // Initialise cd0 such that the mesh will align
822 // in in the x-y directions
823 vector2D cd0(1, 0);
824
825 if (meshControls().relaxOrientation())
826 {
827 // Get the longest edge from the array and use as the primary
828 // direction of the coordinate system of the "square" cell
829 cd0 = edges[edgecd0i];
830 }
831
832 if (meshControls().nearWallAlignedDist() > 0)
833 {
834 pointIndexHit pHit = qSurf_.tree().findNearest
835 (
836 toPoint3D(defVert0),
837 meshControls().nearWallAlignedDist2()
838 );
839
840 if (pHit.hit())
841 {
842 cd0 = toPoint2D(faceNormals[pHit.index()]);
843 }
844 }
845
846 // Rotate by 45deg needed to create an averaging procedure which
847 // encourages the cells to be square
848 cd0 = vector2D(cd0.x() + cd0.y(), cd0.y() - cd0.x());
849
850 // Normalise the primary coordinate direction
851 cd0.normalise();
852
853 // Calculate the orthogonal coordinate direction
854 vector2D cd1(-cd0.y(), cd0.x());
855
856
857 // Restart the circulator
858 ec = ecStart;
859
860 // ... and the counter
861 edgei = 0;
862
863 // Initialise the displacement for the centre and sum-weights
864 vector2D disp = Zero;
865 scalar sumw = 0;
866
867 do
868 {
869 if (!is_infinite(ec))
870 {
871 // Pick up the current edge
872 const vector2D& ei = edges[edgei];
873
874 // Calculate the centre to edge-centre vector
875 vector2D deltai = vertices[edgei] - defVert0;
876
877 // Set the weight for this edge contribution
878 scalar w = 1;
879
880 if (meshControls().squares())
881 {
882 w = magSqr(deltai.x()*ei.y() - deltai.y()*ei.x());
883 // alternative weights
884 //w = mag(deltai.x()*ei.y() - deltai.y()*ei.x());
885 //w = magSqr(ei)*mag(deltai);
886
887 // Use the following for an ~square mesh
888 // Find the coordinate contributions for this edge delta
889 scalar cd0deltai = cd0 & deltai;
890 scalar cd1deltai = cd1 & deltai;
891
892 // Create a "square" displacement
893 if (mag(cd0deltai) > mag(cd1deltai))
894 {
895 disp += (w*cd0deltai)*cd0;
896 }
897 else
898 {
899 disp += (w*cd1deltai)*cd1;
900 }
901 }
902 else
903 {
904 // Use this for a hexagon/pentagon mesh
905 disp += w*deltai;
906 }
907
908 // Sum the weights
909 sumw += w;
910 }
911 else
912 {
913 FatalErrorInFunction
914 << "Infinite triangle found in internal mesh"
915 << exit(FatalError);
916 }
917
918 edgei++;
919
920 } while (++ec != ecStart);
921
922 // Calculate the average displacement
923 disp /= sumw;
924 totalDisp += disp;
925 totalDist += mag(disp);
926
927 // Move the point by a fraction of the average displacement
928 movePoint(vit, defVert0 + relaxation*disp);
929 }
930 }
931
932 Info << "\nTotal displacement = " << totalDisp
933 << " total distance = " << totalDist << endl;
934 */
935}
936
937/*
938void Foam::CV2D::moveInternalPoints(const point2DField& newPoints)
939{
940 label pointi = 0;
941
942 for
943 (
944 Triangulation::Finite_vertices_iterator vit = finite_vertices_begin();
945 vit != finite_vertices_end();
946 ++vit
947 )
948 {
949 if (vit->internalPoint())
950 {
951 movePoint(vit, newPoints[pointi++]);
952 }
953 }
954}
955*/
956
957void Foam::CV2D::write() const
958{
959 if (meshControls().objOutput())
960 {
961 writeFaces("allFaces.obj", false);
962 writeFaces("faces.obj", true);
963 writeTriangles("allTriangles.obj", false);
964 writeTriangles("triangles.obj", true);
965 writePatch("patch.pch");
966 }
967}
968
969
970void Foam::CV2D::write(const word& stage) const
971{
972 if (meshControls().objOutput())
973 {
974 Foam::mkDir(stage + "Faces");
975 Foam::mkDir(stage + "Triangles");
976
977 writeFaces
978 (
979 stage
980 + "Faces/allFaces_"
981 + runTime_.timeName()
982 + ".obj",
983 false
984 );
985
986 writeTriangles
987 (
988 stage
989 + "Triangles/allTriangles_"
990 + runTime_.timeName()
991 + ".obj",
992 false
993 );
994 }
995}
996
997
998// ************************************************************************* //
CGAL::Point_3< K > Point
CGAL::Delaunay_triangulation_3< K, Tds, CompactLocator > Delaunay
scalar y
scalar delta
label n
Conformal-Voronoi 2D automatic mesher with grid or read initial points and point position relaxation ...
Definition: CV2D.H:149
void write() const
const cv2DControls & meshControls() const
Definition: CV2DI.H:119
~CV2D()
Destructor.
void insertSurfacePointPairs()
Insert all surface point-pairs from.
void boundaryConform()
Insert point-pairs where there are protrusions into.
void insertGrid()
Create the initial mesh as a regular grid of points.
void removeSurfacePointPairs()
Remove the point-pairs introduced by insertSurfacePointPairs.
void newPoints()
Move the internal points to the given new locations and update.
void insertPoints(const point2DField &points, const scalar nearness)
Create the initial mesh from the given internal points.
Triangulation::Finite_vertices_iterator Finite_vertices_iterator
Definition: DelaunayMesh.H:73
void normalise()
Normalise the field inplace. See notes in Field.
Type sample01()
Return a sample whose components lie in the range [0,1].
scalar span() const
Return the span.
constant condensation/saturation model.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
patchWriters clear()
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
Namespace for OpenFOAM.
Vector2D< scalar > vector2D
A 2D vector of scalars obtained from the generic Vector2D.
Definition: vector2D.H:51
dimensionedScalar sign(const dimensionedScalar &ds)
PointFrompoint toPoint(const Foam::point &p)
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:515
dimensionedScalar y0(const dimensionedScalar &ds)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
vector2DField point2DField
point2DField is a vector2DField.
messageStream Info
Information stream (stdout output on master, null elsewhere)
vector point
Point is a vector.
Definition: point.H:43
vector2D point2D
Point2D is a vector.
Definition: point2D.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
Field< vector > vectorField
Specialisation of Field<T> for vector.
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
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
messageStream Warning
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
runTime write()
srcOptions insert("case", fileName(rootDirSource/caseDirSource))
labelList f(nPoints)
volScalarField & b
Definition: createFields.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
3D tensor transformation operations.
System unsigned integer.
Random rndGen
Definition: createFields.H:23