derivedFields.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) 2019-2020 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "derivedFields.H"
29#include "volFields.H"
30#include "dictionary.H"
31#include "Time.H"
32#include "mapPolyMesh.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace functionObjects
40{
43}
44}
45
46
47const Foam::Enum
48<
50>
52({
53 { derivedType::NONE , "none" },
54 { derivedType::MASS_FLUX , "rhoU" },
55 { derivedType::TOTAL_PRESSURE , "pTotal" },
56});
57
58
59// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
60
61namespace Foam
62{
63
64static bool calc_rhoU
65(
66 const fvMesh& mesh,
67 const word& derivedName,
68 const scalar rhoRef
69)
70{
71 // rhoU = rho * U
72
73 const auto* rhoPtr = mesh.findObject<volScalarField>("rho");
75
76 volVectorField* result = mesh.getObjectPtr<volVectorField>(derivedName);
77
78 const bool isNew = !result;
79
80 if (!result)
81 {
82 result = new volVectorField
83 (
85 (
86 derivedName,
87 mesh.time().timeName(),
88 mesh,
91 true
92 ),
93 mesh,
95 );
96
97 result->store();
98 }
99
100 if (rhoPtr)
101 {
102 const auto& rho = *rhoPtr;
103
104 *result = (rho * U);
105 }
106 else
107 {
108 const dimensionedScalar rho("rho", dimDensity, rhoRef);
109
110 *result = (rho * U);
111 }
112
113 return isNew;
114}
115
116
117static bool calc_pTotal
118(
119 const fvMesh& mesh,
120 const word& derivedName,
121 const scalar rhoRef
122)
123{
124 // pTotal = p + rho * U^2 / 2
125
126 const auto* rhoPtr = mesh.findObject<volScalarField>("rho");
129
130 volScalarField* result = mesh.getObjectPtr<volScalarField>(derivedName);
131
132 const bool isNew = !result;
133
134 if (!result)
135 {
136 result = new volScalarField
137 (
139 (
140 derivedName,
141 mesh.time().timeName(),
142 mesh,
145 true
146 ),
147 mesh,
149 );
150
151 result->store();
152 }
153
154 if (rhoPtr)
155 {
156 const auto& rho = *rhoPtr;
157
158 *result = (p + 0.5 * rho * magSqr(U));
159 }
160 else
161 {
162 const dimensionedScalar rho("rho", dimDensity, rhoRef);
163
164 *result = (rho * (p + 0.5 * magSqr(U)));
165 }
166
167 return isNew;
168}
169} // End namespace Foam
170
171
172// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
173
175(
176 const word& name,
177 const Time& runTime,
178 const dictionary& dict
179)
180:
182 derivedTypes_(),
183 rhoRef_(1.0)
184{
185 read(dict);
186}
187
188
189// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
190
192{
194
195 rhoRef_ = dict.getOrDefault<scalar>("rhoRef", 1);
196
197 wordList derivedNames(dict.get<wordList>("derived"));
198
199 derivedTypes_.resize(derivedNames.size());
200
201 label nbad = 0, ngood = 0;
202
203 for (const word& key : derivedNames)
204 {
205 derivedTypes_[ngood] = knownNames.lookup(key, derivedType::UNKNOWN);
206
207 switch (derivedTypes_[ngood])
208 {
209 case derivedType::NONE:
210 break;
211
212 case derivedType::UNKNOWN:
213 {
214 derivedNames[nbad++] = key;
215 break;
216 }
217
218 default:
219 {
220 ++ngood;
221 break;
222 }
223 }
224 }
225
226 if (nbad)
227 {
229 << "Ignoring unknown derived names: "
230 << SubList<word>(derivedNames, nbad) << nl;
231 }
232
233 derivedTypes_.resize(ngood);
234
235 // Output the good names
236 forAll(derivedTypes_, i)
237 {
238 derivedNames[i] = knownNames[derivedTypes_[i]];
239 }
240
241 Info<< type() << " derived: "
242 << flatOutput(SubList<word>(derivedNames, ngood)) << nl;
243
244 return true;
245}
246
247
249{
250 Log << type() << " calculating:";
251
252 for (const derivedType category : derivedTypes_)
253 {
254 bool isNew = false;
255
256 switch (category)
257 {
258 case derivedType::MASS_FLUX:
259 {
260 isNew = calc_rhoU(mesh_, knownNames[category], rhoRef_);
261
262 Log << " " << knownNames[category];
263 if (isNew) Log << " (new)";
264 break;
265 }
266
267 case derivedType::TOTAL_PRESSURE:
268 {
269 isNew = calc_pTotal(mesh_, knownNames[category], rhoRef_);
270
271 Log << " " << knownNames[category];
272 if (isNew) Log << " (new)";
273 break;
274 }
275
276 default:
277 break;
278 }
279 }
280
281 Log << nl << endl;
282
283 return true;
284}
285
286
288{
289 for (const derivedType category : derivedTypes_)
290 {
291 switch (category)
292 {
293 case derivedType::NONE:
294 case derivedType::UNKNOWN:
295 break;
296
297 default:
298 {
299 const auto* ioptr =
300 mesh_.cfindObject<regIOobject>(knownNames[category]);
301
302 if (ioptr)
303 {
304 Log << type() << " " << name() << " write:" << nl
305 << " writing field " << ioptr->name() << endl;
306
307 ioptr->write();
308 }
309 break;
310 }
311 }
312 }
313
314 return true;
315}
316
317
319{
320 for (const derivedType category : derivedTypes_)
321 {
322 mesh_.thisDb().checkOut(knownNames[category]);
323 }
324}
325
326
328{
329 if (&mpm.mesh() == &mesh_)
330 {
331 removeDerivedFields();
332 }
333}
334
335
337{
338 if (&m == &mesh_)
339 {
340 removeDerivedFields();
341 }
342}
343
344
345// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
virtual bool read()
Re-read model coefficients if they have changed.
A List obtained as a section of another List.
Definition: SubList.H:70
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base-class for Time/database function objects.
Computes two predefined derived fields, i.e. rhoU, and pTotal, where the defined fields are hard-code...
derivedType
Options for the derived/calculated field type.
virtual bool read(const dictionary &dict)
Read the data.
static const Enum< derivedType > knownNames
Names for derivedType.
void removeDerivedFields()
Remove (checkOut) derived fields from the object registry.
virtual bool execute()
Calculate the derived fields.
virtual bool write()
Write derived fields.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
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
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:363
void movePoints()
Update for new mesh geometry.
void updateMesh()
Update for new mesh topology.
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
const Type & lookupObject(const word &name, const bool recursive=false) const
Type * getObjectPtr(const word &name, const bool recursive=false) const
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:76
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
engineTime & runTime
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
const dimensionSet dimPressure
static bool calc_rhoU(const fvMesh &mesh, const word &derivedName, const scalar rhoRef)
Definition: derivedFields.C:65
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:215
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
const dimensionSet dimDensity
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
static bool calc_pTotal(const fvMesh &mesh, const word &derivedName, const scalar rhoRef)
Info<< "Reading mechanical properties\n"<< endl;IOdictionary mechanicalProperties(IOobject("mechanicalProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));const dictionary &rhoDict(mechanicalProperties.subDict("rho"));word rhoType(rhoDict.get< word >("type"));autoPtr< volScalarField > rhoPtr
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333