Go to the documentation of this file.
45 const GeometricField<Type, fvsPatchField, surfaceMesh>& ssf,
51 const fvMesh&
mesh = ssf.mesh();
53 tmp<GeometricField<GradType, fvPatchField, volMesh>> tgGrad
55 new GeometricField<GradType, fvPatchField, volMesh>
67 extrapolatedCalculatedFvPatchField<GradType>::typeName
70 GeometricField<GradType, fvPatchField, volMesh>& gGrad = tgGrad.ref();
76 Field<GradType>& igGrad = gGrad;
77 const Field<Type>& issf = ssf;
81 GradType Sfssf = Sf[facei]*issf[facei];
83 igGrad[owner[facei]] += Sfssf;
84 igGrad[neighbour[facei]] -= Sfssf;
90 mesh.boundary()[patchi].faceCells();
94 const fvsPatchField<Type>& pssf = ssf.boundaryField()[patchi];
98 igGrad[pFaceCells[facei]] += pSf[facei]*pssf[facei];
104 gGrad.correctBoundaryConditions();
122 const GeometricField<Type, fvPatchField, volMesh>& vsf,
128 tmp<GeometricField<GradType, fvPatchField, volMesh>> tgGrad
132 GeometricField<GradType, fvPatchField, volMesh>& gGrad = tgGrad.ref();
153 >::Boundary& gGradbf = gGrad.boundaryFieldRef();
165 gGradbf[patchi] +=
n *
168 - (
n & gGradbf[patchi])
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
static tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > gradf(const GeometricField< Type, fvsPatchField, surfaceMesh > &, const word &name)
Return the gradient of the given field.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
virtual tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > calcGrad(const GeometricField< Type, fvPatchField, volMesh > &vsf, const word &name) const
Return the gradient of the given field to the gradScheme::grad.
Mesh data needed to do the Finite Volume discretisation.
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
#define forAll(list, i)
Loop across all elements in list.
Field< vector > vectorField
Specialisation of Field<T> for vector.
static void correctBoundaryConditions(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > &)
Correct the boundary values of the gradient using the patchField.
cellMask correctBoundaryConditions()
word name(const complex &c)
Return string representation of complex.
UList< label > labelUList
A UList of labels.
Generic GeometricField class.
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
const Boundary & boundaryField() const
Return const-reference to the boundary field.