dsmcFields.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "dsmcFields.H"
30#include "volFields.H"
31#include "dictionary.H"
32#include "dsmcCloud.H"
33#include "constants.H"
34#include "stringListOps.H"
36
37using namespace Foam::constant;
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
43namespace functionObjects
44{
46
48 (
52 );
53}
54}
55
56
57// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
58
59namespace Foam
60{
61
62static const word& filteredName
63(
64 const word& baseName,
65 const wordList& names,
66 const string& scopePrefix
67)
68{
69 label idx = names.find(baseName);
70
71 if (idx < 0 && !scopePrefix.empty())
72 {
73 // Take the first matching item
74 idx = firstMatchingString(regExp(scopePrefix + baseName), names);
75 }
76
77 if (idx < 0)
78 {
79 return word::null;
80 }
81
82 return names[idx];
83}
84
85} // End namespace Foam
86
87
88// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
89
91(
92 const word& name,
93 const Time& runTime,
94 const dictionary& dict
95)
96:
98{
99 read(dict);
100}
101
102
103// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104
106{
108 return true;
109}
110
111
113{
114 return true;
115}
116
117
119{
120 // This is fairly horrible with too many hard-coded names...
121
122 // Pre-filter names to obtain 'Mean' vol fields
123 const wordList allMeanNames
124 (
125 obr_.sortedNames
126 (
127 regExp("vol.*Field"), // Any vol field type
128 regExp(".+Mean") // Mean field names
129 )
130 );
131
132 // The separator is often ':', but could be something else.
133 // Replace as first char in [..], so that the regex remains valid,
134 // even if the separator happens to be '-'.
135
136 string scopePrefix = ".+[_:]";
137 scopePrefix[3] = IOobject::scopeSeparator;
138
139
140 // Find scoped/unscoped field name and do lookup.
141 // Short-circuit with message if not found (name or field)
142
143 // Note: currently just find a match without and with a scoping prefix
144 // but could refine to pick the longest name etc, or after finding
145 // the first matching field, use the same prefix for all subsequent fields
146
147 #undef doLocalCode
148 #define doLocalCode(Name, FieldType, Member) \
149 \
150 const FieldType* Member##Ptr = nullptr; \
151 { \
152 const word& fldName = \
153 filteredName(Name, allMeanNames, scopePrefix); \
154 \
155 if (!fldName.empty()) \
156 { \
157 Member##Ptr = obr_.cfindObject<FieldType>(fldName); \
158 } \
159 \
160 if (returnReduce(!Member##Ptr, orOp<bool>())) \
161 { \
162 Log << type() << ' ' << name() << " : no " << Name \
163 << " field found - not calculating\n"; \
164 return false; \
165 } \
166 } \
167 /* Define the const reference */ \
168 const FieldType& Member = *Member##Ptr;
169
170
171 // rhoNMean: always required
172 doLocalCode("rhoNMean", volScalarField, rhoNMean);
173
174 // Also check for division by zero
175 {
176 const scalar minval = min(mag(rhoNMean)).value();
177
178 if (minval <= VSMALL)
179 {
180 Log << type() << ' ' << name()
181 << " : Small value (" << minval << ") in rhoNMean field"
182 << " - not calculating to avoid division by zero" << nl;
183 return false;
184 }
185 }
186
187
188 // The other fields
189
190 doLocalCode("rhoMMean", volScalarField, rhoMMean);
191 doLocalCode("momentumMean", volVectorField, momentumMean);
192 doLocalCode("linearKEMean", volScalarField, linearKEMean);
193 doLocalCode("internalEMean", volScalarField, internalEMean);
194 doLocalCode("iDofMean", volScalarField, iDofMean);
195 doLocalCode("fDMean", volVectorField, fDMean);
196 #undef doLocalCode
197
198
199 //
200 // Everything seem to be okay - can execute
201 //
202 {
203 Log << "Calculating dsmcFields." << endl;
204
205 Log << " Calculating UMean field." << nl;
207 (
209 (
210 "UMean",
211 obr_.time().timeName(),
212 obr_,
214 ),
215 momentumMean/rhoMMean
216 );
217
218 Log << " Calculating translationalT field." << endl;
219 volScalarField translationalT
220 (
222 (
223 "translationalT",
224 obr_.time().timeName(),
225 obr_,
227 ),
228
229 2.0/(3.0*physicoChemical::k.value()*rhoNMean)
230 *(linearKEMean - 0.5*rhoMMean*(UMean & UMean))
231 );
232
233 Log << " Calculating internalT field." << endl;
234 volScalarField internalT
235 (
237 (
238 "internalT",
239 obr_.time().timeName(),
240 obr_,
242 ),
243 (2.0/physicoChemical::k.value())*(internalEMean/iDofMean)
244 );
245
246 Log << " Calculating overallT field." << endl;
247 volScalarField overallT
248 (
250 (
251 "overallT",
252 obr_.time().timeName(),
253 obr_,
255 ),
256 2.0/(physicoChemical::k.value()*(3.0*rhoNMean + iDofMean))
257 *(linearKEMean - 0.5*rhoMMean*(UMean & UMean) + internalEMean)
258 );
259
260 Log << " Calculating pressure field." << endl;
262 (
264 (
265 "p",
266 obr_.time().timeName(),
267 obr_,
269 ),
270 physicoChemical::k.value()*rhoNMean*translationalT
271 );
272
273 volScalarField::Boundary& pBf = p.boundaryFieldRef();
274
275 forAll(mesh_.boundaryMesh(), i)
276 {
277 const polyPatch& patch = mesh_.boundaryMesh()[i];
278
279 if (isA<wallPolyPatch>(patch))
280 {
281 pBf[i] =
282 fDMean.boundaryField()[i]
283 & (patch.faceAreas()/mag(patch.faceAreas()));
284 }
285 }
286
287
288 // Report
289
290 Log << " mag(UMean) max/min : "
291 << max(mag(UMean)).value() << token::SPACE
292 << min(mag(UMean)).value() << nl
293
294 << " translationalT max/min : "
295 << max(translationalT).value() << token::SPACE
296 << min(translationalT).value() << nl
297
298 << " internalT max/min : "
299 << max(internalT).value() << token::SPACE
300 << min(internalT).value() << nl
301
302 << " overallT max/min : "
303 << max(overallT).value() << token::SPACE
304 << min(overallT).value() << nl
305
306 << " p max/min : "
307 << max(p).value() << token::SPACE
308 << min(p).value() << endl;
309
310
311 // Write
312 UMean.write();
313
314 translationalT.write();
315
316 internalT.write();
317
318 overallT.write();
319
320 p.write();
321 }
322
323 Log << "dsmcFields written." << nl << endl;
324 return true;
325}
326
327
328// ************************************************************************* //
#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.
volVectorField UMean(UMeanHeader, mesh)
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static char scopeSeparator
Character for scoping object names (':' or '_')
Definition: IOobject.H:300
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
label find(const T &val, label pos=0) const
Find index of the first occurrence of the value.
Definition: UList.C:212
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const Type & value() const
Return const reference to value.
Abstract base-class for Time/database function objects.
This function object calculates and outputs the intensive fields:
Definition: dsmcFields.H:98
virtual bool execute()
Do nothing.
Definition: dsmcFields.C:112
virtual bool write()
Calculate and write the DSMC fields.
Definition: dsmcFields.C:118
virtual bool read(const dictionary &)
Read the dsmcFields data.
Definition: dsmcFields.C:105
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Computes the magnitude of an input field.
Definition: mag.H:153
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
Wrapper around C++11 regular expressions with some additional prefix-handling. The prefix-handling is...
Definition: regExpCxx.H:83
virtual bool write(const bool valid=true) const
Write using setting from DB.
@ SPACE
Space [isspace].
Definition: token.H:125
A class for handling words, derived from Foam::string.
Definition: word.H:68
static const word null
An empty word.
Definition: word.H:80
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
engineTime & runTime
const dimensionedScalar k
Boltzmann constant.
Different types of constants.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
label firstMatchingString(const UnaryMatchPredicate &matcher, const UList< StringType > &input, const bool invert=false)
Find first list item that matches, -1 on failure.
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
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
regExpCxx regExp
Selection of preferred regular expression implementation.
Definition: regExpFwd.H:43
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
static const word & filteredName(const word &baseName, const wordList &names, const string &scopePrefix)
Definition: dsmcFields.C:63
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Operations on lists of strings.
#define doLocalCode(GeoField)