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-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
37 using namespace Foam::constant;
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 namespace functionObjects
44 {
45  defineTypeNameAndDebug(dsmcFields, 0);
46 
48  (
49  functionObject,
50  dsmcFields,
51  dictionary
52  );
53 }
54 }
55 
56 
57 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
58 
59 namespace Foam
60 {
61 
62 static 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 
90 Foam::functionObjects::dsmcFields::dsmcFields
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 
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 (returnReduceOr(!Member##Ptr)) \
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  (
208  IOobject
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  (
221  IOobject
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  (
236  IOobject
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  (
249  IOobject
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  (
263  IOobject
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 // ************************************************************************* //
Different types of constants.
dictionary dict
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
List of names generated by calling name() for each list item and filtered for matches.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label firstMatchingString(const UnaryMatchPredicate &matcher, const UList< StringType > &input, const bool invert=false)
Find first list item that matches, -1 on failure.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
regExpPosix regExp
Selection of preferred regular expression implementation.
Definition: regExpFwd.H:36
volVectorField UMean(UMeanHeader, mesh)
Operations on lists of strings.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
virtual bool write()
Calculate and write the DSMC fields.
Definition: dsmcFields.C:111
Macros for easy insertion into run-time selection tables.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
virtual bool execute()
Do nothing.
Definition: dsmcFields.C:105
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
Space [isspace].
Definition: token.H:131
static const word null
An empty word.
Definition: word.H:84
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
label find(const T &val) const
Find index of the first occurrence of the value.
Definition: UList.C:173
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
defineTypeNameAndDebug(combustionModel, 0)
static char scopeSeparator
Character for scoping object names (&#39;:&#39; or &#39;_&#39;)
Definition: IOobject.H:338
const dimensionedScalar k
Boltzmann constant.
Wrapper around POSIX extended regular expressions with some additional prefix-handling. The prefix-handling is loosely oriented on PCRE regular expressions and provides a simple means of tuning the expressions.
Definition: regExpPosix.H:80
Nothing to be read.
#define Log
Definition: PDRblock.C:28
const std::string patch
OpenFOAM patch number as a std::string.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
volScalarField & p
#define doLocalCode(Name, FieldType, Member)
virtual bool read(const dictionary &)
Read the dsmcFields data.
Definition: dsmcFields.C:98
static const word & filteredName(const word &baseName, const wordList &names, const string &scopePrefix)
Definition: dsmcFields.C:56
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)