38#ifndef interfaceTrackingFvMesh_H
39#define interfaceTrackingFvMesh_H
58class interfaceTrackingFvMesh
60 public dynamicMotionSolverFvMesh
74 wordList nonReflectingFreeSurfacePatches_;
77 wordList pointNormalsCorrectionPatches_;
97 Switch correctContactLineNormals_;
109 mutable std::unique_ptr<areaVectorField> UsPtr_;
113 mutable std::unique_ptr<vectorIOField> controlPointsPtr_;
117 mutable std::unique_ptr<labelList> motionPointsMaskPtr_;
120 mutable std::unique_ptr<vectorField> pointsDisplacementDirPtr_;
123 mutable std::unique_ptr<vectorField> facesDisplacementDirPtr_;
126 mutable std::unique_ptr<areaScalarField> fsNetPhiPtr_;
129 mutable std::unique_ptr<edgeScalarField> phisPtr_;
132 mutable std::unique_ptr<areaScalarField> surfactConcPtr_;
135 mutable std::unique_ptr<volScalarField> bulkSurfactConcPtr_;
138 mutable std::unique_ptr<areaScalarField> surfaceTensionPtr_;
141 mutable std::unique_ptr<surfactantProperties> surfactantPtr_;
144 mutable std::unique_ptr<areaScalarField> contactAnglePtr_;
155 void makeControlPoints();
158 void makeMotionPointsMask();
161 void makeDirections();
164 void makeFsNetPhi()
const;
170 void makeBulkSurfactConc()
const;
173 void makeSurfactConc()
const;
176 void makeSurfaceTension()
const;
179 void makeSurfactant()
const;
182 void makeContactAngle();
185 interfaceTrackingFvMesh(
const interfaceTrackingFvMesh&) =
delete;
188 void operator=(
const interfaceTrackingFvMesh&) =
delete;
191 void initializeData();
194 void updateDisplacementDirections();
197 void initializeControlPointsPosition();
212 void smoothFreeSurfaceMesh();
215 void updateSurfaceFlux();
218 void updateSurfactantConcentration();
221 vector totalPressureForce()
const;
224 vector totalViscousForce()
const;
227 vector totalSurfaceTensionForce()
const;
230 scalar maxCourantNumber();
233 void updateProperties();
236 void correctContactLinePointNormals();
239 void correctPointDisplacement
248 TypeName(
"interfaceTrackingFvMesh");
254 interfaceTrackingFvMesh(
const IOobject&
io,
const bool doInit=
true);
276 virtual bool init(
const bool doInit);
302 return fsPatchIndex_;
308 return pureFreeSurface_;
314 return rigidFreeSurface_;
320 return rigidFreeSurface_;
326 return correctContactLineNormals_;
332 return correctContactLineNormals_;
392 return normalMotionDir_;
428 controlPointsPtr_.reset(
nullptr);
Defines the attributes of an object for which implicit objectRegistry management is supported,...
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Mesh data needed to do the Finite Volume discretisation.
vectorField & controlPoints()
Return control points.
void correctUsBoundaryConditions()
Correct surface velocity boundary conditions.
Switch pureFreeSurface() const noexcept
Pure free-surface.
areaVectorField & Us()
Return free-surface velocity field.
TypeName("interfaceTrackingFvMesh")
Runtime type information.
const areaScalarField & fsNetPhi() const
Return free-surface net flux.
const surfaceScalarField & phi() const
Return constant reference to velocity field.
tmp< scalarField > freeSurfaceSnGradUn()
Return free surface normal derivative of normal velocity comp.
vectorField & pointsDisplacementDir()
Return reference to point displacement direction field.
Switch rigidFreeSurface() const noexcept
Rigid free-surface.
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
const volVectorField & U() const
Return constant reference to velocity field.
Switch correctContactLineNormals() const noexcept
Correct contact line point normals.
void writeVTK() const
Write VTK freeSurface mesh.
void updateUs()
Update free-surface velocity field.
areaScalarField & surfaceTension()
Return surface tension field.
vectorField & facesDisplacementDir()
Return reference to control points displacement direction field.
void writeVTKControlPoints()
Write VTK freeSurface control points.
faMesh & aMesh()
Return reference to finite area mesh.
virtual bool update()
Update the mesh for both mesh motion and topology change.
Switch & correctContactLineNormals() noexcept
Correct contact line point normals.
const label & fsPatchIndex() const
const fvMesh & mesh() const
interfaceTrackingFvMesh(const IOobject &io, pointField &&points, faceList &&faces, labelList &&allOwner, labelList &&allNeighbour, const bool syncPar=true)
* ~interfaceTrackingFvMesh()
Destructor.
const dimensionedScalar & sigma() const
Surface tension coefficient of pure free-surface.
tmp< areaVectorField > surfaceTensionGrad()
Return surface tension gradient.
Switch & rigidFreeSurface() noexcept
Rigid free-surface.
labelList & motionPointsMask()
Return reference to motion points mask field.
volScalarField & bulkSurfactantConcentration()
Return volume surfactant concentration field.
bool normalMotionDir() const
Motion direction swithc.
const volScalarField & p() const
Return constant reference to pressure field.
const surfactantProperties & surfactant() const
Return surfactant properties.
areaScalarField & surfactantConcentration()
Return free-surface surfactant concentration field.
edgeScalarField & Phis()
Return free-surface fluid flux field.
tmp< scalarField > freeSurfacePressureJump()
Return free surface pressure jump.
tmp< vectorField > freeSurfaceSnGradU()
Return free surface normal derivative of velocity.
const faMesh & aMesh() const
Return reference to finite area mesh.
void clearControlPoints()
Clear control points.
virtual const faceList & faces() const
Return raw faces.
A class for managing temporary objects.
List< word > wordList
List of word.
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
List< face > faceList
List of faces.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
GeometricField< vector, faPatchField, areaMesh > areaVectorField
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Field< vector > vectorField
Specialisation of Field<T> for vector.
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.