Go to the documentation of this file.
54 ) <<
"Discretisation scheme not specified"
56 <<
"Valid schemes are :" <<
endl
57 << MeshConstructorTablePtr_->sortedToc()
61 const word schemeName(schemeData);
65 Info<<
"surfaceInterpolationScheme<Type>::New"
66 "(const fvMesh&, Istream&)"
67 " : discretisation scheme = "
72 typename MeshConstructorTable::iterator constructorIter =
73 MeshConstructorTablePtr_->find(schemeName);
75 if (constructorIter == MeshConstructorTablePtr_->end())
80 ) <<
"Unknown discretisation scheme "
81 << schemeName <<
nl <<
nl
82 <<
"Valid schemes are :" <<
endl
83 << MeshConstructorTablePtr_->sortedToc()
87 return constructorIter()(
mesh, schemeData);
100 if (schemeData.
eof())
105 ) <<
"Discretisation scheme not specified"
107 <<
"Valid schemes are :" <<
endl
108 << MeshConstructorTablePtr_->sortedToc()
112 const word schemeName(schemeData);
116 Info<<
"surfaceInterpolationScheme<Type>::New"
117 "(const fvMesh&, const surfaceScalarField&, Istream&)"
118 " : discretisation scheme = "
123 typename MeshFluxConstructorTable::iterator constructorIter =
124 MeshFluxConstructorTablePtr_->find(schemeName);
126 if (constructorIter == MeshFluxConstructorTablePtr_->end())
131 ) <<
"Unknown discretisation scheme "
132 << schemeName <<
nl <<
nl
133 <<
"Valid schemes are :" <<
endl
134 << MeshFluxConstructorTablePtr_->sortedToc()
138 return constructorIter()(
mesh, faceFlux, schemeData);
162 if (surfaceInterpolation::debug)
164 Info<<
"surfaceInterpolationScheme<Type>::uncorrectedInterpolate"
165 "(const GeometricField<Type, fvPatchField, volMesh>&, "
166 "const tmp<surfaceScalarField>&, "
167 "const tmp<surfaceScalarField>&) : "
171 <<
" from cells to faces "
172 "without explicit correction"
193 "interpolate("+vf.name()+
')',
207 sfi[fi] =
lambda[fi]*vfi[P[fi]] +
y[fi]*vfi[
N[fi]];
247 if (surfaceInterpolation::debug)
249 Info<<
"surfaceInterpolationScheme<Type>::interpolate"
250 "(const GeometricField<Type, fvPatchField, volMesh>&, "
251 "const tmp<surfaceScalarField>&) : "
255 <<
" from cells to faces "
256 "without explicit correction"
275 "interpolate("+vf.name()+
')',
289 sfi[fi] =
lambda[fi]*(vfi[P[fi]] - vfi[
N[fi]]) + vfi[
N[fi]];
300 tsf().boundaryField()[
pi] =
325 if (surfaceInterpolation::debug)
327 Info<<
"surfaceInterpolationScheme<Type>::interpolate"
328 "(const GeometricField<Type, fvPatchField, volMesh>&) : "
332 <<
" from cells to faces"
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
A class for handling words, derived from string.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
bool eof() const
Return true if end of input seen.
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.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
virtual ~surfaceInterpolationScheme()
Destructor.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Ostream & endl(Ostream &os)
Add newline and flush stream.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
InternalField & internalField()
Return internal field.
const labelUList & neighbour() const
Internal face neighbour.
Mesh data needed to do the Finite Volume discretisation.
volScalarField sf(fieldObject, mesh)
const labelUList & owner() const
Internal face owner.
errorManipArg< error, int > exit(error &err, const int errNo=1)
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Abstract base class for surface interpolation schemes.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
label size() const
Return the number of elements in the UList.
Generic GeometricField class.
void clear() const
If object pointer points to valid object:
dimensionedScalar lambda(laminarTransport.lookup("lambda"))