Go to the documentation of this file.
40 void Foam::fvMesh::makeSf()
const
52 <<
"face areas already exist"
73 SfPtr_->setOriented();
77 void Foam::fvMesh::makeMagSf()
const
89 <<
"mag face areas already exist"
113 void Foam::fvMesh::makeC()
const
125 <<
"cell centres already exist"
153 void Foam::fvMesh::makeCf()
const
165 <<
"face centres already exist"
197 <<
"Constructing from primitiveMesh::cellVolumes()" <<
endl;
226 <<
"V0 is not available"
239 <<
"V0 is not available"
281 if (!steady() && moving() && time().subCycling())
284 const TimeState& ts0 = time().prevTimeState();
291 if (tFrac < (1 - SMALL))
293 return V0() + tFrac*(V() - V0());
309 if (!steady() && moving() && time().subCycling())
312 const TimeState& ts0 = time().prevTimeState();
322 return V0() + t0Frac*(V() - V0());
410 const labelUList& neighbour = this->neighbour();
414 delta[facei] =
C[neighbour[facei]] -
C[owner[facei]];
417 surfaceVectorField::Boundary& deltabf =
delta.boundaryFieldRef();
421 deltabf[patchi] =
boundary()[patchi].delta();
433 <<
"mesh flux field does not exist, is the mesh actually moving?"
439 if (!time().subCycling() && phiPtr_->timeIndex() != time().timeIndex())
444 phiPtr_->setOriented();
455 <<
"mesh flux field does not exist, is the mesh actually moving?"
int debug
Static debugging option.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
#define InfoInFunction
Report an information message using Foam::Info.
SlicedGeometricField< vector, fvsPatchField, slicedFvsPatchField, surfaceMesh > slicedSurfaceVectorField
SlicedGeometricField< vector, fvPatchField, slicedFvPatchField, volMesh > slicedVolVectorField
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
static constexpr const zero Zero
Global zero.
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycle old-time cell volumes.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
surfaceScalarField & setPhi()
Return cell face motion fluxes.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type & value() const
Return const reference to value.
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
#define forAll(list, i)
Loop across all elements in list.
const DimensionedField< scalar, volMesh > & V0() const
Return old-time cell volumes.
static int debug
Debug switch.
const fileName & pointsInstance() const
Return the current instance directory for points.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
The time value with time-stepping information, user-defined remapping, etc.
const dimensionSet dimArea(sqr(dimLength))
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const volVectorField & C() const
Return cell centres as volVectorField.
errorManip< error > abort(error &err)
const surfaceScalarField & phi() const
Return cell face motion fluxes.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
IOobject(const word &name, const fileName &instance, const objectRegistry ®istry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
scalar deltaTValue() const
Return time step value.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
The internalField of a SlicedGeometricField.
DimensionedField< scalar, volMesh > & setV0()
Return old-time cell volumes.
Graphite solid properties.
const dimensionSet dimVolume(pow3(dimLength))
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
tmp< surfaceVectorField > delta() const
Return face deltas as surfaceVectorField.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const vectorField & faceAreas() const
const surfaceVectorField & Sf() const
Return cell face area vectors.