64 #define DefineContiguousValueLocationType(Type, nullArg) \
66 inline bool contiguous<ValueLocation<Type>>() { return true; }
70 #undef DefineContiguousValueLocationType
79 namespace functionObjects
97 {
"sum",
"average",
"min",
"max",
"minMag",
"maxMag"};
102 void Foam::functionObjects::LagrangianFieldValue::readCoeffs
104 const dictionary&
dict
108 const bool haveFields =
dict.found(
"fields");
109 const bool haveField =
dict.found(
"field");
110 if (haveFields == haveField)
113 <<
"keywords fields and field both "
114 << (haveFields ?
"" :
"un") <<
"defined in "
115 <<
"dictionary " <<
dict.name()
120 dict.lookup(
"fields") >> fields_;
129 const bool haveWeightFields =
dict.found(
"weightFields");
130 const bool haveWeightField =
dict.found(
"weightField");
131 if (haveWeightFields && haveWeightField)
134 <<
"keywords weightFields and weightField both "
135 <<
"defined in dictionary " <<
dict.name()
138 else if (haveWeightFields)
140 dict.lookup(
"weightFields") >> weightFields_;
142 else if (haveWeightField)
145 dict.lookup(
"weightField") >> weightFields_.
first();
150 weightFields_.
clear();
154 dict.readIfPresent<Switch>(
"writeLocation", writeLocation_);
161 void Foam::functionObjects::LagrangianFieldValue::writeName
166 static const direction nComponents = pTraits<Type>::nComponents;
168 const auto& componentNames = pTraits<Type>::componentNames;
172 for (
direction d = 0; d < nComponents; ++ d)
177 word(operationTypeNames_[operation_])
180 + (word(componentNames[d]).empty() ?
"" :
"_")
181 + word(componentNames[d])
189 template<
class Type,
class LocationType>
190 void Foam::functionObjects::LagrangianFieldValue::writeLocationName
193 const word& locationName
196 static const direction nComponents = pTraits<Type>::nComponents;
197 static const direction lnComponents = pTraits<LocationType>::nComponents;
199 const auto& componentNames = pTraits<Type>::componentNames;
200 const auto& lcomponentNames = pTraits<LocationType>::componentNames;
204 for (
direction d = 0; d < nComponents; ++ d)
206 for (
direction ld = 0; ld < lnComponents; ++ ld)
211 word(operationTypeNames_[operation_])
214 + (word(componentNames[d]).empty() ?
"" :
"_")
215 + word(componentNames[d])
219 + (word(lcomponentNames[ld]).empty() ?
"" :
"_")
220 + word(lcomponentNames[ld])
229 void Foam::functionObjects::LagrangianFieldValue::writeNameAndLocationNames
234 writeName<Type>(
name);
240 writeLocationName<Type, label>(
name,
"processor");
242 writeLocationName<Type, label>(
name,
"element");
243 writeLocationName<Type, point>(
name,
"position");
249 void Foam::functionObjects::LagrangianFieldValue::writeValue
256 for (
direction d = 0; d < pTraits<Type>::nComponents; ++ d)
264 template<
class Type,
class LocationType>
265 void Foam::functionObjects::LagrangianFieldValue::writeLocationValue
267 const FixedList<LocationType, pTraits<Type>::nComponents>& value
272 for (
direction d = 0; d < pTraits<Type>::nComponents; ++ d)
274 writeValue(value[d]);
280 template<
class Type,
class Op>
281 void Foam::functionObjects::LagrangianFieldValue::writeValueAndLocationValues
283 const tmp<Field<Type>>& tField,
284 const scalar emptyValue,
288 const Field<Type>& field = tField();
290 ValueLocation<Type> result;
292 for (
direction d = 0; d < pTraits<Type>::nComponents; ++ d)
295 result.proci[d] = -1;
296 result.elementi[d] = -1;
302 FixedList<label, pTraits<Type>::nComponents>& ei = result.elementi;
306 for (
label i = 1; i < field.size(); ++ i)
308 for (
direction d = 0; d < pTraits<Type>::nComponents; ++ d)
317 for (
direction d = 0; d < pTraits<Type>::nComponents; ++ d)
321 result.elementi[d] = ei[d];
322 result.position[d] =
mesh().position(ei[d]);
329 [&op](
const ValueLocation<Type>& a,
const ValueLocation<Type>&
b)
331 ValueLocation<Type> result;
333 for (
direction d = 0; d < pTraits<Type>::nComponents; ++ d)
335 const ValueLocation<Type>& r =
339 result.proci[d] = r.proci[d];
340 result.elementi[d] = r.elementi[d];
341 result.position[d] = r.position[d];
348 writeValue<Type>(result.value);
354 writeLocationValue<Type, label>(result.proci);
356 writeLocationValue<Type, label>(result.elementi);
357 writeLocationValue<Type, point>(result.position);
362 template<
template<
class>
class GeoField>
363 bool Foam::functionObjects::LagrangianFieldValue::multiplyWeight
365 const word& weightFieldName,
369 if (!
mesh().foundObject<GeoField<scalar>>(weightFieldName))
return false;
371 const GeoField<scalar>& w =
380 template<
template<
class>
class GeoField,
class Type>
381 bool Foam::functionObjects::LagrangianFieldValue::writeFieldName
383 const word& fieldName
386 if (!
mesh().foundObject<GeoField<Type>>(fieldName))
return false;
392 writeName<Type>(fieldName);
396 writeNameAndLocationNames<Type>(fieldName);
398 case operationType::minMag:
399 case operationType::maxMag:
400 writeNameAndLocationNames<scalar>(fieldName);
408 template<
template<
class>
class GeoField,
class Type>
409 bool Foam::functionObjects::LagrangianFieldValue::writeFieldValue
412 const word& fieldName
415 if (!
mesh().foundObject<GeoField<Type>>(fieldName))
return false;
417 const typename GeoField<Type>::FieldType& field =
423 writeValue(
gSum(weight*field));
426 writeValue(
gSum(weight*field)/
gSum(weight));
429 writeValueAndLocationValues
437 writeValueAndLocationValues
444 case operationType::minMag:
445 writeValueAndLocationValues
452 case operationType::maxMag:
453 writeValueAndLocationValues
466 void Foam::functionObjects::LagrangianFieldValue::writeFileHeader(
const label)
469 writeCommented(
file(),
"Time");
473 #define WRITE_FIELD_NAME(Type, GeoField) \
474 && !writeFieldName<GeoField, Type>(fields_[fieldi])
484 cannotFindObject(fields_[fieldi]);
487 #undef WRITE_COLUMN_HEADER
507 operation_(operationTypeNames_.
read(
dict.lookup(
"operation"))),
508 writeLocation_(false)
539 result.
append(weightFields_);
564 case operationType::minMag:
565 case operationType::maxMag:
577 forAll(weightFields_, weightFieldi)
579 const word& weightFieldName = weightFields_[weightFieldi];
583 !multiplyWeight<LagrangianField>(weightFieldName, weight)
584 && !multiplyWeight<LagrangianDynamicField>(weightFieldName, weight)
585 && !multiplyWeight<LagrangianInternalField>(weightFieldName, weight)
589 <<
"Weight field " << weightFieldName <<
" was not found"
597 const word& fieldName = fields_[fieldi];
599 #define WRITE_FIELD_VALUE(Type, GeoField) \
600 && !writeFieldValue<GeoField, Type>(weight, fieldName)
610 cannotFindObject(fields_[fieldi]);
613 #undef WRITE_COLUMN_VALUE
#define WRITE_FIELD_NAME(Type, GeoField)
#define WRITE_FIELD_VALUE(Type, GeoField)
#define DefineContiguousValueLocationType(Type, nullArg)
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A 1D vector of objects of type <T> with a fixed size <Size>.
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
void append(const T &)
Append an element at the end of the list.
void resize(const label)
Alias for setSize(const label)
void clear()
Clear the list, i.e. set size to zero.
Initialise the NamedEnum HashTable from the static list of names.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
T & first()
Return the first element of the list.
static bool master(const label communicator=0)
Am I the master process.
static bool & parRun()
Is this a parallel run?
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Abstract base-class for Time/database functionObjects.
Function to log a single reduced quantity generated from the values in a Lagrangian field; e....
virtual wordList fields() const
Return the list of fields required.
LagrangianFieldValue(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
static const NamedEnum< operationType, 6 > operationTypeNames_
Operation type names.
virtual ~LagrangianFieldValue()
Destructor.
virtual bool execute()
Execute. Does nothing.
virtual bool write()
Write the sum.
virtual bool read(const dictionary &)
Read parameters.
operationType
Operation type enumeration.
Base class for function objects relating to a Lagrangian mesh.
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
virtual void resetName(const word &name)
Reset the list of names to a single name entry.
virtual bool write()
Write function.
virtual bool read(const dictionary &)
Read optional controls.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
A class for handling words, derived from string.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
tmp< VolField< Type > > average(const SurfaceField< Type > &ssf)
Area-weighted average a surfaceField creating a volField.
void writeHeader(std::ostream &, const bool isBinary, const std::string &title)
Write header.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Type gSum(const FieldField< Field, Type > &f)
bool read(const char *, int32_t &)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
label & setComponent(label &l, const direction)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void component(LagrangianPatchField< typename LagrangianPatchField< Type >::cmptType > &sf, const LagrangianPatchField< Type > &f, const direction d)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Istream & operator>>(Istream &, pistonPointEdgeData &)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
DimensionedField< Type, LagrangianMesh > LagrangianInternalField
GeometricField< Type, LagrangianMesh, LagrangianPrimitiveDynamicField > LagrangianDynamicField
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
GeometricField< Type, LagrangianMesh > LagrangianField
FixedList< label, pTraits< Type >::nComponents > proci
FixedList< point, pTraits< Type >::nComponents > position
FixedList< label, pTraits< Type >::nComponents > elementi