LagrangianFieldValue.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "LagrangianFieldValue.H"
27 #include "LagrangianFields.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 template<class Type>
37 {
38  Type value;
42 };
43 
44 template<class Type>
46 {
47  return os
48  << vl.value << token::SPACE
49  << vl.proci << token::SPACE
50  << vl.elementi << token::SPACE
51  << vl.position;
52 }
53 
54 template<class Type>
56 {
57  return is
58  >> vl.value
59  >> vl.proci
60  >> vl.elementi
61  >> vl.position;
62 }
63 
64 #define DefineContiguousValueLocationType(Type, nullArg) \
65  template<> \
66  inline bool contiguous<ValueLocation<Type>>() { return true; }
70 #undef DefineContiguousValueLocationType
71 
72 }
73 
74 
75 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
76 
77 namespace Foam
78 {
79 namespace functionObjects
80 {
83  (
87  );
88 }
89 }
90 
91 const Foam::NamedEnum
92 <
94  6
95 >
97 {"sum", "average", "min", "max", "minMag", "maxMag"};
98 
99 
100 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
101 
102 void Foam::functionObjects::LagrangianFieldValue::readCoeffs
103 (
104  const dictionary& dict
105 )
106 {
107  // Read the fields
108  const bool haveFields = dict.found("fields");
109  const bool haveField = dict.found("field");
110  if (haveFields == haveField)
111  {
113  << "keywords fields and field both "
114  << (haveFields ? "" : "un") << "defined in "
115  << "dictionary " << dict.name()
116  << exit(FatalIOError);
117  }
118  else if (haveFields)
119  {
120  dict.lookup("fields") >> fields_;
121  }
122  else if (haveField)
123  {
124  fields_.resize(1);
125  dict.lookup("field") >> fields_.first();
126  }
127 
128  // Read the weight fields
129  const bool haveWeightFields = dict.found("weightFields");
130  const bool haveWeightField = dict.found("weightField");
131  if (haveWeightFields && haveWeightField)
132  {
134  << "keywords weightFields and weightField both "
135  << "defined in dictionary " << dict.name()
136  << exit(FatalIOError);
137  }
138  else if (haveWeightFields)
139  {
140  dict.lookup("weightFields") >> weightFields_;
141  }
142  else if (haveWeightField)
143  {
144  weightFields_.resize(1);
145  dict.lookup("weightField") >> weightFields_.first();
146  }
147  else
148  {
149  // No weights
150  weightFields_.clear();
151  }
152 
153  // Whether or not we are writing the location
154  dict.readIfPresent<Switch>("writeLocation", writeLocation_);
155 
156  resetName(typeName);
157 }
158 
159 
160 template<class Type>
161 void Foam::functionObjects::LagrangianFieldValue::writeName
162 (
163  const word& name
164 )
165 {
166  static const direction nComponents = pTraits<Type>::nComponents;
167 
168  const auto& componentNames = pTraits<Type>::componentNames;
169 
170  if (Pstream::master())
171  {
172  for (direction d = 0; d < nComponents; ++ d)
173  {
174  writeTabbed
175  (
176  file(),
177  word(operationTypeNames_[operation_])
178  + "("
179  + name
180  + (word(componentNames[d]).empty() ? "" : "_")
181  + word(componentNames[d])
182  + ")"
183  );
184  }
185  }
186 }
187 
188 
189 template<class Type, class LocationType>
190 void Foam::functionObjects::LagrangianFieldValue::writeLocationName
191 (
192  const word& name,
193  const word& locationName
194 )
195 {
196  static const direction nComponents = pTraits<Type>::nComponents;
197  static const direction lnComponents = pTraits<LocationType>::nComponents;
198 
199  const auto& componentNames = pTraits<Type>::componentNames;
200  const auto& lcomponentNames = pTraits<LocationType>::componentNames;
201 
202  if (Pstream::master())
203  {
204  for (direction d = 0; d < nComponents; ++ d)
205  {
206  for (direction ld = 0; ld < lnComponents; ++ ld)
207  {
208  writeTabbed
209  (
210  file(),
211  word(operationTypeNames_[operation_])
212  + "("
213  + name
214  + (word(componentNames[d]).empty() ? "" : "_")
215  + word(componentNames[d])
216  + ")"
217  + ":"
218  + locationName
219  + (word(lcomponentNames[ld]).empty() ? "" : "_")
220  + word(lcomponentNames[ld])
221  );
222  }
223  }
224  }
225 }
226 
227 
228 template<class Type>
229 void Foam::functionObjects::LagrangianFieldValue::writeNameAndLocationNames
230 (
231  const word& name
232 )
233 {
234  writeName<Type>(name);
235 
236  if (writeLocation_)
237  {
238  if (Pstream::parRun())
239  {
240  writeLocationName<Type, label>(name, "processor");
241  }
242  writeLocationName<Type, label>(name, "element");
243  writeLocationName<Type, point>(name, "position");
244  }
245 }
246 
247 
248 template<class Type>
249 void Foam::functionObjects::LagrangianFieldValue::writeValue
250 (
251  const Type& value
252 )
253 {
254  if (Pstream::master())
255  {
256  for (direction d = 0; d < pTraits<Type>::nComponents; ++ d)
257  {
258  file() << valueWidth(1) << component(value, d);
259  }
260  }
261 }
262 
263 
264 template<class Type, class LocationType>
265 void Foam::functionObjects::LagrangianFieldValue::writeLocationValue
266 (
267  const FixedList<LocationType, pTraits<Type>::nComponents>& value
268 )
269 {
270  if (Pstream::master())
271  {
272  for (direction d = 0; d < pTraits<Type>::nComponents; ++ d)
273  {
274  writeValue(value[d]);
275  }
276  }
277 }
278 
279 
280 template<class Type, class Op>
281 void Foam::functionObjects::LagrangianFieldValue::writeValueAndLocationValues
282 (
283  const tmp<Field<Type>>& tField,
284  const scalar emptyValue,
285  const Op& op
286 )
287 {
288  const Field<Type>& field = tField();
289 
290  ValueLocation<Type> result;
291 
292  for (direction d = 0; d < pTraits<Type>::nComponents; ++ d)
293  {
294  setComponent(result.value, d) = emptyValue;
295  result.proci[d] = -1;
296  result.elementi[d] = -1;
297  result.position[d] = point::uniform(NaN);
298  }
299 
300  if (field.size())
301  {
302  FixedList<label, pTraits<Type>::nComponents>& ei = result.elementi;
303 
304  ei = Zero;
305 
306  for (label i = 1; i < field.size(); ++ i)
307  {
308  for (direction d = 0; d < pTraits<Type>::nComponents; ++ d)
309  {
310  if (op(component(field[i], d), component(field[ei[d]], d)))
311  {
312  ei[d] = i;
313  }
314  }
315  }
316 
317  for (direction d = 0; d < pTraits<Type>::nComponents; ++ d)
318  {
319  setComponent(result.value, d) = component(field[ei[d]], d);
320  result.proci[d] = Pstream::parRun() ? Pstream::myProcNo() : -1;
321  result.elementi[d] = ei[d];
322  result.position[d] = mesh().position(ei[d]);
323  }
324  }
325 
326  reduce
327  (
328  result,
329  [&op](const ValueLocation<Type>& a, const ValueLocation<Type>& b)
330  {
331  ValueLocation<Type> result;
332 
333  for (direction d = 0; d < pTraits<Type>::nComponents; ++ d)
334  {
335  const ValueLocation<Type>& r =
336  op(component(a.value, d), component(b.value, d)) ? a : b;
337 
338  setComponent(result.value, d) = component(r.value, d);
339  result.proci[d] = r.proci[d];
340  result.elementi[d] = r.elementi[d];
341  result.position[d] = r.position[d];
342  }
343 
344  return result;
345  }
346  );
347 
348  writeValue<Type>(result.value);
349 
350  if (writeLocation_)
351  {
352  if (Pstream::parRun())
353  {
354  writeLocationValue<Type, label>(result.proci);
355  }
356  writeLocationValue<Type, label>(result.elementi);
357  writeLocationValue<Type, point>(result.position);
358  }
359 }
360 
361 
362 template<template<class> class GeoField>
363 bool Foam::functionObjects::LagrangianFieldValue::multiplyWeight
364 (
365  const word& weightFieldName,
366  scalarField& weight
367 ) const
368 {
369  if (!mesh().foundObject<GeoField<scalar>>(weightFieldName)) return false;
370 
371  const GeoField<scalar>& w =
372  mesh().lookupObject<GeoField<scalar>>(weightFieldName);
373 
374  weight *= w;
375 
376  return true;
377 }
378 
379 
380 template<template<class> class GeoField, class Type>
381 bool Foam::functionObjects::LagrangianFieldValue::writeFieldName
382 (
383  const word& fieldName
384 )
385 {
386  if (!mesh().foundObject<GeoField<Type>>(fieldName)) return false;
387 
388  switch (operation_)
389  {
390  case operationType::sum:
392  writeName<Type>(fieldName);
393  break;
394  case operationType::min:
395  case operationType::max:
396  writeNameAndLocationNames<Type>(fieldName);
397  break;
398  case operationType::minMag:
399  case operationType::maxMag:
400  writeNameAndLocationNames<scalar>(fieldName);
401  break;
402  }
403 
404  return true;
405 }
406 
407 
408 template<template<class> class GeoField, class Type>
409 bool Foam::functionObjects::LagrangianFieldValue::writeFieldValue
410 (
411  const scalarField& weight,
412  const word& fieldName
413 )
414 {
415  if (!mesh().foundObject<GeoField<Type>>(fieldName)) return false;
416 
417  const typename GeoField<Type>::FieldType& field =
418  mesh().lookupObject<GeoField<Type>>(fieldName).primitiveField();
419 
420  switch (operation_)
421  {
422  case operationType::sum:
423  writeValue(gSum(weight*field));
424  break;
426  writeValue(gSum(weight*field)/gSum(weight));
427  break;
428  case operationType::min:
429  writeValueAndLocationValues
430  (
431  weight*field,
432  vGreat,
433  lessOp<scalar>()
434  );
435  break;
436  case operationType::max:
437  writeValueAndLocationValues
438  (
439  weight*field,
440  -vGreat,
441  greaterOp<scalar>()
442  );
443  break;
444  case operationType::minMag:
445  writeValueAndLocationValues
446  (
447  mag(weight*field),
448  vGreat,
449  lessOp<scalar>()
450  );
451  break;
452  case operationType::maxMag:
453  writeValueAndLocationValues
454  (
455  mag(weight*field),
456  -vGreat,
457  greaterOp<scalar>()
458  );
459  break;
460  }
461 
462  return true;
463 }
464 
465 
466 void Foam::functionObjects::LagrangianFieldValue::writeFileHeader(const label)
467 {
468  writeHeader(file(), "Lagrangian Sum");
469  writeCommented(file(), "Time");
470 
471  forAll(fields_, fieldi)
472  {
473  #define WRITE_FIELD_NAME(Type, GeoField) \
474  && !writeFieldName<GeoField, Type>(fields_[fieldi])
475 
476  if
477  (
478  true
482  )
483  {
484  cannotFindObject(fields_[fieldi]);
485  }
486 
487  #undef WRITE_COLUMN_HEADER
488  }
489 
490  file().endl();
491 }
492 
493 
494 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
495 
497 (
498  const word& name,
499  const Time& runTime,
500  const dictionary& dict
501 )
502 :
504  logFiles(mesh(), name),
505  fields_(),
506  weightFields_(),
507  operation_(operationTypeNames_.read(dict.lookup("operation"))),
508  writeLocation_(false)
509 {
510  readCoeffs(dict);
511 }
512 
513 
514 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
515 
517 {}
518 
519 
520 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
521 
523 {
525  {
526  readCoeffs(dict);
527  return true;
528  }
529  else
530  {
531  return false;
532  }
533 }
534 
535 
537 {
538  wordList result(fields_);
539  result.append(weightFields_);
540  return result;
541 }
542 
543 
545 {
546  return true;
547 }
548 
549 
551 {
552  logFiles::write();
553 
554  // Filter out operations that don't make sense if there are no elements
555  if (returnReduce(mesh().size(), sumOp<label>()) == 0)
556  {
557  switch (operation_)
558  {
559  case operationType::sum:
560  break;
562  case operationType::min:
563  case operationType::max:
564  case operationType::minMag:
565  case operationType::maxMag:
566  return true;
567  }
568  }
569 
570  if (Pstream::master())
571  {
572  writeTime(file());
573  }
574 
575  // Construct the weights
576  scalarField weight(mesh().size(), 1);
577  forAll(weightFields_, weightFieldi)
578  {
579  const word& weightFieldName = weightFields_[weightFieldi];
580 
581  if
582  (
583  !multiplyWeight<LagrangianField>(weightFieldName, weight)
584  && !multiplyWeight<LagrangianDynamicField>(weightFieldName, weight)
585  && !multiplyWeight<LagrangianInternalField>(weightFieldName, weight)
586  )
587  {
589  << "Weight field " << weightFieldName << " was not found"
590  << exit(FatalError);
591  }
592  }
593 
594  // Write the field values
595  forAll(fields_, fieldi)
596  {
597  const word& fieldName = fields_[fieldi];
598 
599  #define WRITE_FIELD_VALUE(Type, GeoField) \
600  && !writeFieldValue<GeoField, Type>(weight, fieldName)
601 
602  if
603  (
604  true
608  )
609  {
610  cannotFindObject(fields_[fieldi]);
611  }
612 
613  #undef WRITE_COLUMN_VALUE
614  }
615 
616  if (Pstream::master())
617  {
618  file().endl();
619  }
620 
621  return true;
622 }
623 
624 
625 // ************************************************************************* //
#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.
Definition: UList.H:433
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>.
Definition: FixedList.H:78
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void resize(const label)
Alias for setSize(const label)
Definition: ListI.H:138
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
T & first()
Return the first element of the list.
Definition: UListI.H:114
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
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 bool execute()
Execute. Does nothing.
virtual bool read(const dictionary &)
Read parameters.
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...
Definition: logFiles.H:60
virtual void resetName(const word &name)
Reset the list of names to a single name entry.
Definition: logFiles.C:66
virtual bool write()
Write function.
Definition: logFiles.C:173
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.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
volScalarField & b
Definition: createFields.H:27
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
tmp< VolField< Type > > average(const SurfaceField< Type > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
void writeHeader(std::ostream &, const bool isBinary, const std::string &title)
Write header.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
Type gSum(const FieldField< Field, Type > &f)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
label & setComponent(label &l, const direction)
Definition: label.H:86
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)
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
uint8_t direction
Definition: direction.H:45
GeometricField< Type, LagrangianMesh > LagrangianField
#define Op(opName, op)
Definition: ops.H:100
dictionary dict
FixedList< label, pTraits< Type >::nComponents > proci
FixedList< point, pTraits< Type >::nComponents > position
FixedList< label, pTraits< Type >::nComponents > elementi