Calculates the face fluxes. More...
#include <cutFaceAdvect.H>


Public Member Functions | |
| cutFaceAdvect (const fvMesh &mesh, const volScalarField &alpha1) | |
| Construct from fvMesh and a scalarField. | |
| label | calcSubFace (const label faceI, const vector &normal, const vector &base) |
| Calculates cut centre and cut area (plicReconstruction). | |
| label | calcSubFace (const face &f, const pointField &points, const scalarField &val, const scalar cutValue) |
| Calculates cut centre and cut area (iso reconstruction). | |
| label | calcSubFace (const label faceI, const label sign, const scalar cutValue) |
| Calculates cut centre and cut area (iso reconstruction). | |
| scalar | timeIntegratedFaceFlux (const label faceI, const vector &x0, const vector &n0, const scalar Un0, const scalar dt, const scalar phi, const scalar magSf) |
| Calculate time integrated flux for a face. | |
| scalar | timeIntegratedFaceFlux (const label faceI, const scalarField ×, const scalar Un0, const scalar dt, const scalar phi, const scalar magSf) |
| Calculate time integrated flux for a face. | |
| scalar | timeIntegratedArea (const label faceI, const scalar dt, const scalar magSf, const scalar Un0) |
| Calculate time integrated area for a face. | |
| scalar | timeIntegratedArea (const pointField &fPts, const scalarField &pTimes, const scalar dt, const scalar magSf, const scalar Un0) |
| Calculate time integrated area for a face. | |
| void | quadAreaCoeffs (const DynamicPointList &pf0, const DynamicPointList &pf1, scalar &quadArea, scalar &intQuadArea) const |
| For face with vertices p calculate its area and integrated area. | |
| void | cutPoints (const label faceI, const scalar f0, DynamicList< point > &cutPoints) |
| Get cutPoints of face. | |
| void | cutPoints (const pointField &pts, const scalarField &f, const scalar f0, DynamicList< point > &cutPoints) |
| Get cutPoints of face. | |
| const point & | subFaceCentre () const noexcept |
| Returns centre of cutted face. | |
| const vector & | subFaceArea () const noexcept |
| Returns area vector of cutted face. | |
| const DynamicList< point > & | subFacePoints () const noexcept |
| Returns the cut edge of the cutted face. | |
| const DynamicList< point > & | surfacePoints () const noexcept |
| Returns point of the face in sorted of cutted face. | |
| void | clearStorage () |
| Resets internal variables. | |
| Public Member Functions inherited from cutFace | |
| cutFace (const fvMesh &mesh) | |
| Construct from fvMesh. | |
Additional Inherited Members | |
| Static Public Attributes inherited from cutFace | |
| static int | debug = 0 |
| Protected Member Functions inherited from cutFace | |
| void | calcSubFace (const label faceI, const scalarList &pointStatus, label firstFullySubmergedPoint, DynamicList< point > &subFacePoints, DynamicList< point > &surfacePoints, label &faceStatus, vector &subFaceCentre, vector &subFaceArea) |
| Calculate cut points along edges of face with pointStatus, pointfield and computes geometric information and the face status where: | |
| void | calcSubFace (const label faceI, const scalarList &pointStatus, const scalarList &weights, label firstFullySubmergedPoint, DynamicList< point > &subFacePoints, DynamicList< point > &surfacePoints, label &faceStatus, vector &subFaceCentre, vector &subFaceArea) |
| Calculate cut points along edges of face with pointStatus, pointfield, weights and computes geometric information and the face status where: | |
| void | calcSubFaceCentreAndArea (DynamicList< point > &subFacePoints, vector &subFaceCentre, vector &subFaceArea) |
| Calculates centre and normal of the face. | |
| void | calcSubFace (const face &f, const pointField &points, const scalarList &pointStatus, label firstFullySubmergedPoint, DynamicList< point > &subFacePoints, DynamicList< point > &surfacePoints, label &faceStatus, vector &subFaceCentre, vector &subFaceArea) |
| Calculate cut points along edges of face with pointStatus, pointfield and computes geometric information and the face status where: | |
Calculates the face fluxes.
Reference:
Roenby, J., Bredmose, H. and Jasak, H. (2016).
A computational method for sharp interface advection
Royal Society Open Science, 3
doi 10.1098/rsos.160405
Original code supplied by Johan Roenby, DHI (2016)
Definition at line 62 of file cutFaceAdvect.H.
| cutFaceAdvect | ( | const fvMesh & | mesh, |
| const volScalarField & | alpha1 ) |
Construct from fvMesh and a scalarField.
Length of scalarField should equal number of mesh points
Definition at line 28 of file cutFaceAdvect.C.
References alpha1, clearStorage(), cutFace::cutFace(), mesh, and Foam::Zero.

Calculates cut centre and cut area (plicReconstruction).
Definition at line 52 of file cutFaceAdvect.C.
References cutFace::calcSubFace(), clearStorage(), f(), forAll, Foam::mag(), and Foam::Zero.
Referenced by timeIntegratedArea(), timeIntegratedArea(), and timeIntegratedFaceFlux().


| Foam::label calcSubFace | ( | const face & | f, |
| const pointField & | points, | ||
| const scalarField & | val, | ||
| const scalar | cutValue ) |
Calculates cut centre and cut area (iso reconstruction).
Definition at line 116 of file cutFaceAdvect.C.
References cutFace::calcSubFace(), clearStorage(), f(), forAll, Foam::mag(), points, and Foam::Zero.

| Foam::label calcSubFace | ( | const label | faceI, |
| const label | sign, | ||
| const scalar | cutValue ) |
Calculates cut centre and cut area (iso reconstruction).
Definition at line 569 of file cutFaceAdvect.C.
References cutFace::calcSubFace(), f(), forAll, Foam::mag(), Foam::sign(), and Foam::Zero.

| Foam::scalar timeIntegratedFaceFlux | ( | const label | faceI, |
| const vector & | x0, | ||
| const vector & | n0, | ||
| const scalar | Un0, | ||
| const scalar | dt, | ||
| const scalar | phi, | ||
| const scalar | magSf ) |
Calculate time integrated flux for a face.
Definition at line 180 of file cutFaceAdvect.C.
References calcSubFace(), clearStorage(), cutFace::debug, e, Foam::endl(), f(), forAll, Foam::mag(), nPoints, phi, pi(), Foam::sign(), subFaceArea(), timeIntegratedArea(), and WarningInFunction.

| Foam::scalar timeIntegratedFaceFlux | ( | const label | faceI, |
| const scalarField & | times, | ||
| const scalar | Un0, | ||
| const scalar | dt, | ||
| const scalar | phi, | ||
| const scalar | magSf ) |
Calculate time integrated flux for a face.
Definition at line 327 of file cutFaceAdvect.C.
References clearStorage(), forAll, nPoints, phi, pi(), Foam::sign(), UList< T >::size(), and timeIntegratedArea().

| Foam::scalar timeIntegratedArea | ( | const label | faceI, |
| const scalar | dt, | ||
| const scalar | magSf, | ||
| const scalar | Un0 ) |
Calculate time integrated area for a face.
Definition at line 632 of file cutFaceAdvect.C.
References alpha, DynamicList< T, SizeMin >::append(), beta(), calcSubFace(), cutPoints(), e, UList< T >::first(), UList< T >::last(), Foam::mag(), Foam::max(), Foam::neg(), Foam::pos0(), quadAreaCoeffs(), Foam::sign(), Foam::sortedOrder(), and subFaceArea().
Referenced by timeIntegratedFaceFlux(), and timeIntegratedFaceFlux().


| Foam::scalar timeIntegratedArea | ( | const pointField & | fPts, |
| const scalarField & | pTimes, | ||
| const scalar | dt, | ||
| const scalar | magSf, | ||
| const scalar | Un0 ) |
Calculate time integrated area for a face.
Definition at line 374 of file cutFaceAdvect.C.
References alpha, DynamicList< T, SizeMin >::append(), beta(), calcSubFace(), cutPoints(), e, UList< T >::first(), Foam::identity(), UList< T >::last(), Foam::mag(), Foam::max(), Foam::neg(), Foam::pos0(), quadAreaCoeffs(), Foam::sign(), UList< T >::size(), Foam::sortedOrder(), and subFaceArea().

| void quadAreaCoeffs | ( | const DynamicPointList & | pf0, |
| const DynamicPointList & | pf1, | ||
| scalar & | quadArea, | ||
| scalar & | intQuadArea ) const |
For face with vertices p calculate its area and integrated area.
between isocutting lines with isovalues f0 and f1 given vertex values f
Definition at line 772 of file cutFaceAdvect.C.
References A, alpha, B, beta(), C, D, Foam::endl(), Foam::mag(), UList< T >::size(), WarningInFunction, and Foam::Zero.
Referenced by timeIntegratedArea(), and timeIntegratedArea().


| void cutPoints | ( | const label | faceI, |
| const scalar | f0, | ||
| DynamicList< point > & | cutPoints ) |
Get cutPoints of face.
Definition at line 874 of file cutFaceAdvect.C.
References cutPoints(), Foam::endl(), f(), forAll, Foam::mag(), nPoints, pi(), s(), and WarningInFunction.
Referenced by cutPoints(), cutPoints(), timeIntegratedArea(), and timeIntegratedArea().


| void cutPoints | ( | const pointField & | pts, |
| const scalarField & | f, | ||
| const scalar | f0, | ||
| DynamicList< point > & | cutPoints ) |
Get cutPoints of face.
Definition at line 929 of file cutFaceAdvect.C.
References cutPoints(), Foam::endl(), f(), forAll, Foam::mag(), nPoints, pi(), pts, s(), and WarningInFunction.

|
inlinenoexcept |
Returns centre of cutted face.
Definition at line 283 of file cutFaceAdvect.H.
References Foam::noexcept.
|
inlinenoexcept |
Returns area vector of cutted face.
Definition at line 291 of file cutFaceAdvect.H.
References Foam::noexcept.
Referenced by timeIntegratedArea(), timeIntegratedArea(), and timeIntegratedFaceFlux().

|
inlinenoexcept |
Returns the cut edge of the cutted face.
Definition at line 299 of file cutFaceAdvect.H.
References Foam::noexcept.
|
inlinenoexcept |
Returns point of the face in sorted of cutted face.
Definition at line 307 of file cutFaceAdvect.H.
References Foam::noexcept.
| void clearStorage | ( | ) |
Resets internal variables.
Definition at line 979 of file cutFaceAdvect.C.
References Foam::Zero.
Referenced by calcSubFace(), calcSubFace(), cutFaceAdvect(), timeIntegratedFaceFlux(), and timeIntegratedFaceFlux().
