layerAverage.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) 2011-2024 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 "layerAverage.H"
27 #include "FaceCellWave.H"
28 #include "layerInfo.H"
29 #include "regionSplit.H"
30 #include "syncTools.H"
31 #include "volFields.H"
32 #include "writeFile.H"
33 #include "polyTopoChangeMap.H"
34 #include "polyMeshMap.H"
35 #include "polyDistributionMap.H"
36 #include "OSspecific.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  namespace functionObjects
44  {
47  }
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 void Foam::functionObjects::layerAverage::calcLayers()
54 {
55  // Initialise the faces from which the layers extrude
56  DynamicList<label> startFaces;
57  DynamicList<layerInfo> startFacesInfo;
58  forAll(patchIndices_, i)
59  {
60  const polyPatch& pp = mesh_.boundaryMesh()[patchIndices_[i]];
61  forAll(pp, j)
62  {
63  startFaces.append(pp.start() + j);
64  startFacesInfo.append(layerInfo(0, -1));
65  }
66  }
67  forAll(zoneIndices_, i)
68  {
69  const faceZone& fz = mesh_.faceZones()[zoneIndices_[i]];
70  forAll(fz, j)
71  {
72  startFaces.append(fz[j]);
73  startFacesInfo.append(layerInfo(0, fz.flipMap()[j] ? -1 : +1));
74  }
75  }
76 
77  // Wave to generate layer indices
78  List<layerInfo> allFaceLayerInfo(mesh_.nFaces());
79  List<layerInfo> allCellLayerInfo(mesh_.nCells());
80  FaceCellWave<layerInfo> wave
81  (
82  mesh_,
83  startFaces,
84  startFacesInfo,
85  allFaceLayerInfo,
86  allCellLayerInfo,
88  );
89 
90  // Copy indices out of the wave and determine the total number of layers
91  nLayers_ = 0;
92  cellLayer_ = labelList(mesh_.nCells(), -1);
93  forAll(cellLayer_, celli)
94  {
95  if (allCellLayerInfo[celli].valid(wave.data()))
96  {
97  const label layeri = allCellLayerInfo[celli].cellLayer();
98  nLayers_ = max(nLayers_, layeri + 1);
99  cellLayer_[celli] = layeri;
100  }
101  }
102  reduce(nLayers_, maxOp<label>());
103 
104  // Report
105  if (nLayers_ != 0)
106  {
107  Info<< " Detected " << nLayers_ << " layers" << nl << endl;
108  }
109  else
110  {
111  WarningInFunction<< "No layers detected" << endl;
112  }
113 
114  // Write the indices for debugging
115  if (debug)
116  {
117  tmp<volScalarField> cellLayer =
119  (
120  "cellLayer",
121  mesh_,
123  );
124  cellLayer.ref().primitiveFieldRef() = List<scalar>(cellLayer_);
125  cellLayer.ref().write();
126  }
127 
128  // Sum number of entries per layer
129  layerVolume_ = sum<scalar>(mesh_.V());
130 
131  // Average the cell centres
132  layerCentre_ = sum<vector>(mesh_.V()*mesh_.C())/layerVolume_;
133 
134  // If symmetric, keep only half of the coordinates
135  if (symmetric_)
136  {
137  layerCentre_.setSize(nLayers_/2);
138  }
139 }
140 
141 
143 Foam::functionObjects::layerAverage::weight() const
144 {
145  if (weightFields_.empty())
146  {
147  return tmp<VolInternalField<scalar>>();
148  }
149 
150  tmp<VolInternalField<scalar>> tresult =
152  (
153  "weight",
154  mesh_,
155  dimensionedScalar(dimless, scalar(1))
156  );
157 
158  forAll(weightFields_, i)
159  {
160  const VolInternalField<scalar>& weightField =
161  mesh_.lookupObject<VolInternalField<scalar>>(weightFields_[i]);
162 
163  tresult.ref() *= weightField;
164  }
165 
166  return tresult;
167 }
168 
169 
170 template<>
172 Foam::functionObjects::layerAverage::symmetricCoeff<Foam::vector>() const
173 {
174  direction i = -1;
175 
176  switch (axis_)
177  {
181  i = label(axis_) - label(coordSet::axisType::X);
182  break;
187  << "Symmetric layer average requested with "
188  << coordSet::axisTypeNames_[axis_] << " axis. Symmetric "
189  << "averaging is only possible along coordinate axes."
190  << exit(FatalError);
191  break;
192  }
193 
194  vector result = vector::one;
195  result[i] = -1;
196  return result;
197 }
198 
199 
200 template<>
202 Foam::functionObjects::layerAverage::symmetricCoeff<Foam::symmTensor>() const
203 {
204  return sqr(symmetricCoeff<vector>());
205 }
206 
207 
208 template<>
210 Foam::functionObjects::layerAverage::symmetricCoeff<Foam::tensor>() const
211 {
212  return symmetricCoeff<symmTensor>();
213 }
214 
215 
216 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
217 
219 (
220  const word& name,
221  const Time& runTime,
222  const dictionary& dict
223 )
224 :
225  fvMeshFunctionObject(name, runTime, dict)
226 {
227  read(dict);
228 }
229 
230 
231 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
232 
234 {}
235 
236 
237 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
238 
240 {
241  Info<< type() << " " << name() << ":" << nl;
242 
243  patchIndices_ = patchSet(dict, true).toc();
244 
245  zoneIndices_ =
247  (
248  dict.lookupOrDefault<wordReList>("zones", wordReList()),
249  mesh_.faceZones().toc()
250  );
251 
252  if (patchIndices_.empty() && zoneIndices_.empty())
253  {
255  << "No patches or zones specified" << endl;
256  }
257 
258  symmetric_ = dict.lookupOrDefault<bool>("symmetric", false);
259 
260  axis_ =
262  [
263  dict.lookupOrDefault<word>
264  (
265  "axis",
267  )
268  ];
269 
270  fields_ = dict.lookup<wordList>("fields");
271 
272  weightFields_ =
273  dict.found("weightFields")
274  ? dict.lookup<wordList>("weightFields")
275  : dict.found("weightField")
276  ? wordList(1, dict.lookup<word>("weightField"))
277  : wordList();
278 
279  formatter_ = setWriter::New(dict.lookup("setFormat"), dict);
280 
281  calcLayers();
282 
283  return true;
284 }
285 
286 
288 {
289  wordList result(fields_);
290  result.append(weightFields_);
291  return result;
292 }
293 
294 
296 {
297  return true;
298 }
299 
300 
302 {
303  // Get the weights
304  tmp<VolInternalField<scalar>> weight(this->weight());
305  tmp<Field<scalar>> layerWeight
306  (
307  weight.valid() ? sum<scalar>(mesh_.V()*weight) : tmp<Field<scalar>>()
308  );
309 
310  // Create list of available fields
312  forAll(fields_, fieldi)
313  {
314  if
315  (
316  false
317  #define FoundTypeField(Type, nullArg) \
318  || foundObject<VolInternalField<Type>>(fields_[fieldi])
320  #undef FoundTypeField
321  )
322  {
323  fieldNames.append(fields_[fieldi]);
324  }
325  else
326  {
327  cannotFindObject(fields_[fieldi]);
328  }
329  }
330 
331  // Collapse the fields
332  #define DeclareTypeValueSets(Type, nullArg) \
333  PtrList<Field<Type>> Type##ValueSets(fieldNames.size());
335  #undef DeclareTypeValueSets
336  forAll(fieldNames, fieldi)
337  {
338  const word& fieldName = fieldNames[fieldi];
339 
340  #define CollapseTypeFields(Type, nullArg) \
341  if (mesh_.foundObject<VolInternalField<Type>>(fieldName)) \
342  { \
343  const VolInternalField<Type>& field = \
344  mesh_.lookupObject<VolInternalField<Type>>(fieldName); \
345  \
346  Type##ValueSets.set \
347  ( \
348  fieldi, \
349  average<Type>(weight, layerWeight, field) \
350  ); \
351  }
353  #undef CollapseTypeFields
354  }
355 
356  // Write
357  if (Pstream::master() && layerCentre_.size())
358  {
359  // Make output directory
360  const fileName outputPath =
361  time_.globalPath()
363  /(mesh_.name() != polyMesh::defaultRegion ? mesh_.name() : word())
364  /name()
365  /time_.name();
366  mkDir(outputPath);
367 
368  scalarField layerDistance(layerCentre_.size(), 0);
369  for (label i = 1; i < layerCentre_.size(); ++ i)
370  {
371  layerDistance[i] =
372  layerDistance[i-1] + mag(layerCentre_[i] - layerCentre_[i-1]);
373  }
374 
375  formatter_->write
376  (
377  outputPath,
378  typeName,
379  coordSet
380  (
381  identityMap(layerCentre_.size()),
382  word::null,
383  layerCentre_,
385  layerDistance,
387  ),
388  fieldNames
389  #define TypeValueSetsParameter(Type, nullArg) , Type##ValueSets
392  );
393  }
394 
395  return true;
396 }
397 
398 
400 {
401  if (&mesh == &mesh_)
402  {
403  Info<< type() << " " << name() << ":" << nl;
404  calcLayers();
405  }
406 }
407 
408 
410 (
411  const polyTopoChangeMap& map
412 )
413 {
414  if (&map.mesh() == &mesh_)
415  {
416  Info<< type() << " " << name() << ":" << nl;
417  calcLayers();
418  }
419 }
420 
421 
423 {
424  if (&map.mesh() == &mesh_)
425  {
426  Info<< type() << " " << name() << ":" << nl;
427  calcLayers();
428  }
429 }
430 
431 
433 (
434  const polyDistributionMap& map
435 )
436 {
437  if (&map.mesh() == &mesh_)
438  {
439  Info<< type() << " " << name() << ":" << nl;
440  calcLayers();
441  }
442 }
443 
444 
445 // ************************************************************************* //
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void setSize(const label)
Reset size of List.
Definition: List.C:281
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static const Form one
Definition: VectorSpace.H:119
Holds list of sampling positions.
Definition: coordSet.H:51
static const NamedEnum< axisType, 6 > axisTypeNames_
String representation of axis enums.
Definition: coordSet.H:69
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
A class for handling file names.
Definition: fileName.H:82
Abstract base-class for Time/database functionObjects.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
virtual wordList fields() const
Return the list of fields required.
Definition: layerAverage.C:287
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: layerAverage.C:410
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: layerAverage.C:433
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: layerAverage.C:422
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: layerAverage.C:399
layerAverage(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: layerAverage.C:219
virtual bool execute()
Do nothing.
Definition: layerAverage.C:295
virtual bool write()
Calculate and write the graphs.
Definition: layerAverage.C:301
virtual ~layerAverage()
Destructor.
Definition: layerAverage.C:233
virtual bool read(const dictionary &)
Read the field average data.
Definition: layerAverage.C:239
static const word outputPrefix
Directory prefix.
Definition: writeFile.H:72
const volVectorField & C() const
Return cell centres.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
label nTotalCells() const
Return total number of cells in decomposed mesh.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
const polyMesh & mesh() const
Return polyMesh.
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
const polyMesh & mesh() const
Return polyMesh.
Definition: polyMeshMap.H:75
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:267
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1515
const faceZoneList & faceZones() const
Return face zones.
Definition: polyMesh.H:446
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const polyMesh & mesh() const
Return polyMesh.
label nCells() const
label nFaces() const
static autoPtr< setWriter > New(const word &writeType, const IOstream::streamFormat writeFormat=IOstream::ASCII, const IOstream::compressionType writeCompression=IOstream::UNCOMPRESSED)
Select given write options.
Definition: setWriter.C:286
A class for managing temporary objects.
Definition: tmp.H:55
bool valid() const
Is this temporary object valid,.
Definition: tmpI.H:167
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
static List< word > fieldNames
Definition: globalFoam.H:46
#define FoundTypeField(Type, nullArg)
#define DeclareTypeValueSets(Type, nullArg)
#define TypeValueSetsParameter(Type, nullArg)
#define CollapseTypeFields(Type, nullArg)
#define WarningInFunction
Report a warning using Foam::Warning.
bool valid(const PtrList< ModelType > &l)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
label wave(const fvMesh &mesh, const List< labelPair > &changedPatchAndFaces, const label nCorrections, GeometricField< scalar, PatchField, GeoMesh > &distance, TrackingData &td, GeometricField< DataType, PatchField, GeoMesh > &... data)
Wave distance (and maybe additional) data from faces. If nCorrections is.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< word > wordList
A List of words.
Definition: fileName.H:54
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
dimensionedSymmTensor sqr(const dimensionedVector &dv)
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
messageStream Info
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensioned< scalar > mag(const dimensioned< Type > &)
typename VolField< Type >::Internal VolInternalField
Definition: volFieldsFwd.H:61
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
static const char nl
Definition: Ostream.H:266
uint8_t direction
Definition: direction.H:45
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict