Go to the documentation of this file.
30 #include "surfaceInterpolate.H"
45 namespace patchDistMethods
54 Foam::patchDistMethods::advectionDiffusion::advectionDiffusion
65 coeffs_(
dict.subDict(
type() +
"Coeffs")),
75 epsilon_(coeffs_.getOrDefault<scalar>(
"epsilon", 0.1)),
76 tolerance_(coeffs_.getOrDefault<scalar>(
"tolerance", 1
e-3)),
77 maxIter_(coeffs_.getOrDefault<
int>(
"maxIter", 10)),
98 pdmPredictor_->correct(
y);
107 mesh_.time().timeName(),
115 patchTypes<vector>(mesh_, patchIDs_)
119 volVectorField::Boundary& nybf = ny.boundaryFieldRef();
121 for (
const label patchi : patchIDs_)
123 nybf[patchi] == -
patches[patchi].nf();
127 scalar initialResidual = 0;
132 ny /= (
mag(ny) + SMALL);
135 nf /= (
mag(nf) + SMALL);
149 initialResidual = yEqn.solve().initialResidual();
151 }
while (initialResidual > tolerance_ && ++iter < maxIter_);
Defines the attributes of an object for which implicit objectRegistry management is supported,...
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Calculation of approximate distance to nearest patch for all cells and boundary by solving the Eikona...
static constexpr const zero Zero
Global zero (0)
Calculate the divergence of the given field.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
defineTypeNameAndDebug(advectionDiffusion, 0)
Calculate the matrix for the divergence of the given field and flux.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Calculate the matrix for the laplacian of the field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
bool notNull(const T *ptr)
True if ptr is not a pointer (of type T) to the nullObject.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Specialisation of patchDist for wall distance calculation.
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Calculate the matrix for implicit and explicit sources.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
addToRunTimeSelectionTable(patchDistMethod, advectionDiffusion, dictionary)
Calculate the gradient of the given field.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const dimensionedScalar e
Elementary charge.
const polyBoundaryMesh & patches
A special matrix type and solver, designed for finite volume solutions of scalar equations....
static const GeometricField< vector, fvPatchField, volMesh > & null()
Return a null geometric field.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
const dimensionSet dimless
Dimensionless.