34 const MPLICcellStorage& cellInfo,
41 const scalar cellAlpha = cellInfo.cellAlpha();
43 const UIndirectList<scalar>& cellMagSfs = cellInfo.magSf();
58 if (
mag(cellInfo.cellAlphaMax() - cellInfo.cellAlphaMin()) < vSmall)
64 findPointAlphaBounds(cellInfo, tetDecom);
67 if (
mag(pCubicAlphas_.a() - pCubicAlphas_.d()) < vSmall)
73 calcPointAlphaInterior(cellInfo, tetDecom);
76 const FixedList<scalar, 4> coeffs(solveVanderMatrix());
79 findRoots(cellInfo, coeffs, tetDecom);
83 if (
mag(cutAlpha_) > rootSmall && (1 -
mag(cellAlpha/cutAlpha_)) > 0.1)
89 unweighted_ ? calcAlphaf(cellMagSfs) : calcAlphaUf();
95 void Foam::MPLICcell::findPointAlphaBounds
97 const MPLICcellStorage& cellInfo,
101 const scalar cellAlpha = cellInfo.cellAlpha();
104 UIndirectList<scalar>(cellInfo.pointsAlpha(), cellInfo.cellPoints());
108 cellPointsAlpha_.
append(cellAlpha);
111 sort(cellPointsAlpha_);
115 pointsAlpha_.clear();
116 pointsAlpha_.append(cellPointsAlpha_[0]);
117 for (
label i=1; i<cellPointsAlpha_.
size(); i++)
119 if (
mag(cellPointsAlpha_[i-1] - cellPointsAlpha_[i]) > vSmall)
121 pointsAlpha_.append(cellPointsAlpha_[i]);
125 const label nAlphas = pointsAlpha_.size();
131 pCubicAlphas_.a() = -1;
132 pCubicAlphas_.d() = -1;
133 cCubicAlphas_.a() = -1;
134 cCubicAlphas_.d() = -1;
139 else if (nAlphas == 2)
141 pCubicAlphas_.a() = pointsAlpha_.first();
142 pCubicAlphas_.d() = pointsAlpha_.last();
143 cCubicAlphas_.a() = 1;
144 cCubicAlphas_.d() = 0;
153 label index =
label(round((nAlphas)/2.0)) - 1;
156 scalar target = pointsAlpha_[index];
157 scalar cutAlpha = calcAlpha(cellInfo, target, tetDecom);
158 scalar targetOld = target;
161 for (
label i = 1; i < nAlphas-1; ++i)
163 (cutAlpha >= cellAlpha) ? index++ : index--;
167 if (index == nAlphas - 1)
170 cCubicAlphas_.d() = 0;
171 pCubicAlphas_.a() = target;
172 pCubicAlphas_.d() = pointsAlpha_[index];
180 cCubicAlphas_.a() = 1;
182 pCubicAlphas_.a() = pointsAlpha_[index];
183 pCubicAlphas_.d() = target;
188 target = pointsAlpha_[index];
189 cutAlpha = calcAlpha(cellInfo, target, tetDecom);
191 if (cutAlphaOld > cellAlpha && cutAlpha < cellAlpha)
193 cCubicAlphas_.a() = cutAlphaOld;
195 pCubicAlphas_.a() = targetOld;
196 pCubicAlphas_.d() = target;
200 else if (cutAlphaOld < cellAlpha && cutAlpha > cellAlpha)
203 cCubicAlphas_.d() = cutAlphaOld;
204 pCubicAlphas_.a() = target;
205 pCubicAlphas_.d() = targetOld;
217 void Foam::MPLICcell::calcPointAlphaInterior
219 const MPLICcellStorage& cellInfo,
223 for (
label i=1; i<=2; i++)
226 pCubicAlphas_.a() + (pCubicAlphas_.d() - pCubicAlphas_.a())*(i/3.0);
228 cCubicAlphas_[i] = calcAlpha(cellInfo, pCubicAlphas_[i], tetDecom);
259 const vector4& b = cCubicAlphas_;
261 return FixedList<scalar, 4>
263 scalar(- 4.5*b[0] + 13.5*b[1] - 13.5*b[2] + 4.5*b[3]),
264 scalar(9.0*b[0] - 22.5*b[1] + 18.0*b[2] - 4.5*b[3]),
265 scalar(- 5.5*b[0] + 9.0*b[1] - 4.5*b[2] + 1.0*b[3]),
271 void Foam::MPLICcell::findRoots
273 const MPLICcellStorage& cellInfo,
274 const FixedList<scalar, 4>& coeff,
278 const scalar cellAlpha = cellInfo.cellAlpha();
281 const Roots<3> roots =
282 cubicEqn(coeff[0], coeff[1], coeff[2], coeff[3] - cellAlpha).roots();
285 scalar rootOld =
SMALL;
289 Roots<3> selectedRoots;
291 const scalar pMax =
cmptMax(pCubicAlphas_);
292 const scalar pMin =
cmptMin(pCubicAlphas_);
293 const scalar sFactor = pCubicAlphas_.d() - pCubicAlphas_.a();
298 const scalar root = roots[rooti]*sFactor + pCubicAlphas_.a();
301 if (root < pMax && root > pMin && rootOld != root)
303 target = (target == 0) ? root : target;
304 selectedRoots[nRoots++] = root;
311 cutAlpha_ = calcAlpha(cellInfo, target, tetDecom);
315 scalar error =
mag(cutAlpha_ - cellAlpha);
318 if (nRoots > 0 && error > 1
e-3)
320 scalar minError = error;
323 for (
label rooti=1; rooti<nRoots; rooti++)
325 const scalar targeti =
326 calcAlpha(cellInfo, selectedRoots[rooti], tetDecom);
328 error =
mag(cellAlpha - targeti);
330 if (error < minError)
337 cutAlpha_ = calcAlpha(cellInfo, selectedRoots[minIndex], tetDecom);
342 Foam::scalar Foam::MPLICcell::calcAlpha
344 const MPLICcellStorage& cellInfo,
351 return calcCutCellVolumeAlpha(cellInfo, target);
355 return calcTetCutCellVolumeAlpha(cellInfo, target);
360 void Foam::MPLICcell::calcSubCellVolume()
362 vector cEst = subFaceCentres_[0];
363 for(
label i = 1; i < subFaceCentres_.size(); i++)
365 cEst += subFaceCentres_[i];
367 cEst /= subFaceCentres_.size();
372 subCellVolume_ += subFaceAreas_[i] & (subFaceCentres_[i] - cEst);
374 subCellVolume_ /= 3.0;
378 Foam::scalar Foam::MPLICcell::calcCutCellVolumeAlpha
380 const MPLICcellStorage& cellInfo,
384 const scalar V = cellInfo.V();
387 if (cellInfo.cellAlphaMax() > target && cellInfo.cellAlphaMin() < target)
390 const bool status = singleCutCell(cellInfo, target);
391 if (!status && multiCut_)
393 multiCutCell(cellInfo, target);
400 if (subFaceCentres_.size() != 0)
406 if (subCellVolume_ <= 0)
408 resetFaceFields(cellInfo.size());
413 return min(subCellVolume_, V)/V;
415 else if (target <= cellInfo.cellAlphaMin())
419 subFaceMagSf_ = cellInfo.magSf();
431 resetFaceFields(cellInfo.size());
439 Foam::scalar Foam::MPLICcell::calcTetCutCellVolumeAlpha
441 const MPLICcellStorage& cellInfo,
446 resetFaceFields(cellInfo.size());
449 pointsAlpha_ = cellInfo.pointsAlpha();
450 pointsAlpha_.append(target);
453 scalar cellVolume = 0;
455 if (
min(pointsAlpha_) < target &&
max(pointsAlpha_) > target)
458 const vector& a = cellInfo.C();
462 tetPointsAlpha_[0] = cellInfo.cellAlpha();
466 tetPointsU_[0] = cellInfo.cellU();
473 face f = cellInfo.faces()[cellInfo.cellFaces()[facei]];
476 if (!cellInfo.isOwner()[facei])
481 const label& bL = f[0];
482 const point& b = cellInfo.points()[bL];
483 tetPointsAlpha_[1] = cellInfo.pointsAlpha()[bL];
486 tetPointsU_[1] = cellInfo.pointsU()[bL];
490 for (
label i = 1; i < f.size()-1; ++i)
493 const label cL = f[i];
494 const label dL = f[i + 1];
497 const point& c = cellInfo.points()[cL];
498 const point& d = cellInfo.points()[dL];
501 tetPointsAlpha_[2] = cellInfo.pointsAlpha()[cL];
502 tetPointsAlpha_[3] = cellInfo.pointsAlpha()[dL];
506 tetPointsU_[2] = cellInfo.pointsU()[cL];
507 tetPointsU_[3] = cellInfo.pointsU()[dL];
511 const scalar tetMax =
max(tetPointsAlpha_);
512 const scalar tetMin =
min(tetPointsAlpha_);
518 cellVolume += cellTet.mag();
521 if (tetMin < target && tetMax > target)
524 tetPoints_ = {a,
b,
c, d};
525 tetSf_[0] = cellTet.Sa();
526 tetSf_[1] = cellTet.Sb();
527 tetSf_[2] = cellTet.Sc();
528 tetSf_[3] = cellTet.Sd();
535 const bool ow = cellInfo.isOwner()[facei];
538 cutTetCell(target, facei,ow);
540 if (subFaceCentres_.size() > 0)
547 else if (tetMin >= target)
549 subCellVolume_ += cellTet.mag();
552 subFaceMagSf_[facei] +=
mag(cellTet.Sa());
566 if (cellInfo.isOwner()[facei])
568 alphaPhiU_[facei] += phiU;
572 alphaPhiU_[facei] -= phiU;
583 if (subCellVolume_ <= 0)
585 resetFaceFields(cellInfo.size());
589 return min(subCellVolume_, cellVolume)/cellVolume;
591 else if (target <=
min(pointsAlpha_))
595 subFaceMagSf_ = cellInfo.magSf();
601 subCellVolume_ = cellVolume;
606 resetFaceFields(cellInfo.size());
613 bool Foam::MPLICcell::singleCutCell
615 const MPLICcellStorage& cellInfo,
620 resetFaceFields(cellInfo.size());
626 bool moreCutsPerFace = 0;
632 if (cellInfo.facesAlphaMin()[facei] >= target)
636 cellInfo.Sf()[facei],
637 cellInfo.Cf()[facei],
638 cellInfo.magSf()[facei],
639 cellInfo.isOwner()[facei]
644 subFaceMagSf_[facei] = cellInfo.magSf()[facei];
648 alphaPhiU_[facei] = phiU_[facei];
652 else if (cellInfo.facesAlphaMax()[facei] < target)
660 cellInfo.faces()[cellInfo.cellFaces()[facei]],
662 cellInfo.pointsAlpha(),
665 cellInfo.isOwner()[facei]
674 else if (cutType == 1)
677 cutPoints_.append(faceCutter_.
cutPoints());
683 const vector Cf = faceCutter_.
Cf(Sf);
684 const scalar magSf =
mag(Sf);
685 appendSfCf(Sf, Cf, magSf);
689 subFaceMagSf_[facei] += magSf;
693 alphaPhiU_[facei] += faceCutter_.
alphaPhiU();
700 bool cutOrientationDiffers = 0;
701 if (cutPoints_.size() > 2)
703 cutOrientationDiffers = cutStatusCalcSf();
704 const vector Cf = calcCutCf(cutSf_);
705 appendSfCf(cutSf_, Cf,
mag(cutSf_));
709 if (cutOrientationDiffers || moreCutsPerFace)
722 bool Foam::MPLICcell::multiCutCell
724 const MPLICcellStorage& cellInfo,
729 resetFaceFields(cellInfo.size());
732 if (!addressingCalculated_)
734 calcAddressing(cellInfo);
738 boolList isEdgeCutOld(cellInfo.cellEdges().size(),
false);
739 boolList isEdgeCut(cellInfo.cellEdges().size(),
false);
742 boolList submerged(cellInfo.size(),
false);
745 label facei, nextFace, faceEdgei, status;
758 while (j < cellInfo.size())
760 facei = (status == 0) ? j : nextFace;
763 if (cellInfo.facesAlphaMin()[facei] >= target && !submerged[facei])
765 submerged[facei] =
true;
768 cellInfo.Sf()[facei],
769 cellInfo.Cf()[facei],
770 cellInfo.magSf()[facei],
771 cellInfo.isOwner()[facei]
777 subFaceMagSf_[facei] = cellInfo.magSf()[facei];
781 alphaPhiU_[facei] = phiU_[facei];
788 cellInfo.faces()[cellInfo.cellFaces()[facei]],
789 localFaceEdges_[facei],
794 cellInfo.pointsAlpha(),
798 cellInfo.isOwner()[facei]
804 const label edgei = localFaceEdges_[facei][faceEdgei];
805 const labelList& edgeFaces = localEdgeFaces_[edgei];
806 nextFace = edgeFaces[edgeFaces[0] == facei];
807 faceEdgei =
findIndex(localFaceEdges_[nextFace], edgei);
811 cutPoints_.append(faceCutter_.
cutPoints());
812 cutEdges_.append(faceCutter_.
cutEdges());
815 if (faceCutter_.
subPoints().size() > 2 && !submerged[facei])
818 const vector Cf = faceCutter_.
Cf(Sf);
819 const scalar magSf =
mag(Sf);
822 appendSfCf(Sf, Cf, magSf);
826 subFaceMagSf_[facei] += magSf;
830 alphaPhiU_[facei] += faceCutter_.
alphaPhiU();
835 if (cutEdges_.size() > 0 && cutEdges_.first() == cutEdges_.last())
845 isEdgeCutOld = isEdgeCut;
847 if (cutPoints_.size() == 0)
854 const vector Sf = calcCutSf();
855 const vector Cf = calcCutCf(Sf);
856 appendSfCf(Sf, Cf,
mag(Sf));
868 bool Foam::MPLICcell::cutTetCell
871 const label faceOrig,
881 const face& f = tetFaces_[facei];
890 tetPointsAlpha_[f[0]],
891 tetPointsAlpha_[f[1]]
893 tetPointsAlpha_[f[2]]
897 const vector& Sf = tetSf_[facei];
898 const vector& Cf = tetCf_[facei];
899 appendSfCf(Sf, Cf,
mag(Sf));
901 if (unweighted_ && facei == 0)
903 subFaceMagSf_[faceOrig] +=
mag(Sf);
905 else if (!unweighted_ && facei == 0)
907 const face& f0 = tetFaces_[0];
919 alphaPhiU_[faceOrig] += phiU;
923 alphaPhiU_[faceOrig] -= phiU;
934 tetPointsAlpha_[f[0]],
935 tetPointsAlpha_[f[1]]
937 tetPointsAlpha_[f[2]]
956 cutPoints_.append(faceCutter_.
cutPoints());
962 const vector Cf = faceCutter_.
Cf(Sf);
963 const scalar magSf =
mag(Sf);
964 appendSfCf(Sf, Cf, magSf);
967 if (unweighted_ && facei == 0)
969 subFaceMagSf_[faceOrig] += magSf;
973 else if (!unweighted_ && facei == 0)
977 alphaPhiU_[faceOrig] += faceCutter_.
alphaPhiU();
981 alphaPhiU_[faceOrig] -= faceCutter_.
alphaPhiU();
988 if (cutPoints_.size() > 2)
990 const vector Sf = calcCutSf();
991 const vector Cf = calcCutCf(Sf);
992 appendSfCf(Sf, Cf,
mag(Sf));
1004 unweighted_(unweighted),
1005 multiCut_(multiCut),
1006 faceCutter_(unweighted),
1010 subFaceCentres_(10),
1032 addressingCalculated_ =
false;
1035 label status = calcMatchAlphaCutCell(cellInfo);
1038 if (status == 0 && multiCut_)
1040 status = calcMatchAlphaCutCell(cellInfo,
true);
1043 if (status == 0 || status == -1)
const vector Sf() const
Return subface surface area vector.
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
tetrahedron< point, const point & > tetPointRef
#define forAll(list, i)
Loop across all elements in list.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A 1D vector of objects of type <T> with a fixed size <Size>.
const vector Cf(const vector &area) const
Return subface centre.
void size(const label)
Override size to be inconsistent with allocated storage.
const DynamicList< point > & cutPoints() const
Access to cut points.
const DynamicList< label > & cutEdges() const
Access to cut edges.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Holds information (coordinate and normal) regarding nearest wall point.
static const scalar SMALL
Vector< scalar > vector
A scalar version of the templated Vector.
Provides local cell addressing for geometry and data for MPLIC class.
const dimensionedScalar c
Speed of light in a vacuum.
List< bool > boolList
Bool container classes.
scalar cutAlpha() const
Return volume fraction corresponding to the cut.
Form normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
List< label > labelList
A List of labels.
MPLICcell(const bool unweighted=true, const bool multiCut=true)
Construct for given interpolation and PLIC type.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
vector point
Point is a vector.
label cutFace(const labelList &f, const labelList &faceEdges, const pointField &points, const boolList &isEdgeCutOld, boolList &isEdgeCut, label &faceEdgei, const UList< scalar > &pointsAlpha, const UList< vector > &pointsU, const label facei, const scalar target, const bool ow)
Function to cut for multi cut.
bool matchAlpha(const MPLICcellStorage &cellInfo)
Match cell volume ratios.
triangle< point, const point & > triPointRef
dimensioned< scalar > mag(const dimensioned< Type > &)
const doubleScalar e
Elementary charge.
const DynamicList< point > & subPoints() const
Access to submerged face points.
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
scalar alphaPhiU() const
Calculate and return alphaPhiU.