66 const string& scopePrefix
69 label idx =
names.find(baseName);
71 if (idx < 0 && !scopePrefix.empty())
90Foam::functionObjects::dsmcFields::dsmcFields
136 string scopePrefix =
".+[_:]";
148 #define doLocalCode(Name, FieldType, Member) \
150 const FieldType* Member##Ptr = nullptr; \
152 const word& fldName = \
153 filteredName(Name, allMeanNames, scopePrefix); \
155 if (!fldName.empty()) \
157 Member##Ptr = obr_.cfindObject<FieldType>(fldName); \
160 if (returnReduceOr(!Member##Ptr)) \
162 Log << type() << ' ' << name() << " : no " << Name \
163 << " field found - not calculating\n"; \
168 const FieldType& Member = *Member##Ptr;
176 const scalar minval =
min(
mag(rhoNMean)).value();
178 if (minval <= VSMALL)
181 <<
" : Small value (" << minval <<
") in rhoNMean field"
182 <<
" - not calculating to avoid division by zero" <<
nl;
203 Log <<
"Calculating dsmcFields." <<
endl;
205 Log <<
" Calculating UMean field." <<
nl;
211 obr_.time().timeName(),
215 momentumMean/rhoMMean
218 Log <<
" Calculating translationalT field." <<
endl;
224 obr_.time().timeName(),
230 *(linearKEMean - 0.5*rhoMMean*(
UMean &
UMean))
233 Log <<
" Calculating internalT field." <<
endl;
239 obr_.time().timeName(),
246 Log <<
" Calculating overallT field." <<
endl;
252 obr_.time().timeName(),
257 *(linearKEMean - 0.5*rhoMMean*(
UMean &
UMean) + internalEMean)
260 Log <<
" Calculating pressure field." <<
endl;
266 obr_.time().timeName(),
275 forAll(mesh_.boundaryMesh(), i)
277 const polyPatch& patch = mesh_.boundaryMesh()[i];
282 fDMean.boundaryField()[i]
283 & (patch.faceAreas()/
mag(patch.faceAreas()));
290 Log <<
" mag(UMean) max/min : "
294 <<
" translationalT max/min : "
296 <<
min(translationalT).value() <<
nl
298 <<
" internalT max/min : "
300 <<
min(internalT).value() <<
nl
302 <<
" overallT max/min : "
304 <<
min(overallT).value() <<
nl
314 translationalT.write();
323 Log <<
"dsmcFields written." <<
nl <<
endl;
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)
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
@ NO_READ
Nothing to be read.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static char scopeSeparator
Character for scoping object names (':' or '_').
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Abstract base-class for Time/database function objects.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
This function object calculates and outputs the intensive fields:
virtual bool execute()
Execute the function-object operations, currently does nothing.
virtual bool write()
Write the function-object results.
virtual bool read(const dictionary &)
Read the function-object dictionary.
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
const objectRegistry & obr_
Reference to the region objectRegistry.
A patch is a list of labels that address the faces in the global face list.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
A class for handling words, derived from Foam::string.
static const word null
An empty word.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define doLocalCode(FieldType, Variable)
const dimensionedScalar k
Boltzmann constant.
Different types of constants.
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
List< word > wordList
List of word.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
GeometricField< vector, fvPatchField, volMesh > volVectorField
label firstMatchingString(const UnaryMatchPredicate &matcher, const UList< StringType > &input, const bool invert=false)
Find first list item that matches, -1 on failure.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
regExpCxx regExp
Selection of preferred regular expression implementation.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
constexpr char nl
The newline '\n' character (0x0a).
static const word & filteredName(const word &baseName, const wordList &names, const string &scopePrefix)
#define forAll(list, i)
Loop across all elements in list.
Operations on lists of strings.