36template<
class Type,
class Limiter,
template<
class>
class LimitFunc>
37void Foam::LimitedScheme<Type, Limiter, LimitFunc>::calcLimiter
39 const GeometricField<Type, fvPatchField, volMesh>&
phi,
40 surfaceScalarField& limiterField
43 typedef GeometricField<typename Limiter::phiType, fvPatchField, volMesh>
46 typedef GeometricField<typename Limiter::gradPhiType, fvPatchField, volMesh>
51 tmp<VolFieldType> tlPhi = LimitFunc<Type>()(
phi);
52 const VolFieldType& lPhi = tlPhi();
54 tmp<GradVolFieldType> tgradc(fvc::grad(lPhi));
55 const GradVolFieldType& gradc = tgradc();
64 scalarField& pLim = limiterField.primitiveFieldRef();
68 label own = owner[face];
69 label nei = neighbour[face];
71 pLim[face] = Limiter::limiter
74 this->faceFlux_[face],
83 surfaceScalarField::Boundary& bLim = limiterField.boundaryFieldRef();
91 const scalarField& pCDweights = CDweights.boundaryField()[patchi];
93 this->faceFlux_.boundaryField()[patchi];
95 const Field<typename Limiter::phiType> plPhiP
97 lPhi.boundaryField()[patchi].patchInternalField()
99 const Field<typename Limiter::phiType> plPhiN
101 lPhi.boundaryField()[patchi].patchNeighbourField()
103 const Field<typename Limiter::gradPhiType> pGradcP
105 gradc.boundaryField()[patchi].patchInternalField()
107 const Field<typename Limiter::gradPhiType> pGradcN
109 gradc.boundaryField()[patchi].patchNeighbourField()
113 vectorField pd(CDweights.boundaryField()[patchi].patch().delta());
117 pLim[face] = Limiter::limiter
135 limiterField.setOriented();
141template<
class Type,
class Limiter,
template<
class>
class LimitFunc>
150 const word limiterFieldName(
type() +
"Limiter(" +
phi.name() +
')');
152 if (this->
mesh().cache(
"limiter"))
163 mesh.time().timeName(),
175 auto& limiterField = *fldptr;
177 calcLimiter(
phi, limiterField);
187 auto tlimiterField = surfaceScalarField::New
194 calcLimiter(
phi, tlimiterField.ref());
196 return tlimiterField;
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
Generic GeometricField class.
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvsPatchField< scalar >::calculatedType())
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
virtual tmp< surfaceScalarField > limiter(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
Mesh data needed to do the Finite Volume discretisation.
bool store()
Register object with its registry and transfer ownership to the registry.
const fvMesh & mesh() const
Return mesh reference.
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
A class for handling words, derived from Foam::string.
Calculate the gradient of the given field.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
UList< label > labelUList
A UList of labels.
#define forAll(list, i)
Loop across all elements in list.