48 <<
"Discretisation scheme not specified\n\n"
50 << MeshConstructorTablePtr_->sortedToc()
54 const word schemeName(schemeData);
56 if (surfaceInterpolation::debug || surfaceInterpolationScheme<Type>::debug)
61 auto* ctorPtr = MeshConstructorTable(schemeName);
70 *MeshConstructorTablePtr_
74 return ctorPtr(
mesh, schemeData);
90 <<
"Discretisation scheme not specified"
92 <<
"Valid schemes are :" <<
endl
93 << MeshConstructorTablePtr_->sortedToc()
97 const word schemeName(schemeData);
99 if (surfaceInterpolation::debug || surfaceInterpolationScheme<Type>::debug)
101 InfoInFunction <<
"Discretisation scheme = " << schemeName << endl;
104 auto* ctorPtr = MeshFluxConstructorTable(schemeName);
113 *MeshFluxConstructorTablePtr_
114 ) <<
exit(FatalIOError);
117 return ctorPtr(
mesh, faceFlux, schemeData);
132 if (surfaceInterpolation::debug)
138 <<
" from cells to faces without explicit correction"
159 "interpolate("+vf.
name()+
')',
171 for (label fi=0; fi<P.
size(); fi++)
173 sfi[fi] =
lambda[fi]*vfi[P[fi]] +
y[fi]*vfi[
N[fi]];
206template<
class SFType>
223 if (surfaceInterpolation::debug)
229 <<
" from cells to faces without explicit correction"
251 "interpolate("+vf.
name()+
')',
263 const typename SFType::Internal& Sfi = Sf.internalField();
265 for (label fi=0; fi<P.
size(); fi++)
270 sfi[fi] = Sfi[fi] & (
lambda[fi]*(vfi[P[fi]] - vfi[
N[fi]]) + vfi[
N[fi]]);
281 const typename SFType::Patch& pSf = Sf.boundaryField()[
pi];
337 if (surfaceInterpolation::debug)
343 <<
" from cells to faces"
357 tsf.ref().oriented() = Sf.
oriented();
395 return tSfDotinterpVf;
406 if (surfaceInterpolation::debug)
412 <<
" from cells to faces"
416 tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tsf
421 tsf.ref() += correction(vf);
constexpr scalar pi(M_PI)
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
orientedType oriented() const noexcept
Return oriented type.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Generic GeometricField class.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & name() const noexcept
Return the object name.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
const fileName & instance() const noexcept
Read access to instance path component.
bool eof() const noexcept
True if end of input seen.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
void size(const label n)
Older name for setAddressableSize.
Mesh data needed to do the Finite Volume discretisation.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) -2 >::type type
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
virtual tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &) const =0
Return the interpolation weighting factors for the given field.
static tmp< GeometricField< typename innerProduct< typename SFType::value_type, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const SFType &Sf, const GeometricField< Type, fvPatchField, volMesh > &vf, const tmp< surfaceScalarField > &tlambdas)
Return the face-interpolate of the given cell field.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
Return the face-interpolate of the given cell field.
const fvMesh & mesh() const
Return mesh reference.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the face-interpolate.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
Mesh data needed to do the Finite Volume discretisation.
A class for managing temporary objects.
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
A class for handling words, derived from Foam::string.
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define InfoInFunction
Report an information message using Foam::Info.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
errorManipArg< error, int > exit(error &err, const int errNo=1)
UList< label > labelUList
A UList of labels.
fvsPatchField< scalar > fvsPatchScalarField
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
const Vector< label > N(dict.get< Vector< label > >("N"))