Go to the documentation of this file.
81 oversetPatch_(ptf.oversetPatch_)
93 oversetPatch_(ptf.oversetPatch_)
105 if (this->oversetPatch_.master())
108 const fvMesh&
mesh = this->internalField().mesh();
110 const word& fldName = this->internalField().name();
112 if (&
mesh.lduAddr() != &
mesh.fvMesh::lduAddr())
119 Info<<
"Skipping overset interpolation for solved-for field "
123 else if (!fvSchemes.
found(
"oversetInterpolation"))
126 <<
"Missing required dictionary entry"
127 <<
" 'oversetInterpolation'"
128 <<
". Skipping overset interpolation for field "
131 else if (fvSchemes.
found(
"oversetInterpolationRequired"))
136 if (fvSchemes.
found(
"oversetInterpolationSuppressed"))
139 <<
"Cannot have both dictionary entry"
140 <<
" 'oversetInterpolationSuppresed' and "
141 <<
" 'oversetInterpolationRequired' for field "
144 const dictionary& intDict = fvSchemes.
subDict
146 "oversetInterpolationRequired"
148 if (intDict.found(fldName))
152 Info<<
"Interpolating field " << fldName <<
endl;
162 const_cast<Field<Type>&
>
164 this->primitiveField()
170 Info<<
"Skipping overset interpolation for field "
176 const dictionary* dictPtr
178 fvSchemes.
findDict(
"oversetInterpolationSuppressed")
184 bool skipInterpolate = suppress.found(fldName);
190 || dictPtr->found(fldName);
197 Info<<
"Skipping suppressed overset interpolation"
198 <<
" for field " << fldName <<
endl;
205 Info<<
"Interpolating non-suppressed field " << fldName
210 const_cast<Field<Type>&
>
212 this->primitiveField()
226 this->writeEntry(
"value",
os);
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
A class for handling words, derived from Foam::string.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Ostream & endl(Ostream &os)
A HashTable with keys but without contents that is similar to std::unordered_set.
oversetFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
virtual void write(Ostream &os) const
const oversetFvPatch & oversetPatch_
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
OBJstream os(runTime.globalPath()/outputName)
Mesh data needed to do the Finite Volume discretisation.
errorManipArg< error, int > exit(error &err, const int errNo=1)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Patch for indicating interpolated boundaries (in overset meshes).
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
HashSet< word, Hash< word > > wordHashSet
A HashSet with word keys and string hasher.
#define FatalIOErrorInFunction(ios)
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
#define IOWarningInFunction(ios)
virtual void initEvaluate(const Pstream::commsTypes commsType)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...