Go to the documentation of this file.
40 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
41 GeometricField<Type, PatchField, GeoMesh>* variablesSet::allocateNamedField
45 const word& solverName
48 typedef GeometricField<Type, PatchField, GeoMesh> fieldType;
51 localIOdictionary
dict
65 dictionary& bField(
dict.subDict(
"boundaryField"));
69 for (entry& dEntry : bField)
73 dEntry.dict().add<word>(
"solverName",
solverName,
true);
79 return (
new fieldType(io,
mesh,
dict));
83 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
84 bool variablesSet::readFieldOK
86 autoPtr<GeometricField<Type, PatchField, GeoMesh>>& fieldPtr,
89 const word& solverName,
90 const bool useSolverNameForFields
93 typedef GeometricField<Type, PatchField, GeoMesh> fieldType;
96 IOobject headerCustomName
99 mesh.time().timeName(),
105 IOobject headerBaseName
108 mesh.time().timeName(),
114 bool fieldFound(
false);
119 headerCustomName.typeHeaderOk<fieldType>(
false)
125 allocateNamedField<Type, PatchField, GeoMesh>
135 else if (headerBaseName.typeHeaderOk<fieldType>(
false))
139 allocateNamedField<Type, PatchField, GeoMesh>
150 Info<<
"Field " << customName <<
" not found" <<
endl;
151 Info<<
"Reading base field " << baseName <<
" and renaming ... "
153 fieldPtr.ref().rename(customName);
164 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
165 autoPtr<GeometricField<Type, PatchField, GeoMesh>>
189 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
225 const word& baseName,
226 const word& solverName,
227 const bool useSolverNameForFields
239 useSolverNameForFields
247 <<
"Could not read field with custom ("
248 <<
word(baseName + solverName) <<
") "
249 <<
"or base (" << baseName <<
") name"
259 const word& baseName,
260 const word& solverName,
261 const bool useSolverNameForFields
267 setField(fieldPtr,
mesh, baseName, solverName, useSolverNameForFields);
277 const word& solverName
282 typedef typename VolFieldType::Boundary Boundary;
285 const word baseName = baseField.name();
286 const word customName = baseName + solverName;
290 baseField.rename(customName);
319 if (headerCustomName.typeHeaderOk<VolFieldType>(
true))
321 Info<<
"Reading custom turbulence field " << customName
322 <<
" and replacing " << baseName <<
nl <<
endl;
323 VolFieldType customField(headerCustomName,
mesh);
332 Boundary& customBoundary = customField.boundaryFieldRef();
333 forAll(baseBoundary, patchI)
338 customBoundary[patchI].clone(baseField.
ref())
345 template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
353 if (
field.nOldTimes())
bool useSolverNameForFields() const
Append solver name to fields?
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & solverName() const
Return solver name.
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
tmp< GeometricField< Type, fvPatchField, volMesh > > allocateField(const fvMesh &mesh, const word &baseName, const word &solverName, const bool useSolverNameForFields)
static word timeName(const scalar t, const int precision=precision_)
Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));if(U.nOldTimes()){ volVectorField *Uold=&U.oldTime();volScalarField *Kold=&K.oldTime();*Kold==0.5 *magSqr(*Uold);while(Uold->nOldTimes()) { Uold=&Uold-> oldTime()
void swapAndRename(autoPtr< GeometricField< Type, PatchField, GeoMesh >> &p1, autoPtr< GeometricField< Type, PatchField, GeoMesh >> &p2)
Swap autoPtrs and rename managed fields.
Ostream & endl(Ostream &os)
Add newline and flush stream.
#define forAll(list, i)
Loop across all elements in list.
static void nullifyField(GeometricField< Type, PatchField, GeoMesh > &fieldPtr)
Nullify field and old times, if present.
messageStream Info
Information stream (stdout output on master, null elsewhere)
surfacesMesh setField(triSurfaceToAgglom)
T * ptr() noexcept
Same as release().
static void setField(autoPtr< GeometricField< Type, fvPatchField, volMesh >> &fieldPtr, const fvMesh &mesh, const word &baseName, const word &solverName, const bool useSolverNameForFields)
Read vol fields.
Generic dimensioned Type class.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Mesh data needed to do the Finite Volume discretisation.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
autoPtr< GeometricField< Type, PatchField, GeoMesh > > allocateRenamedField(const autoPtr< GeometricField< Type, PatchField, GeoMesh >> &bf)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define DebugInfo
Report an information message using Foam::Info.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
const Time & time() const
Return the top-level database.
void renameTurbulenceField(GeometricField< Type, fvPatchField, volMesh > &baseField, const word &solverName)
Generic GeometricField class.