cutFaceAdvect.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) 2016-2017 DHI
9 Copyright (C) 2018-2019 Johan Roenby
10 Copyright (C) 2019-2020 DLR
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "cutFaceAdvect.H"
31#include "OFstream.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
36(
37 const fvMesh& mesh,
39)
40:
42 mesh_(mesh),
43 //alpha1_(alpha1),
44 subFaceCentre_(Zero),
45 subFaceArea_(Zero),
46 subFacePoints_(10),
47 surfacePoints_(4),
48 pointStatus_(10),
49 weight_(10),
50 pTimes_(10),
51 faceStatus_(-1)
52{
54}
55
56
57// * * * * * * * * * * * Public Member Functions * * * * * * * * * * * * * //
58
60(
61 const label faceI,
62 const vector& normal,
63 const vector& base
64)
65{
66 clearStorage();
67
68 const face& f = mesh_.faces()[faceI];
69
70 label inLiquid = 0;
71 label firstFullySubmergedPoint = -1;
72
73 // Loop face and calculate pointStatus
74 forAll(f, i)
75 {
76 scalar value = (mesh_.points()[f[i]] - base) & normal;
77 if (mag(value) < 10 * SMALL)
78 {
79 value = 0;
80 }
81 pointStatus_.append(value);
82 if (pointStatus_[i] > 10 * SMALL)
83 {
84 inLiquid++;
85 if (firstFullySubmergedPoint == -1)
86 {
87 firstFullySubmergedPoint = i;
88 }
89 }
90 }
91
92 if (inLiquid == f.size()) // fluid face
93 {
94 faceStatus_ = -1;
95 subFaceCentre_ = mesh_.faceCentres()[faceI];
96 subFaceArea_ = mesh_.faceAreas()[faceI];
97 return faceStatus_;
98 }
99 else if (inLiquid == 0) // gas face
100 {
101 faceStatus_ = 1;
102 subFaceCentre_ = Zero;
103 subFaceArea_ = Zero;
104 return faceStatus_;
105 }
106
108 (
109 faceI,
110 pointStatus_,
111 firstFullySubmergedPoint,
112 subFacePoints_,
113 surfacePoints_,
114 faceStatus_,
115 subFaceCentre_,
116 subFaceArea_
117 );
118
119 return faceStatus_;
120}
121
122
124(
125 const face& f,
126 const pointField& points,
127 const scalarField& val,
128 const scalar cutValue
129)
130{
131 clearStorage();
132
133 label inLiquid = 0;
134 label firstFullySubmergedPoint = -1;
135 scalarList pointStatus(f.size());
136
137 // Loop face and calculate pointStatus
138 forAll(f, i)
139 {
140 pointStatus[i] = val[f[i]] - cutValue;
141 if (mag(pointStatus[i]) < 10 * SMALL)
142 {
143 pointStatus[i] = 0;
144 }
145 if (pointStatus[i] > 10 * SMALL)
146 {
147 ++inLiquid;
148 if (firstFullySubmergedPoint == -1)
149 {
150 firstFullySubmergedPoint = i;
151 }
152 }
153 }
154
155 if (inLiquid == f.size()) // fluid face
156 {
157 faceStatus_ = -1;
158 subFaceCentre_ = f.centre(points);
159 subFaceArea_ = f.areaNormal(points);
160 return faceStatus_;
161 }
162 else if (inLiquid == 0) // gas face
163 {
164 faceStatus_ = 1;
165 subFaceCentre_ = Zero;
166 subFaceArea_ = Zero;
167 return faceStatus_;
168 }
169
171 (
172 f,
173 points,
174 pointStatus,
175 firstFullySubmergedPoint,
176 subFacePoints_,
177 surfacePoints_,
178 faceStatus_,
179 subFaceCentre_,
180 subFaceArea_
181 );
182
183 return faceStatus_;
184}
185
186
188(
189 const label faceI,
190 const vector& x0,
191 const vector& n0, // has to has the length of 1
192 const scalar Un0,
193 const scalar dt,
194 const scalar phi,
195 const scalar magSf
196)
197{
198 clearStorage();
199
200/* Temporarily taken out
201 // Treating rare cases where isoface normal is not calculated properly
202 if (mag(n0) < 0.5)
203 {
204 scalar alphaf = 0.0;
205 scalar waterInUpwindCell = 0.0;
206
207 if (phi > 0 || !mesh_.isInternalFace(faceI))
208 {
209 const label upwindCell = mesh_.faceOwner()[faceI];
210 alphaf = alpha1_[upwindCell];
211 waterInUpwindCell = alphaf * mesh_.V()[upwindCell];
212 }
213 else
214 {
215 const label upwindCell = mesh_.faceNeighbour()[faceI];
216 alphaf = alpha1_[upwindCell];
217 waterInUpwindCell = alphaf * mesh_.V()[upwindCell];
218 }
219
220 return min(alphaf * phi * dt, waterInUpwindCell);
221 }*/
222
223 // Find sorted list of times where the isoFace will arrive at face points
224 // given initial position x0 and velocity Un0*n0
225
226 // Get points for this face
227 const face& f = mesh_.faces()[faceI];
228 const label nPoints = f.size();
229
230 if (mag(Un0) > 1e-12) // Note: tolerances
231 {
232 // Here we estimate time of arrival to the face points from their normal
233 // distance to the initial surface and the surface normal velocity
234
235 for (const scalar fi : f)
236 {
237 scalar value = ((mesh_.points()[fi] - x0) & n0) / Un0;
238 if (mag(value) < 10 * SMALL)
239 {
240 value = 0;
241 }
242 pTimes_.append(value);
243 }
244
245 scalar dVf = 0;
246
247 // Check if pTimes changes direction more than twice when looping face
248 label nShifts = 0;
249 forAll(pTimes_, pi)
250 {
251 const label oldEdgeSign =
252 sign(pTimes_[(pi + 1) % nPoints] - pTimes_[pi]);
253 const label newEdgeSign =
254 sign(pTimes_[(pi + 2) % nPoints] - pTimes_[(pi + 1) % nPoints]);
255
256 if (newEdgeSign != oldEdgeSign)
257 {
258 nShifts++;
259 }
260 }
261
262 if (nShifts == 2 || nShifts == 0)
263 {
264 dVf = phi / magSf * timeIntegratedArea(faceI, dt, magSf, Un0);
265 }
266 else if (nShifts > 2) // triangle decompose the non planar face
267 {
268 const pointField fPts(f.points(mesh_.points()));
269 pointField fPts_tri(3);
270 scalarField pTimes_tri(3);
271 fPts_tri[0] = mesh_.faceCentres()[faceI];
272 pTimes_tri[0] = ((fPts_tri[0] - x0) & n0)/Un0;
273 for (label pi = 0; pi < nPoints; ++pi)
274 {
275 fPts_tri[1] = fPts[pi];
276 pTimes_tri[1] = pTimes_[pi];
277 fPts_tri[2] = fPts[(pi + 1) % nPoints];
278 pTimes_tri[2] = pTimes_[(pi + 1) % nPoints];
279 const scalar magSf_tri =
280 mag
281 (
282 0.5
283 *(fPts_tri[2] - fPts_tri[0])
284 ^(fPts_tri[1] - fPts_tri[0])
285 );
286
287 const scalar phi_tri = phi*magSf_tri/magSf;
288 dVf +=
289 phi_tri/magSf_tri
290 *timeIntegratedArea
291 (
292 fPts_tri,
293 pTimes_tri,
294 dt,
295 magSf_tri,
296 Un0
297 );
298 }
299 }
300 else
301 {
302 if (debug)
303 {
305 << "Warning: nShifts = " << nShifts << " on face "
306 << faceI << " with pTimes = " << pTimes_
307 << " owned by cell " << mesh_.faceOwner()[faceI]
308 << endl;
309 }
310 }
311
312 return dVf;
313 }
314 else
315 {
316 // Un0 is almost zero and isoFace is treated as stationary
317 calcSubFace(faceI, -n0, x0);
318 const scalar alphaf = mag(subFaceArea() / magSf);
319
320 if (debug)
321 {
323 << "Un0 is almost zero (" << Un0
324 << ") - calculating dVf on face " << faceI
325 << " using subFaceFraction giving alphaf = " << alphaf
326 << endl;
327 }
328
329 return phi * dt * alphaf;
330 }
331}
332
333
335(
336 const label faceI,
337 const scalarField& times,
338 const scalar Un0,
339 const scalar dt,
340 const scalar phi,
341 const scalar magSf
342)
343{
344 clearStorage();
345
346 label nPoints = times.size();
347
348 {
349 // Here we estimate time of arrival to the face points from their normal
350 // distance to the initial surface and the surface normal velocity
351
352 pTimes_.append(times);
353
354 scalar dVf = 0;
355
356 // Check if pTimes changes direction more than twice when looping face
357 label nShifts = 0;
358 forAll(pTimes_, pi)
359 {
360 const label oldEdgeSign =
361 sign(pTimes_[(pi + 1) % nPoints] - pTimes_[pi]);
362 const label newEdgeSign =
363 sign(pTimes_[(pi + 2) % nPoints] - pTimes_[(pi + 1) % nPoints]);
364
365 if (newEdgeSign != oldEdgeSign)
366 {
367 ++nShifts;
368 }
369 }
370
371 if (nShifts == 2)
372 {
373 dVf = phi/magSf*timeIntegratedArea(faceI, dt, magSf, Un0);
374 }
375 // not possible to decompose face
376 return dVf;
377 }
378}
379
380
382(
383 const pointField& fPts,
384 const scalarField& pTimes,
385 const scalar dt,
386 const scalar magSf,
387 const scalar Un0
388)
389{
390 // Initialise time integrated area returned by this function
391 scalar tIntArea = 0.0;
392
393 // Finding ordering of vertex points
394 const labelList order(Foam::sortedOrder(pTimes));
395 const scalar firstTime = pTimes[order.first()];
396 const scalar lastTime = pTimes[order.last()];
397
398 // Dealing with case where face is not cut by surface during time interval
399 // [0,dt] because face was already passed by surface
400 if (lastTime <= 0)
401 {
402 // If all face cuttings were in the past and cell is filling up (Un0>0)
403 // then face must be full during whole time interval
404 tIntArea = magSf * dt * pos0(Un0);
405 return tIntArea;
406 }
407
408 // Dealing with case where face is not cut by surface during time interval
409 // [0, dt] because dt is too small for surface to reach closest face point
410 if (firstTime >= dt)
411 {
412 // If all cuttings are in the future but non of them within [0,dt] then
413 // if cell is filling up (Un0 > 0) face must be empty during whole time
414 // interval
415 tIntArea = magSf * dt * (1 - pos0(Un0));
416 return tIntArea;
417 }
418
419 // If we reach this point in the code some part of the face will be swept
420 // during [tSmall, dt-tSmall]. However, it may be the case that there are no
421 // vertex times within the interval. This will happen sometimes for small
422 // time steps where both the initial and the final face-interface
423 // intersection line (FIIL) will be along the same two edges.
424
425 // Face-interface intersection line (FIIL) to be swept across face
426 DynamicList<point> FIIL(3);
427 // Submerged area at beginning of each sub time interval time
428 scalar initialArea = 0.0;
429 // Running time keeper variable for the integration process
430 scalar time = 0.0;
431
432 // Special treatment of first sub time interval
433 if (firstTime > 0)
434 {
435 // If firstTime > 0 the face is uncut in the time interval
436 // [0, firstTime] and hence fully submerged in fluid A or B.
437 // If Un0 > 0 cell is filling up and it must initially be empty.
438 // If Un0 < 0 cell must initially be full(y immersed in fluid A).
439 time = firstTime;
440 initialArea = magSf * (1.0 - pos0(Un0));
441 tIntArea = initialArea * time;
442 cutPoints(fPts, pTimes, time, FIIL);
443 }
444 else
445 {
446 // If firstTime <= 0 then face is initially cut and we must
447 // calculate the initial submerged area and FIIL:
448 time = 0.0;
449 // Note: calcSubFace assumes well-defined 2-point FIIL!!!!
450 // calcSubFace(fPts, -sign(Un0)*pTimes, time);
451 // calcSubFace(fPts, -sign(Un0)*pTimes, time)
452 calcSubFace(face(identity(pTimes.size())), fPts, pTimes, time);
453 initialArea = mag(subFaceArea());
454 cutPoints(fPts, pTimes, time, FIIL);
455 }
456
457 // Making sorted array of all vertex times that are between max(0,firstTime)
458 // and dt and further than tSmall from the previous time.
459 DynamicList<scalar> sortedTimes(pTimes.size());
460 {
461 scalar prevTime = time;
462 const scalar tSmall = max(1e-6*dt, 10*SMALL);
463
464 for (const scalar timeI : order)
465 {
466 if (timeI > prevTime + tSmall && timeI <= dt)
467 {
468 sortedTimes.append(timeI);
469 prevTime = timeI;
470 }
471 }
472 }
473
474 // Sweeping all quadrilaterals corresponding to the intervals defined above
475 for (const scalar newTime : sortedTimes)
476 {
477 // New face-interface intersection line
478 DynamicList<point> newFIIL(3);
479 cutPoints(fPts, pTimes, newTime, newFIIL);
480
481 // quadrilateral area coefficients
482 scalar alpha = 0, beta = 0;
483 quadAreaCoeffs(FIIL, newFIIL, alpha, beta);
484 // Integration of area(t) = A*t^2+B*t from t = 0 to 1
485 tIntArea += (newTime - time) *
486 (initialArea + sign(Un0) * (alpha/3.0 + 0.5*beta));
487 // Adding quad area to submerged area
488 initialArea += sign(Un0)*(alpha + beta);
489
490 FIIL = newFIIL;
491 time = newTime;
492 }
493
494 // Special treatment of last time interval
495 if (lastTime > dt)
496 {
497 // FIIL will end up cutting the face at dt
498 // New face-interface intersection line
499 DynamicList<point> newFIIL(3);
500 cutPoints(fPts, pTimes, dt, newFIIL);
501
502 // quadrilateral area coefficients
503 scalar alpha = 0, beta = 0;
504 quadAreaCoeffs(FIIL, newFIIL, alpha, beta);
505 // Integration of area(t) = A*t^2+B*t from t = 0 to 1
506 tIntArea += (dt - time) *
507 (initialArea + sign(Un0)*(alpha / 3.0 + 0.5 * beta));
508 }
509 else
510 {
511 // FIIL will leave the face at lastTime and face will be fully in fluid
512 // A or fluid B in the time interval from lastTime to dt.
513 tIntArea += magSf*(dt - lastTime)*pos0(Un0);
514 }
515
516 return tIntArea;
517}
518
519
520void Foam::cutFaceAdvect::isoFacesToFile
521(
522 const DynamicList<List<point>>& faces,
523 const word& filNam,
524 const word& filDir
525) const
526{
527 // Writing isofaces to vtk file for inspection in paraview
528
529 fileName outputFile(filDir/(filNam + ".vtk"));
530
531 mkDir(filDir);
532 Info<< "Writing file: " << outputFile << endl;
533
534 OFstream os(outputFile);
535 os << "# vtk DataFile Version 2.0" << nl
536 << filNam << nl
537 << "ASCII" << nl
538 << "DATASET POLYDATA" << nl;
539
540 label nPoints{0};
541 for (const List<point>& f : faces)
542 {
543 nPoints += f.size();
544 }
545
546 os << "POINTS " << nPoints << " float" << nl;
547 for (const List<point>& f : faces)
548 {
549 for (const point& p : f)
550 {
551 os << p.x() << ' ' << p.y() << ' ' << p.z() << nl;
552 }
553 }
554
555 os << "POLYGONS "
556 << faces.size() << ' ' << (nPoints + faces.size()) << nl;
557
558 label pointi = 0;
559 for (const List<point>& f : faces)
560 {
561 label endp = f.size();
562 os << endp;
563
564 endp += pointi;
565
566 while (pointi < endp)
567 {
568 os << ' ' << pointi;
569 ++pointi;
570 }
571 os << nl;
572 }
573}
574
575
577(
578 const label faceI,
579 const label sign,
580 const scalar cutValue
581)
582{
583 // clearStorage();
584 const face& f = mesh_.faces()[faceI];
585 label inLiquid = 0;
586 label firstFullySubmergedPoint = -1;
587
588 // Loop face and calculate pointStatus
589 forAll(f, i)
590 {
591 scalar value = (sign * pTimes_[i] - cutValue);
592
593 if (mag(value) < 10 * SMALL)
594 {
595 value = 0;
596 }
597 pointStatus_.append(value);
598 if (pointStatus_[i] > 10 * SMALL)
599 {
600 inLiquid++;
601 if (firstFullySubmergedPoint == -1)
602 {
603 firstFullySubmergedPoint = i;
604 }
605 }
606 }
607
608 if (inLiquid == f.size()) // fluid face
609 {
610 faceStatus_ = -1;
611 subFaceCentre_ = mesh_.faceCentres()[faceI];
612 subFaceArea_ = mesh_.faceAreas()[faceI];
613 return faceStatus_;
614 }
615 else if (inLiquid == 0) // gas face
616 {
617 faceStatus_ = 1;
618 subFaceCentre_ = Zero;
619 subFaceArea_ = Zero;
620 return faceStatus_;
621 }
622
624 (
625 faceI,
626 pointStatus_,
627 firstFullySubmergedPoint,
628 subFacePoints_,
629 surfacePoints_,
630 faceStatus_,
631 subFaceCentre_,
632 subFaceArea_
633 );
634
635 return faceStatus_;
636}
637
638
640(
641 const label faceI,
642 const scalar dt,
643 const scalar magSf,
644 const scalar Un0
645)
646{
647 // Initialise time integrated area returned by this function
648 scalar tIntArea = 0.0;
649
650 // Finding ordering of vertex points
651 const labelList order(Foam::sortedOrder(pTimes_));
652 const scalar firstTime = pTimes_[order.first()];
653 const scalar lastTime = pTimes_[order.last()];
654
655 // Dealing with case where face is not cut by surface during time interval
656 // [0,dt] because face was already passed by surface
657 if (lastTime <= 0)
658 {
659 // If all face cuttings were in the past and cell is filling up (Un0>0)
660 // then face must be full during whole time interval
661 tIntArea = magSf* dt * pos0(Un0);
662 return tIntArea;
663 }
664
665 // Dealing with case where face is not cut by surface during time interval
666 // [0, dt] because dt is too small for surface to reach closest face point
667 if (firstTime >= dt)
668 {
669 // If all cuttings are in the future but non of them within [0,dt] then
670 // if cell is filling up (Un0 > 0) face must be empty during whole time
671 // interval
672 tIntArea = magSf * dt * (1 - pos0(Un0));
673 return tIntArea;
674 }
675
676 // If we reach this point in the code some part of the face will be swept
677 // during [tSmall, dt-tSmall]. However, it may be the case that there are no
678 // vertex times within the interval. This will happen sometimes for small
679 // time steps where both the initial and the final face-interface
680 // intersection line (FIIL) will be along the same two edges.
681
682 // Face-interface intersection line (FIIL) to be swept across face
683 DynamicList<point> FIIL(3);
684 // Submerged area at beginning of each sub time interval time
685 scalar initialArea = 0.0;
686 // Running time keeper variable for the integration process
687 scalar time = 0.0;
688
689 // Special treatment of first sub time interval
690 if (firstTime > 0)
691 {
692 // If firstTime > 0 the face is uncut in the time interval
693 // [0, firstTime] and hence fully submerged in fluid A or B.
694 // If Un0 > 0 cell is filling up and it must initially be empty.
695 // If Un0 < 0 cell must initially be full(y immersed in fluid A).
696 time = firstTime;
697 initialArea = magSf * (1.0 - pos0(Un0));
698 tIntArea = initialArea * time;
699 cutPoints(faceI, time, FIIL);
700 }
701 else
702 {
703 // If firstTime <= 0 then face is initially cut and we must
704 // calculate the initial submerged area and FIIL:
705 time = 0.0;
706 // Note: calcSubFace assumes well-defined 2-point FIIL!!!!
707 calcSubFace(faceI, -sign(Un0), time);
708 initialArea = mag(subFaceArea());
709 cutPoints(faceI, time, FIIL);
710 }
711
712 // Making sorted array of all vertex times that are between max(0,firstTime)
713 // and dt and further than tSmall from the previous time.
714 DynamicList<scalar> sortedTimes(pTimes_.size());
715 {
716 scalar prevTime = time;
717 const scalar tSmall = max(1e-6*dt, 10*SMALL);
718 for (const label oI : order)
719 {
720 const scalar timeI = pTimes_[oI];
721 if (timeI > prevTime + tSmall && timeI <= dt)
722 {
723 sortedTimes.append(timeI);
724 prevTime = timeI;
725 }
726 }
727 }
728
729 // Sweeping all quadrilaterals corresponding to the intervals defined above
730 for (const scalar newTime : sortedTimes)
731 {
732 // New face-interface intersection line
733 DynamicList<point> newFIIL(3);
734 cutPoints(faceI, newTime, newFIIL);
735
736 // quadrilateral area coefficients
737 scalar alpha = 0, beta = 0;
738
739 quadAreaCoeffs(FIIL, newFIIL, alpha, beta);
740 // Integration of area(t) = A*t^2+B*t from t = 0 to 1
741 tIntArea +=
742 (newTime - time)
743 * (initialArea + sign(Un0)
744 * (alpha / 3.0 + 0.5 * beta));
745 // Adding quad area to submerged area
746 initialArea += sign(Un0) * (alpha + beta);
747
748 FIIL = newFIIL;
749 time = newTime;
750 }
751
752 // Special treatment of last time interval
753 if (lastTime > dt)
754 {
755 // FIIL will end up cutting the face at dt
756 // New face-interface intersection line
757 DynamicList<point> newFIIL(3);
758 cutPoints(faceI, dt, newFIIL);
759
760 // quadrilateral area coefficients
761 scalar alpha = 0, beta = 0;
762 quadAreaCoeffs(FIIL, newFIIL, alpha, beta);
763 // Integration of area(t) = A*t^2+B*t from t = 0 to 1
764 tIntArea +=
765 (dt - time)
766 * (initialArea + sign(Un0) * (alpha / 3.0 + 0.5 * beta));
767 }
768 else
769 {
770 // FIIL will leave the face at lastTime and face will be fully in fluid
771 // A or fluid B in the time interval from lastTime to dt.
772 tIntArea += magSf * (dt - lastTime) * pos0(Un0);
773 }
774
775 return tIntArea;
776}
777
778
780(
781 const DynamicList<point>& pf0,
782 const DynamicList<point>& pf1,
783 scalar& alpha,
784 scalar& beta
785) const
786{
787 // Number of points in provided face-interface intersection lines
788 const label np0 = pf0.size();
789 const label np1 = pf1.size();
790
791 // quad area coeffs such that area(t) = alpha*t^2 + beta*t.
792 // With time interval normalised, we have full quadArea = alpha + beta
793 // and time integrated quad area = alpha/3 + beta/2;
794 alpha = 0.0;
795 beta = 0.0;
796
797 if (np0 && np1)
798 {
799 // Initialising quadrilateral vertices A, B, C and D
800 vector A(pf0[0]);
801 vector C(pf1[0]);
802 vector B(pf0[0]);
803 vector D(pf1[0]);
804
805 if (np0 == 2)
806 {
807 B = pf0[1];
808 }
809 else if (np0 > 2)
810 {
811 WarningInFunction << "Vertex face was cut at pf0 = " << pf0 << endl;
812 }
813
814 if (np1 == 2)
815 {
816 D = pf1[1];
817 }
818 else if (np1 > 2)
819 {
820 WarningInFunction << "Vertex face was cut at pf1 = " << pf1 << endl;
821 }
822
823 // Swapping pf1 points if pf0 and pf1 point in same general direction
824 // (because we want a quadrilateral ABCD where pf0 = AB and pf1 = CD)
825 if (((B - A) & (D - C)) > 0)
826 {
827 vector tmp = D;
828 D = C;
829 C = tmp;
830 }
831
832 // Defining local coordinates (xhat, yhat) for area integration of swept
833 // quadrilateral ABCD such that A = (0,0), B = (Bx,0), C = (Cx,Cy) and
834 // D = (Dx,Dy) with Cy = 0 and Dy > 0.
835
836 const scalar Bx = mag(B - A);
837
838 vector xhat(Zero);
839 if (Bx > 10 * SMALL)
840 {
841 // If |AB| > 0 ABCD we use AB to define xhat
842 xhat = (B - A) / mag(B - A);
843 }
844 else if (mag(C - D) > 10 * SMALL)
845 {
846 // If |AB| ~ 0 ABCD is a triangle ACD and we use CD for xhat
847 xhat = (C - D) / mag(C - D);
848 }
849 else
850 {
851 return;
852 }
853
854 // Defining vertical axis in local coordinates
855 vector yhat = D - A;
856 yhat -= ((yhat & xhat) * xhat);
857
858 if (mag(yhat) > 10 * SMALL)
859 {
860 yhat /= mag(yhat);
861
862 const scalar Cx = (C - A) & xhat;
863 const scalar Cy = mag((C - A) & yhat);
864 const scalar Dx = (D - A) & xhat;
865 const scalar Dy = mag((D - A) & yhat);
866
867 // area = ((Cx - Bx)*Dy - Dx*Cy)/6.0 + 0.25*Bx*(Dy + Cy);
868 alpha = 0.5 * ((Cx - Bx) * Dy - Dx * Cy);
869 beta = 0.5 * Bx * (Dy + Cy);
870 }
871 }
872 else
873 {
875 << "Vertex face was cut at " << pf0 << " and at "
876 << pf1 << endl;
877 }
878}
879
880
882(
883 const label faceI,
884 const scalar f0,
885 DynamicList<point>& cutPoints
886)
887{
888 const face& f = mesh_.faces()[faceI];
889 const label nPoints = f.size();
890 scalar f1(pTimes_[0]);
891
892 // Snapping vertex value to f0 if very close (needed for 2D cases)
893 if (mag(f1 - f0) < 10 * SMALL)
894 {
895 f1 = f0;
896 }
897
898 forAll(f, pi)
899 {
900 label pi2 = (pi + 1) % nPoints;
901 scalar f2 = pTimes_[pi2];
902
903 // Snapping vertex value
904 if (mag(f2 - f0) < 10 * SMALL)
905 {
906 f2 = f0;
907 }
908
909 if ((f1 < f0 && f2 > f0) || (f1 > f0 && f2 < f0))
910 {
911 const scalar s = (f0 - f1) / (f2 - f1);
912 cutPoints.append
913 (
914 mesh_.points()[f[pi]]
915 + s*(mesh_.points()[f[pi2]] - mesh_.points()[f[pi]])
916 );
917 }
918 else if (f1 == f0)
919 {
920 cutPoints.append(mesh_.points()[f[pi]]);
921 }
922 f1 = f2;
923 }
924
925 if (cutPoints.size() > 2)
926 {
928 << "cutPoints = " << cutPoints
929 << " for pts = " << f.points(mesh_.points())
930 << ", f - f0 = " << f - f0 << " and f0 = " << f0
931 << endl;
932 }
933}
934
935
937(
938 const pointField& pts,
939 const scalarField& f,
940 const scalar f0,
941 DynamicList<point>& cutPoints
942)
943{
944 const label nPoints = pts.size();
945 scalar f1(f[0]);
946
947 // Snapping vertex value to f0 if very close (needed for 2D cases)
948 if (mag(f1 - f0) < 10 * SMALL)
949 {
950 f1 = f0;
951 }
952
953 forAll(pts, pi)
954 {
955 label pi2 = (pi + 1) % nPoints;
956 scalar f2 = f[pi2];
957
958 // Snapping vertex value
959 if (mag(f2 - f0) < 10 * SMALL)
960 {
961 f2 = f0;
962 }
963
964 if ((f1 < f0 && f2 > f0) || (f1 > f0 && f2 < f0))
965 {
966 const scalar s = (f0 - f1) / (f2 - f1);
967 cutPoints.append(pts[pi] + s * (pts[pi2] - pts[pi]));
968 }
969 else if (f1 == f0)
970 {
971 cutPoints.append(pts[pi]);
972 }
973 f1 = f2;
974 }
975
976 if (cutPoints.size() > 2)
977 {
979 << "cutPoints = " << cutPoints << " for pts = " << pts
980 << ", f - f0 = " << f - f0 << " and f0 = " << f0
981 << endl;
982 }
983}
984
985
987{
988 subFaceCentre_ = Zero;
989 subFaceArea_ = Zero;
990 subFacePoints_.clear();
991 surfacePoints_.clear();
992 pointStatus_.clear();
993 pTimes_.clear();
994 weight_.clear();
995 faceStatus_ = -1;
996}
997
998
999// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
surfaceScalarField & phi
const volScalarField & alpha1
Graphite solid properties.
Definition: C.H:53
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
Output to file stream, using an OSstream.
Definition: OFstream.H:57
T & first()
Return the first element of the list.
Definition: UListI.H:202
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
T & last()
Return the last element of the list.
Definition: UListI.H:216
Calculates the face fluxes.
Definition: cutFaceAdvect.H:69
label calcSubFace(const label faceI, const vector &normal, const vector &base)
Calculates cut centre and cut area (plicReconstruction)
Definition: cutFaceAdvect.C:60
void clearStorage()
Resets internal variables.
scalar timeIntegratedArea(const label faceI, const scalar dt, const scalar magSf, const scalar Un0)
Calculate time integrated area for a face.
scalar timeIntegratedFaceFlux(const label faceI, const vector &x0, const vector &n0, const scalar Un0, const scalar dt, const scalar phi, const scalar magSf)
Calculate time integrated flux for a face.
void quadAreaCoeffs(const DynamicPointList &pf0, const DynamicPointList &pf1, scalar &quadArea, scalar &intQuadArea) const
For face with vertices p calculate its area and integrated area.
Base class for cutting a face, faceI, of an fvMesh, mesh_, at its intersections.
Definition: cutFace.H:60
void calcSubFace(const label faceI, const scalarList &pointStatus, label firstFullySubmergedPoint, DynamicList< point > &subFacePoints, DynamicList< point > &surfacePoints, label &faceStatus, vector &subFaceCentre, vector &subFaceArea)
Definition: cutFace.C:39
const pointField & cutPoints() const
Points for combined set of faces.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
A class for handling file names.
Definition: fileName.H:76
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
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))
#define WarningInFunction
Report a warning using Foam::Warning.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
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
messageStream Info
Information stream (stdout output on master, null elsewhere)
vector point
Point is a vector.
Definition: point.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
volScalarField & alpha
const dimensionedScalar & D
volScalarField & e
Definition: createFields.H:11
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333