40namespace reconstruction
50void Foam::reconstruction::plicRDF::interpolateNormal()
55 leastSquareGrad<scalar> lsGrad(
"polyDegree1",mesh_.
geometricD());
78 DynamicField<vector> cellCentre(100);
79 DynamicField<scalar> alphaValues(100);
81 DynamicList<vector> foundNormals(30);
91 for (
const label gblIdx : stencil[celli])
94 exchangeFields.getValue(
normal_, mapNormal, gblIdx);
100 exchangeFields.getValue(
centre_, mapCentre, gblIdx);
102 estimatedNormal +=
n /
max(
mag(distanceToIntSeg), SMALL);
103 weight += 1/
max(
mag(distanceToIntSeg), SMALL);
104 foundNormals.append(
n);
108 if (weight != 0 &&
mag(estimatedNormal) != 0)
110 estimatedNormal /= weight;
111 estimatedNormal /=
mag(estimatedNormal);
114 bool tooCoarse =
false;
116 if (foundNormals.size() > 1 &&
mag(estimatedNormal) != 0)
122 if ((estimatedNormal & foundNormals[i]) <= 0.98)
135 if (
mag(estimatedNormal) != 0 && !tooCoarse)
137 interfaceNormal_[i] = estimatedNormal;
146 const label gblIdx = stencil[celli][i];
149 exchangeFields.getValue(mesh_.
C(), mapCC, gblIdx)
153 exchangeFields.getValue(
alpha1_, mapAlpha, gblIdx)
156 cellCentre -= mesh_.
C()[celli];
157 interfaceNormal_[i] = lsGrad.grad(cellCentre, alphaValues);
166 leastSquareGrad<scalar> lsGrad(
"polyDegree1", mesh_.geometricD());
169 exchangeFields.setUpCommforZone(interfaceCell_,
false);
173 exchangeFields.getDatafromOtherProc(interfaceCell_, mesh_.C())
177 exchangeFields.getDatafromOtherProc(interfaceCell_,
phi)
180 DynamicField<vector> cellCentre(100);
181 DynamicField<scalar> phiValues(100);
185 forAll(interfaceLabels_, i)
187 const label celli = interfaceLabels_[i];
192 for (
const label gblIdx : stencil[celli])
196 exchangeFields.getValue(mesh_.C(), mapCC, gblIdx)
200 exchangeFields.getValue(
phi, mapPhi, gblIdx)
204 cellCentre -= mesh_.C()[celli];
205 interfaceNormal_[i] = lsGrad.grad(cellCentre, phiValues);
210void Foam::reconstruction::plicRDF::setInitNormals(
bool interpolate)
215 interfaceLabels_.clear();
219 if (sIterPLIC_.isASurfaceCell(alpha1_[celli]))
221 interfaceCell_[celli] =
true;
222 interfaceLabels_.append(celli);
225 interfaceNormal_.setSize(interfaceLabels_.size());
227 RDF_.markCellsNearSurf(interfaceCell_, 1);
228 const boolList& nextToInterface_ = RDF_.nextToInterface();
229 exchangeFields.updateStencil(nextToInterface_);
242void Foam::reconstruction::plicRDF::calcResidual
244 List<normalRes>& normalResidual
249 exchangeFields.setUpCommforZone(interfaceCell_,
false);
251 Map<vector> mapNormal
253 exchangeFields.getDatafromOtherProc(interfaceCell_, normal_)
258 forAll(interfaceLabels_, i)
260 const label celli = interfaceLabels_[i];
261 if (
mag(normal_[celli]) == 0 ||
mag(interfaceNormal_[i]) == 0)
263 normalResidual[i].celli = celli;
264 normalResidual[i].normalResidual = 0;
265 normalResidual[i].avgAngle = 0;
269 scalar avgDiffNormal = 0;
270 scalar maxDiffNormal = GREAT;
272 const vector cellNormal = normal_[celli]/
mag(normal_[celli]);
275 const label gblIdx = stencil[celli][j];
277 exchangeFields.getValue(normal_, mapNormal, gblIdx);
279 if (
mag(normal) != 0 && j != 0)
282 scalar cosAngle =
max(
min((cellNormal &
n), 1), -1);
283 avgDiffNormal +=
acos(cosAngle) *
mag(normal);
284 weight +=
mag(normal);
285 if (cosAngle < maxDiffNormal)
287 maxDiffNormal = cosAngle;
294 avgDiffNormal /= weight;
303 scalar normalRes = (1 - (cellNormal & newCellNormal));
304 normalResidual[i].celli = celli;
305 normalResidual[i].normalResidual = normalRes;
306 normalResidual[i].avgAngle = avgDiffNormal;
331 interfaceNormal_(0.2*mesh_.nCells()),
333 isoFaceTol_(modelDict().getOrDefault<scalar>(
"isoFaceTol", 1
e-8)),
334 surfCellTol_(modelDict().getOrDefault<scalar>(
"surfCellTol", 1
e-8)),
335 tol_(modelDict().getOrDefault<scalar>(
"tol", 1
e-6)),
336 relTol_(modelDict().getOrDefault<scalar>(
"relTol", 0.1)),
337 iteration_(modelDict().getOrDefault<label>(
"iterations", 5)),
338 interpolateNormal_(modelDict().getOrDefault(
"interpolateNormal", true)),
340 sIterPLIC_(mesh_,surfCellTol_)
342 setInitNormals(
false);
350 if (
mag(interfaceNormal_[i]) == 0)
388 const bool uptodate = alreadyReconstructed(forceUpdate);
390 if (uptodate && !forceUpdate)
395 if (mesh_.topoChanging())
398 if (interfaceCell_.size() != mesh_.nCells())
400 interfaceCell_.resize(mesh_.nCells());
403 interfaceCell_ =
false;
406 setInitNormals(interpolateNormal_);
412 const boolList& nextToInterface_ = RDF_.nextToInterface();
414 bitSet tooCoarse(mesh_.nCells(),
false);
416 for (label iter=0; iter<iteration_; ++iter)
418 forAll(interfaceLabels_, i)
420 const label celli = interfaceLabels_[i];
421 if (
mag(interfaceNormal_[i]) == 0 || tooCoarse.
test(celli))
425 sIterPLIC_.vofCutCell
434 if (sIterPLIC_.cellStatus() == 0)
437 normal_[celli] = sIterPLIC_.surfaceArea();
438 centre_[celli] = sIterPLIC_.surfaceCentre();
439 if (
mag(normal_[celli]) == 0)
441 normal_[celli] =
Zero;
442 centre_[celli] =
Zero;
447 normal_[celli] =
Zero;
448 centre_[celli] =
Zero;
452 normal_.correctBoundaryConditions();
453 centre_.correctBoundaryConditions();
457 nHatb *= 1/(mesh_.magSf().boundaryField());
468 RDF_.updateContactAngle(alpha1_, U_, nHatb);
470 calcResidual(normalResidual);
473 label resCounter = 0;
475 scalar avgNormRes = 0;
480 const label celli = normalResidual[i].celli;
481 const scalar normalRes= normalResidual[i].normalResidual;
482 const scalar avgA = normalResidual[i].avgAngle;
484 if (avgA > 0.26 && iter > 0)
486 tooCoarse.
set(celli);
492 scalar discreteError = 0.01*
sqr(avgA);
493 if (discreteError != 0)
495 normRes= normalRes/
max(discreteError, tol_);
499 normRes= normalRes/tol_;
501 avgNormRes += normRes;
522 <<
"initial residual absolute = "
523 << avgRes/resCounter <<
nl
524 <<
"initial residual normalized = "
525 << avgNormRes/resCounter <<
nl;
531 (avgNormRes/resCounter < relTol_ || avgRes/resCounter < tol_)
534 || iter + 1 == iteration_
538 <<
"iterations = " << iter <<
nl
539 <<
"final residual absolute = "
540 << avgRes/resCounter <<
nl
541 <<
"final residual normalized = " << avgNormRes/resCounter
559 if (
mag(normal_[celli]) != 0)
561 vector n = normal_[celli]/
mag(normal_[celli]);
562 scalar cutValue = (centre_[celli] - mesh_.C()[celli]) & (
n);
569 alpha1_[celli] =
cutCell.VolumeOfFluid();
573 alpha1_.correctBoundaryConditions();
574 alpha1_.oldTime () = alpha1_;
575 alpha1_.oldTime().correctBoundaryConditions();
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const volScalarField & alpha1
scalar deltaTValue() const noexcept
Return time step value.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
void set(const bitSet &bitset)
Set specified bits from another bitset.
bool test(const label pos) const
Test value at specified position, never auto-vivify entries.
Class for cutting a cell, cellI, of an fvMesh, mesh_, at its intersection with an surface defined by ...
Service routines for cutting a cell, celli, of an fvMesh, mesh_, at its intersection with a surface.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
const volVectorField & C() const
Return cell centres as volVectorField.
const Time & time() const
Return the top-level database.
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Original code supplied by Henning Scheufler, DLR (2019)
boolList interfaceCell_
Is interface cell?
volVectorField normal_
Interface area normals.
const volVectorField & U_
Reference to the velocity field.
const volVectorField & centre() const noexcept
const-Reference to interface centres
DynamicField< label > interfaceLabels_
List of cell labels that have an interface.
volScalarField & alpha1_
Reference to the VoF Field.
volVectorField centre_
Interface centres.
Reconstructs an interface (centre and normal vector) consisting of planes to match the internal fluid...
virtual void reconstruct(bool forceUpdate=true)
Reconstruct interface.
virtual void mapAlphaField() const
Map VoF Field in case of refinement.
const point & surfaceCentre() const
The centre of cutting isosurface.
label cellStatus()
The cellStatus.
label vofCutCell(const label celli, const scalar alpha1, const scalar tol, const label maxIter, vector normal)
Finds matching cutValue for the given value fraction.
const vector & surfaceArea() const
The area vector of cutting isosurface.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Class for parallel communication in a narrow band. It either provides a Map with the neighbouring val...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define DebugInfo
Report an information message using Foam::Info.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
vector point
Point is a vector.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
List< labelList > labelListList
A List of labelList.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
List< bool > boolList
A List of bools.
static constexpr const zero Zero
Global zero (0)
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
dimensionedScalar acos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a)
#define addProfilingInFunction(name)
#define forAll(list, i)
Loop across all elements in list.