78 for (
const label gblIdx : stencil[celli])
82 exchangeFields.
getValue(mesh_.
C(), mapCC, gblIdx)
90 cellCentre -= mesh_.
C()[celli];
91 interfaceNormal_[i] = lsGrad.grad(cellCentre, phiValues);
98Foam::reconstruction::gradAlpha::gradAlpha
116 isoFaceTol_(modelDict().getOrDefault<scalar>(
"isoFaceTol", 1
e-8)),
117 surfCellTol_(modelDict().getOrDefault<scalar>(
"surfCellTol", 1
e-8)),
118 sIterPLIC_(mesh_,surfCellTol_)
129 const bool uptodate = alreadyReconstructed(forceUpdate);
131 if (uptodate && !forceUpdate)
136 if (mesh_.topoChanging())
139 interfaceCell_.resize_nocopy(mesh_.nCells());
141 interfaceCell_ =
false;
143 interfaceLabels_.clear();
147 if (sIterPLIC_.isASurfaceCell(alpha1_[celli]))
149 interfaceCell_[celli] =
true;
150 interfaceLabels_.append(celli);
153 interfaceNormal_.resize(interfaceLabels_.size());
159 forAll(interfaceLabels_, i)
161 const label celli = interfaceLabels_[i];
162 if (
mag(interfaceNormal_[i]) == 0)
167 sIterPLIC_.vofCutCell
176 if (sIterPLIC_.cellStatus() == 0)
178 normal_[celli] = sIterPLIC_.surfaceArea();
179 centre_[celli] = sIterPLIC_.surfaceCentre();
180 if (
mag(normal_[celli]) == 0)
182 normal_[celli] =
Zero;
183 centre_[celli] =
Zero;
200 cutCellPLIC cutCell(mesh_);
204 if (
mag(normal_[celli]) != 0)
206 vector n = normal_[celli]/
mag(normal_[celli]);
207 scalar cutValue = (centre_[celli] - mesh_.C()[celli]) & (
n);
214 alpha1_[celli] = cutCell.VolumeOfFluid();
218 alpha1_.correctBoundaryConditions();
219 alpha1_.oldTime () = alpha1_;
220 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.
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
A HashTable to objects of type <T> with a label key.
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.
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.
bool alreadyReconstructed(bool forceUpdate=true) const
Is the interface already up-to-date?
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
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...
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.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Namespace of functions to calculate explicit derivatives.
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< labelList > labelListList
List of labelList.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
#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.