Go to the documentation of this file.
42 namespace patchDistMethods
59 coeffs_(
dict.subDict(
type() +
"Coeffs")),
69 epsilon_(coeffs_.lookupOrDefault<scalar>(
"epsilon", 0.1)),
70 tolerance_(coeffs_.lookupOrDefault<scalar>(
"tolerance", 1
e-3)),
71 maxIter_(coeffs_.lookupOrDefault<
int>(
"maxIter", 10)),
92 pdmPredictor_->correct(
y);
101 mesh_.time().timeName(),
109 patchTypes<vector>(mesh_, patchIDs_)
121 scalar initialResidual = 0;
126 ny /= (
mag(ny) + SMALL);
129 nf /= (
mag(nf) + SMALL);
143 initialResidual = yEqn.
solve().initialResidual();
145 }
while (initialResidual > tolerance_ && ++iter < maxIter_);
advectionDiffusion(const advectionDiffusion &)
Disallow default bitwise copy construct.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Calculation of approximate distance to nearest patch for all cells and boundary by solving the Eikona...
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Calculate the divergence of the given field.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(advectionDiffusion, 0)
Calculate the matrix for the divergence of the given field and flux.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Calculate the matrix for the laplacian of the field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
A list of keyword definitions, which are a keyword followed by any number of values (e....
Specialisation of patchDist for wall distance calculation.
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Calculate the matrix for implicit and explicit sources.
addToRunTimeSelectionTable(patchDistMethod, advectionDiffusion, dictionary)
Calculate the gradient of the given field.
const dimensionedScalar e
Elementary charge.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)