Go to the documentation of this file.
49 const GeometricField<Type, fvPatchField, volMesh>& vsf,
55 const fvMesh&
mesh = vsf.mesh();
57 tmp<GeometricField<GradType, fvPatchField, volMesh> > tlsGrad
59 new GeometricField<GradType, fvPatchField, volMesh>
74 pTraits<GradType>::zero
76 zeroGradientFvPatchField<GradType>::typeName
79 GeometricField<GradType, fvPatchField, volMesh>& lsGrad = tlsGrad();
92 label ownFaceI = own[facei];
93 label neiFaceI = nei[facei];
95 Type deltaVsf = vsf[neiFaceI] - vsf[ownFaceI];
97 lsGrad[ownFaceI] += ownLs[facei]*deltaVsf;
98 lsGrad[neiFaceI] -= neiLs[facei]*deltaVsf;
107 lsGrad.boundaryField()[
patchi].patch().faceCells();
109 if (vsf.boundaryField()[
patchi].coupled())
111 const Field<Type> neiVsf
113 vsf.boundaryField()[
patchi].patchNeighbourField()
116 forAll(neiVsf, patchFaceI)
118 lsGrad[faceCells[patchFaceI]] +=
119 patchOwnLs[patchFaceI]
120 *(neiVsf[patchFaceI] - vsf[faceCells[patchFaceI]]);
125 const fvPatchField<Type>& patchVsf = vsf.boundaryField()[
patchi];
127 forAll(patchVsf, patchFaceI)
129 lsGrad[faceCells[patchFaceI]] +=
130 patchOwnLs[patchFaceI]
131 *(patchVsf[patchFaceI] - vsf[faceCells[patchFaceI]]);
137 lsGrad.correctBoundaryConditions();
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
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.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
Mesh data needed to do the Finite Volume discretisation.
typeOfRank< typename pTraits< arg1 >::cmptType, int(pTraits< arg1 >::rank)+int(pTraits< arg2 >::rank) >::type 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
U correctBoundaryConditions()
UList< label > labelUList
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
word name(const complex &)
Return a string representation of a complex.