variablesSetTemplates.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2007-2019 PCOpt/NTUA
9 Copyright (C) 2013-2019 FOSS GP
10 Copyright (C) 2019 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "localIOdictionary.H"
31#include "FieldField.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37
38// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39
40template<class Type, template<class> class PatchField, class GeoMesh>
41GeometricField<Type, PatchField, GeoMesh>* variablesSet::allocateNamedField
42(
43 const fvMesh& mesh,
44 const IOobject& io,
45 const word& solverName
46)
47{
48 typedef GeometricField<Type, PatchField, GeoMesh> fieldType;
49
50 // Read-in boundary conditions from given IOobject
51 localIOdictionary dict
52 (
53 IOobject
54 (
55 io.name(),
56 io.instance(),
57 io.local(),
58 io.db(),
61 false
62 ),
63 fieldType::typeName
64 );
65 dictionary& bField(dict.subDict("boundaryField"));
66
67 // Add solverName to all patch entries.
68 // Reduntant if not adjoint fields, but overhead should be small
69 for (entry& dEntry : bField)
70 {
71 if (dEntry.isDict())
72 {
73 dEntry.dict().add<word>("solverName", solverName, true);
74 }
75 }
77 << bField << endl;
78
79 return (new fieldType(io, mesh, dict));
80}
81
82
83template<class Type, template<class> class PatchField, class GeoMesh>
84bool variablesSet::readFieldOK
85(
86 autoPtr<GeometricField<Type, PatchField, GeoMesh>>& fieldPtr,
87 const fvMesh& mesh,
88 const word& baseName,
89 const word& solverName,
90 const bool useSolverNameForFields
91)
92{
93 typedef GeometricField<Type, PatchField, GeoMesh> fieldType;
94
95 word customName = baseName + solverName;
96 IOobject headerCustomName
97 (
98 customName,
99 mesh.time().timeName(),
100 mesh,
103 );
104
105 IOobject headerBaseName
106 (
107 baseName,
108 mesh.time().timeName(),
109 mesh,
112 );
113
114 bool fieldFound(false);
115
116 // Read field with full name (i.e. baseName plus solverName) if present
117 if
118 (
119 headerCustomName.typeHeaderOk<fieldType>(false)
121 )
122 {
123 fieldPtr.reset
124 (
125 allocateNamedField<Type, PatchField, GeoMesh>
126 (
127 mesh,
128 headerCustomName,
130 )
131 );
132 fieldFound = true;
133 }
134 // else, see whether the base field exists
135 else if (headerBaseName.typeHeaderOk<fieldType>(false))
136 {
137 fieldPtr.reset
138 (
139 allocateNamedField<Type, PatchField, GeoMesh>
140 (
141 mesh,
142 headerBaseName,
144 )
145 );
146
147 // Rename field if necessary
149 {
150 Info<< "Field " << customName << " not found" << endl;
151 Info<< "Reading base field " << baseName << " and renaming ... "
152 << endl;
153 fieldPtr.ref().rename(customName);
154 }
155 fieldFound = true;
156 }
157
158 return fieldFound;
159}
160
161
162// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
163
164template<class Type, template<class> class PatchField, class GeoMesh>
165autoPtr<GeometricField<Type, PatchField, GeoMesh>>
167(
169)
170{
172 autoPtr<fieldType> returnField(nullptr);
173 if (bf.valid())
174 {
175 const word timeName = bf().mesh().time().timeName();
176 returnField.reset
177 (
178 new fieldType
179 (
180 bf().name() + timeName,
181 bf()
182 )
183 );
184 }
185 return returnField;
186}
187
188
189template<class Type, template<class> class PatchField, class GeoMesh>
191(
194)
195{
196 // Swapping pointers is OK for the mean flow fields known by the
197 // variablesSet (and, in essence, by the solver).
198 // The problem is that turbulence models know references to U and phi
199 // which cannot be swapped.
200 /*
201 const word name1 = p1().name();
202 const word name2 = p2().name();
203 p1.swap(p2);
204
205 p2().rename("temp");
206 p1().rename(name1);
207 p2().rename(name2);
208 */
209
210 // Copy back-up fields to original instead. Slower but there seems to be
211 // no other way
213 p1() == p2();
214 p2() == temp;
215}
216
217
218// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
219
220template<class Type>
222(
224 const fvMesh& mesh,
225 const word& baseName,
226 const word& solverName,
227 const bool useSolverNameForFields
228)
229{
230 // Try to read in field with custom or base name
231 bool fieldFound
232 (
233 readFieldOK
234 (
235 fieldPtr,
236 mesh,
237 baseName,
240 )
241 );
242
243 // No base or custom field found. This is fatal
244 if (!fieldFound)
245 {
247 << "Could not read field with custom ("
248 << word(baseName + solverName) << ") "
249 << "or base (" << baseName << ") name"
250 << exit(FatalError);
251 }
252}
253
254
255template<class Type>
257(
258 const fvMesh& mesh,
259 const word& baseName,
260 const word& solverName,
261 const bool useSolverNameForFields
262)
263{
265
266 autoPtr<VolFieldType> fieldPtr(nullptr);
267 setField(fieldPtr, mesh, baseName, solverName, useSolverNameForFields);
268
269 return tmp<VolFieldType>(fieldPtr.ptr());
270}
271
272
273template<class Type>
275(
277 const word& solverName
278)
279{
280 // typedefs
282 typedef typename VolFieldType::Boundary Boundary;
283
284 // Name of custom field, to be potentially read in
285 const word baseName = baseField.name();
286 const word customName = baseName + solverName;
287 const fvMesh& mesh = baseField.mesh();
288
289 // Renaming of the base field
290 baseField.rename(customName);
291
292 // Create field with baseName and write it, to enable continuation
293 // Note: gives problems for multi-point runs since we end up with
294 // multiple db entries with the same name (one from here and one from
295 // the solver that will construct a turbulenceModel).
296 // Handled through solver.write() for now
297 /*
298 if (!mesh.foundObject<VolFieldType>(baseName))
299 {
300 autoPtr<VolFieldType> baseCopy(new VolFieldType(baseField));
301 baseCopy().IOobject::writeOpt(baseField.writeOpt());
302 baseCopy().rename(baseName);
303 regIOobject::store(baseCopy);
304 }
305 */
306
307 // Check whether a field with the custom name exists, read it in and
308 // set supplied base field to that
309 IOobject headerCustomName
310 (
311 customName,
312 mesh.time().timeName(),
313 mesh,
316 false // Do not register
317 );
318
319 if (headerCustomName.typeHeaderOk<VolFieldType>(true))
320 {
321 Info<< "Reading custom turbulence field " << customName
322 << " and replacing " << baseName << nl << endl;
323 VolFieldType customField(headerCustomName, mesh);
324
325 // Copy internalfield
326 baseField.primitiveFieldRef() = customField.primitiveField();
327
328 // We might apply different boundary conditions per operating point
329 // We need to read them from the custom files and substitute the ones
330 // known by the turbulence model field
331 Boundary& baseBoundary = baseField.boundaryFieldRef();
332 Boundary& customBoundary = customField.boundaryFieldRef();
333 forAll(baseBoundary, patchI)
334 {
335 baseBoundary.set
336 (
337 patchI,
338 customBoundary[patchI].clone(baseField.ref())
339 );
340 }
341 }
342}
343
344
345template<class Type, template<class> class PatchField, class GeoMesh>
347(
349)
350{
352 field == dimensioned<Type>(field.dimensions(), Zero);
353 if (field.nOldTimes())
354 {
355 fieldType& oldTime = field.oldTime();
357 }
358}
359
360
361// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
362
363} // End namespace Foam
364
365// ************************************************************************* //
const Mesh & mesh() const
Return mesh.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:49
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (uses typeFilePath to find file) and check its info.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:500
const fileName & local() const noexcept
Read access to local path component.
Definition: IOobjectI.H:208
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:196
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
void reset(autoPtr< T > &&other) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:117
T * ptr() noexcept
Same as release().
Definition: autoPtr.H:172
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
virtual void rename(const word &newName)
Rename.
Definition: regIOobject.C:417
A class for managing temporary objects.
Definition: tmp.H:65
void swapAndRename(autoPtr< GeometricField< Type, PatchField, GeoMesh > > &p1, autoPtr< GeometricField< Type, PatchField, GeoMesh > > &p2)
Swap autoPtrs and rename managed fields.
autoPtr< GeometricField< Type, PatchField, GeoMesh > > allocateRenamedField(const autoPtr< GeometricField< Type, PatchField, GeoMesh > > &bf)
void renameTurbulenceField(GeometricField< Type, fvPatchField, volMesh > &baseField, const word &solverName)
static void setField(autoPtr< GeometricField< Type, fvPatchField, volMesh > > &fieldPtr, const fvMesh &mesh, const word &baseName, const word &solverName, const bool useSolverNameForFields)
Read vol fields.
const word & solverName() const
Return solver name.
Definition: variablesSet.C:84
bool useSolverNameForFields() const
Append solver name to fields?
Definition: variablesSet.C:90
virtual autoPtr< variablesSet > clone() const
Clone the variablesSet.
Definition: variablesSet.C:75
tmp< GeometricField< Type, fvPatchField, volMesh > > allocateField(const fvMesh &mesh, const word &baseName, const word &solverName, const bool useSolverNameForFields)
static void nullifyField(GeometricField< Type, PatchField, GeoMesh > &fieldPtr)
Nullify field and old times, if present.
A class for handling words, derived from Foam::string.
Definition: word.H:68
rDeltaTY field()
dynamicFvMesh & mesh
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()
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
word timeName
Definition: getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333