50void Foam::reconstruction::plicRDF::interpolateNormal()
91 for (
const label gblIdx : stencil[celli])
95 point p = mesh_.C()[celli]-
U_[celli]*dt;
110 exchangeFields.
getValue(mesh_.C(), mapCC, gblIdx)
115 estimatedNormal +=
n /
max(
mag(distanceToIntSeg), SMALL);
116 weight += 1/
max(
mag(distanceToIntSeg), SMALL);
117 foundNormals.append(
n);
121 if (weight != 0 &&
mag(estimatedNormal) != 0)
123 estimatedNormal /= weight;
124 estimatedNormal /=
mag(estimatedNormal);
127 bool tooCoarse =
false;
129 if (foundNormals.size() > 1 &&
mag(estimatedNormal) != 0)
135 if ((estimatedNormal & foundNormals[i]) <= 0.98)
148 if (
mag(estimatedNormal) != 0 && !tooCoarse)
150 interfaceNormal_[i] = estimatedNormal;
159 const label gblIdx = stencil[celli][i];
164 mesh_.C(), mapCC, gblIdx,
169 exchangeFields.
getValue(mesh_.C(), mapCC, gblIdx)
178 cellCentre -= mesh_.C()[celli];
179 interfaceNormal_[i] = lsGrad.grad(cellCentre, alphaValues);
188 leastSquareGrad<scalar> lsGrad(
"polyDegree1", mesh_.geometricD());
191 exchangeFields.setUpCommforZone(interfaceCell_,
false);
195 exchangeFields.getDatafromOtherProc(interfaceCell_, mesh_.C())
199 exchangeFields.getDatafromOtherProc(interfaceCell_,
phi)
202 DynamicField<vector> cellCentre(100);
203 DynamicField<scalar> phiValues(100);
207 forAll(interfaceLabels_, i)
209 const label celli = interfaceLabels_[i];
214 for (
const label gblIdx : stencil[celli])
218 exchangeFields.getPosition
220 mesh_.C(), mapCC, gblIdx,
221 exchangeFields.getCyclicPatches
225 exchangeFields.getValue(mesh_.C(), mapCC, gblIdx)
231 exchangeFields.getValue(
phi, mapPhi, gblIdx)
235 cellCentre -= mesh_.C()[celli];
236 interfaceNormal_[i] = lsGrad.grad(cellCentre, phiValues);
241void Foam::reconstruction::plicRDF::setInitNormals(
bool interpolate)
246 interfaceLabels_.clear();
250 if (sIterPLIC_.isASurfaceCell(alpha1_[celli]))
252 interfaceCell_[celli] =
true;
253 interfaceLabels_.append(celli);
256 interfaceNormal_.setSize(interfaceLabels_.size());
258 RDF_.markCellsNearSurf(interfaceCell_, 1);
259 const boolList& nextToInterface_ = RDF_.nextToInterface();
260 exchangeFields.updateStencil(nextToInterface_);
273void Foam::reconstruction::plicRDF::calcResidual
275 List<normalRes>& normalResidual
280 exchangeFields.setUpCommforZone(interfaceCell_,
false);
282 Map<vector> mapNormal
284 exchangeFields.getDatafromOtherProc(interfaceCell_, normal_)
289 forAll(interfaceLabels_, i)
291 const label celli = interfaceLabels_[i];
292 if (
mag(normal_[celli]) == 0 ||
mag(interfaceNormal_[i]) == 0)
294 normalResidual[i].celli = celli;
295 normalResidual[i].normalResidual = 0;
296 normalResidual[i].avgAngle = 0;
300 scalar avgDiffNormal = 0;
301 scalar maxDiffNormal = GREAT;
303 const vector cellNormal = normal_[celli]/
mag(normal_[celli]);
306 const label gblIdx = stencil[celli][j];
308 exchangeFields.getValue(normal_, mapNormal, gblIdx);
310 if (
mag(normal) != 0 && j != 0)
313 scalar cosAngle =
clamp((cellNormal &
n), -1, 1);
314 avgDiffNormal +=
acos(cosAngle) *
mag(normal);
315 weight +=
mag(normal);
316 if (cosAngle < maxDiffNormal)
318 maxDiffNormal = cosAngle;
325 avgDiffNormal /= weight;
334 scalar normalRes = (1 - (cellNormal & newCellNormal));
335 normalResidual[i].celli = celli;
336 normalResidual[i].normalResidual = normalRes;
337 normalResidual[i].avgAngle = avgDiffNormal;
344Foam::reconstruction::plicRDF::plicRDF
362 interfaceNormal_(0.2*mesh_.nCells()),
364 isoFaceTol_(modelDict().getOrDefault<scalar>(
"isoFaceTol", 1
e-8)),
365 surfCellTol_(modelDict().getOrDefault<scalar>(
"surfCellTol", 1
e-8)),
366 tol_(modelDict().getOrDefault<scalar>(
"tol", 1
e-6)),
367 relTol_(modelDict().getOrDefault<scalar>(
"relTol", 0.1)),
368 iteration_(modelDict().getOrDefault<label>(
"iterations", 5)),
369 interpolateNormal_(modelDict().getOrDefault(
"interpolateNormal", true)),
371 sIterPLIC_(mesh_,surfCellTol_)
373 setInitNormals(
false);
381 if (
mag(interfaceNormal_[i]) == 0)
419 const bool uptodate = alreadyReconstructed(forceUpdate);
421 if (uptodate && !forceUpdate)
426 if (mesh_.topoChanging())
429 if (interfaceCell_.size() != mesh_.nCells())
431 interfaceCell_.resize(mesh_.nCells());
434 interfaceCell_ =
false;
437 setInitNormals(interpolateNormal_);
443 const boolList& nextToInterface_ = RDF_.nextToInterface();
445 bitSet tooCoarse(mesh_.nCells(),
false);
447 for (label iter=0; iter<iteration_; ++iter)
449 forAll(interfaceLabels_, i)
451 const label celli = interfaceLabels_[i];
452 if (
mag(interfaceNormal_[i]) == 0 || tooCoarse.test(celli))
456 sIterPLIC_.vofCutCell
465 if (sIterPLIC_.cellStatus() == 0)
468 normal_[celli] = sIterPLIC_.surfaceArea();
469 centre_[celli] = sIterPLIC_.surfaceCentre();
470 if (
mag(normal_[celli]) == 0)
472 normal_[celli] =
Zero;
473 centre_[celli] =
Zero;
478 normal_[celli] =
Zero;
479 centre_[celli] =
Zero;
483 normal_.correctBoundaryConditions();
484 centre_.correctBoundaryConditions();
485 List<normalRes> normalResidual(interfaceLabels_.size());
488 nHatb *= 1/(mesh_.magSf().boundaryField());
499 RDF_.updateContactAngle(alpha1_, U_, nHatb);
501 calcResidual(normalResidual);
504 label resCounter = 0;
506 scalar avgNormRes = 0;
511 const label celli = normalResidual[i].celli;
512 const scalar normalRes= normalResidual[i].normalResidual;
513 const scalar avgA = normalResidual[i].avgAngle;
515 if (avgA > 0.26 && iter > 0)
517 tooCoarse.set(celli);
523 scalar discreteError = 0.01*
sqr(avgA);
524 if (discreteError != 0)
526 normRes= normalRes/
max(discreteError, tol_);
530 normRes= normalRes/tol_;
532 avgNormRes += normRes;
538 reduce(avgRes,sumOp<scalar>());
539 reduce(avgNormRes,sumOp<scalar>());
540 reduce(resCounter,sumOp<label>());
553 <<
"initial residual absolute = "
554 << avgRes/resCounter <<
nl
555 <<
"initial residual normalized = "
556 << avgNormRes/resCounter <<
nl;
562 (avgNormRes/resCounter < relTol_ || avgRes/resCounter < tol_)
565 || iter + 1 == iteration_
569 <<
"iterations = " << iter <<
nl
570 <<
"final residual absolute = "
571 << avgRes/resCounter <<
nl
572 <<
"final residual normalized = " << avgNormRes/resCounter
586 cutCellPLIC cutCell(mesh_);
590 if (
mag(normal_[celli]) != 0)
592 vector n = normal_[celli]/
mag(normal_[celli]);
593 scalar cutValue = (centre_[celli] - mesh_.C()[celli]) & (
n);
600 alpha1_[celli] = cutCell.VolumeOfFluid();
604 alpha1_.correctBoundaryConditions();
605 alpha1_.oldTime () = alpha1_;
606 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
Dynamically sized Field. Similar to DynamicList, but inheriting from a Field instead of a List.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void append(const T &val)
Append an element at the end of the list.
void clear()
Clear the list, i.e. set size to zero.
A HashTable to objects of type <T> with a label key.
scalar deltaTValue() const noexcept
Return time step value.
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 for True 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,...
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
const volVectorField & C() const
Return cell centres as volVectorField.
const Time & time() const
Return the top-level database.
Estimates the gradient with a least square scheme in a cell.
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.
bool alreadyReconstructed(bool forceUpdate=true) const
Is the interface already up-to-date?
const volVectorField & centre() const noexcept
const-Reference to interface centres
DynamicField< label > interfaceLabels_
List of cell labels that have an interface.
reconstructionSchemes(const reconstructionSchemes &)=delete
No copy construct.
volScalarField & alpha1_
Reference to the VoF Field.
volVectorField centre_
Interface centres.
dictionary & modelDict() noexcept
Access to the model dictionary.
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...
List< label > getCyclicPatches(const label celli, const label globalIdx, const vector globalIdxCellCentre) const
Finds and returns list of all cyclic patch labels to which celli's.
void setUpCommforZone(const boolList &zone, bool updateStencil=true)
Update stencil with boolList the size has to match mesh nCells.
Type getValue(const VolumeField< Type > &phi, const Map< Type > &valuesFromOtherProc, const label gblIdx) const
Gives patchNumber and patchFaceNumber for a given Geometric volume field.
static zoneDistribute & New(const fvMesh &)
Selector.
const labelListList & getStencil() noexcept
Stencil reference.
Map< Type > getDatafromOtherProc(const boolList &zone, const VolumeField< Type > &phi)
Returns stencil and provides a Map with globalNumbering.
vector getPosition(const VolumeField< vector > &positions, const Map< vector > &valuesFromOtherProc, const label gblIdx, const List< label > cyclicPatchID=List< label >()) const
#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.
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< labelList > labelListList
List of labelList.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
vector point
Point is a vector.
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 profiling trigger with specified name and description corresponding to the compiler-defined fu...
#define forAll(list, i)
Loop across all elements in list.