Go to the documentation of this file.
41 MeshObject<fvMesh,
Foam::MoveableMeshObject, leastSquaresVectors>(
mesh),
47 mesh_.pointsInstance(),
61 mesh_.pointsInstance(),
87 Info<<
"leastSquaresVectors::calcLeastSquaresVectors() :"
88 <<
"Calculating least square gradient vectors"
92 const fvMesh&
mesh = mesh_;
96 const labelUList& neighbour = mesh_.neighbour();
108 label own = owner[facei];
109 label nei = neighbour[facei];
114 dd[own] += (1 -
w[facei])*wdd;
115 dd[nei] +=
w[facei]*wdd;
119 surfaceVectorField::GeometricBoundaryField& blsP =
120 pVectors_.boundaryField();
127 const fvPatch&
p = pw.patch();
128 const labelUList& faceCells =
p.patch().faceCells();
137 const vector& d = pd[patchFacei];
139 dd[faceCells[patchFacei]] +=
140 ((1 - pw[patchFacei])*pMagSf[patchFacei]/
magSqr(d))*
sqr(d);
147 const vector& d = pd[patchFacei];
149 dd[faceCells[patchFacei]] +=
163 label own = owner[facei];
164 label nei = neighbour[facei];
167 scalar magSfByMagSqrd = magSf[facei]/
magSqr(d);
169 pVectors_[facei] = (1 -
w[facei])*magSfByMagSqrd*(invDd[own] & d);
170 nVectors_[facei] = -
w[facei]*magSfByMagSqrd*(invDd[nei] & d);
180 const fvPatch&
p = pw.patch();
190 const vector& d = pd[patchFacei];
192 patchLsP[patchFacei] =
193 ((1 - pw[patchFacei])*pMagSf[patchFacei]/
magSqr(d))
194 *(invDd[faceCells[patchFacei]] & d);
201 const vector& d = pd[patchFacei];
203 patchLsP[patchFacei] =
204 pMagSf[patchFacei]*(1.0/
magSqr(d))
205 *(invDd[faceCells[patchFacei]] & d);
212 Info<<
"leastSquaresVectors::calcLeastSquaresVectors() :"
213 <<
"Finished calculating least square gradient vectors"
221 calcLeastSquaresVectors();
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
virtual ~leastSquaresVectors()
Destructor.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
#define forAll(list, i)
Loop across all elements in list.
fvsPatchField< scalar > fvsPatchScalarField
void calcLeastSquaresVectors()
Calculate Least-squares gradient vectors.
Ostream & endl(Ostream &os)
Add newline and flush stream.
virtual bool movePoints()
Delete the least square vectors when the mesh moves.
void calcLeastSquaresVectors()
Construct Least-squares gradient vectors.
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
fvsPatchField< vector > fvsPatchVectorField
leastSquaresVectors(const fvMesh &)
Construct given an fvMesh.
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
SymmTensor< scalar > symmTensor
SymmTensor of scalars.
static const SymmTensor zero
Vector< scalar > vector
A scalar version of the templated Vector.
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
UList< label > labelUList
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< scalar > magSqr(const dimensioned< Type > &)