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-2023 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(patchIDs_, i)
59  {
60  const polyPatch& pp = mesh_.boundaryMesh()[patchIDs_[i]];
61  forAll(pp, j)
62  {
63  startFaces.append(pp.start() + j);
64  startFacesInfo.append(layerInfo(0, -1));
65  }
66  }
67  forAll(zoneIDs_, i)
68  {
69  const faceZone& fz = mesh_.faceZones()[zoneIDs_[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  layerCount_ = sum(scalarField(mesh_.nCells(), 1));
130 
131  // Average the cell centres
132  layerCentre_ = sum(mesh_.cellCentres())/layerCount_;
133 
134  // If symmetric, keep only half of the coordinates
135  if (symmetric_)
136  {
137  layerCentre_.setSize(nLayers_/2);
138  }
139 }
140 
141 
142 template<>
144 Foam::functionObjects::layerAverage::symmetricCoeff<Foam::vector>() const
145 {
146  direction i = -1;
147 
148  switch (axis_)
149  {
153  i = label(axis_) - label(coordSet::axisType::X);
154  break;
159  << "Symmetric layer average requested with "
160  << coordSet::axisTypeNames_[axis_] << " axis. Symmetric "
161  << "averaging is only possible along coordinate axes."
162  << exit(FatalError);
163  break;
164  }
165 
166  vector result = vector::one;
167  result[i] = -1;
168  return result;
169 }
170 
171 
172 template<>
174 Foam::functionObjects::layerAverage::symmetricCoeff<Foam::symmTensor>() const
175 {
176  return sqr(symmetricCoeff<vector>());
177 }
178 
179 
180 template<>
182 Foam::functionObjects::layerAverage::symmetricCoeff<Foam::tensor>() const
183 {
184  return symmetricCoeff<symmTensor>();
185 }
186 
187 
188 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
189 
191 (
192  const word& name,
193  const Time& runTime,
194  const dictionary& dict
195 )
196 :
197  fvMeshFunctionObject(name, runTime, dict)
198 {
199  read(dict);
200 }
201 
202 
203 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
204 
206 {}
207 
208 
209 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
210 
212 {
213  Info<< type() << " " << name() << ":" << nl;
214 
215  patchIDs_ =
216  mesh_.boundaryMesh().patchSet
217  (
218  dict.lookupOrDefault<wordReList>("patches", wordReList())
219  ).toc();
220 
221  zoneIDs_ =
223  (
224  dict.lookupOrDefault<wordReList>("zones", wordReList()),
225  mesh_.faceZones().names()
226  );
227 
228  if (patchIDs_.empty() && zoneIDs_.empty())
229  {
231  << "No patches or zones specified" << endl;
232  }
233 
234  symmetric_ = dict.lookupOrDefault<bool>("symmetric", false);
235 
236  axis_ =
238  [
239  dict.lookupOrDefault<word>
240  (
241  "axis",
243  )
244  ];
245 
246  fields_ = dict.lookup<wordList>("fields");
247 
248  formatter_ = setWriter::New(dict.lookup("setFormat"), dict);
249 
250  calcLayers();
251 
252  return true;
253 }
254 
255 
257 {
258  return fields_;
259 }
260 
261 
263 {
264  return true;
265 }
266 
267 
269 {
270  // Create list of available fields
272  forAll(fields_, fieldi)
273  {
274  if
275  (
276  false
277  #define FoundTypeField(Type, nullArg) \
278  || foundObject<VolField<Type>>(fields_[fieldi])
280  #undef FoundTypeField
281  )
282  {
283  fieldNames.append(fields_[fieldi]);
284  }
285  else
286  {
287  cannotFindObject(fields_[fieldi]);
288  }
289  }
290 
291  // Collapse the fields
292  #define DeclareTypeValueSets(Type, nullArg) \
293  PtrList<Field<Type>> Type##ValueSets(fieldNames.size());
295  #undef DeclareTypeValueSets
296  forAll(fieldNames, fieldi)
297  {
298  #define CollapseTypeFields(Type, nullArg) \
299  if (mesh_.foundObject<VolField<Type>>(fieldNames[fieldi])) \
300  { \
301  const VolField<Type>& field = \
302  mesh_.lookupObject<VolField<Type>>(fieldNames[fieldi]); \
303  \
304  Type##ValueSets.set \
305  ( \
306  fieldi, \
307  average(field.primitiveField()) \
308  ); \
309  }
311  #undef CollapseTypeFields
312  }
313 
314  // Write
315  if (Pstream::master() && layerCentre_.size())
316  {
317  // Make output directory
318  const fileName outputPath =
319  time_.globalPath()
321  /(mesh_.name() != polyMesh::defaultRegion ? mesh_.name() : word())
322  /name()
323  /time_.name();
324  mkDir(outputPath);
325 
326  scalarField layerDistance(layerCentre_.size(), 0);
327  for (label i = 1; i < layerCentre_.size(); ++ i)
328  {
329  layerDistance[i] =
330  layerDistance[i-1] + mag(layerCentre_[i] - layerCentre_[i-1]);
331  }
332 
333  formatter_->write
334  (
335  outputPath,
336  typeName,
337  coordSet
338  (
339  identityMap(layerCentre_.size()),
340  word::null,
341  layerCentre_,
343  layerDistance,
345  ),
346  fieldNames
347  #define TypeValueSetsParameter(Type, nullArg) , Type##ValueSets
350  );
351  }
352 
353  return true;
354 }
355 
356 
358 {
359  if (&mesh == &mesh_)
360  {
361  Info<< type() << " " << name() << ":" << nl;
362  calcLayers();
363  }
364 }
365 
366 
368 (
369  const polyTopoChangeMap& map
370 )
371 {
372  if (&map.mesh() == &mesh_)
373  {
374  Info<< type() << " " << name() << ":" << nl;
375  calcLayers();
376  }
377 }
378 
379 
381 {
382  if (&map.mesh() == &mesh_)
383  {
384  Info<< type() << " " << name() << ":" << nl;
385  calcLayers();
386  }
387 }
388 
389 
391 (
392  const polyDistributionMap& map
393 )
394 {
395  if (&map.mesh() == &mesh_)
396  {
397  Info<< type() << " " << name() << ":" << nl;
398  calcLayers();
399  }
400 }
401 
402 
403 // ************************************************************************* //
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.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
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:114
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:160
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:256
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: layerAverage.C:368
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: layerAverage.C:391
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: layerAverage.C:380
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: layerAverage.C:357
layerAverage(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: layerAverage.C:191
virtual bool execute()
Do nothing.
Definition: layerAverage.C:262
virtual bool write()
Calculate and write the graphs.
Definition: layerAverage.C:268
virtual ~layerAverage()
Destructor.
Definition: layerAverage.C:205
virtual bool read(const dictionary &)
Read the field average data.
Definition: layerAverage.C:211
static const word outputPrefix
Directory prefix.
Definition: writeFile.H:72
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
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:445
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1563
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:266
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:403
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const polyMesh & mesh() const
Return polyMesh.
label nCells() const
const vectorField & cellCentres() 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 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:306
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)
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
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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:251
const dimensionSet dimless
messageStream Info
FOR_ALL_FIELD_TYPES(DefineContiguousFvWallLocationDataType)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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 > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
uint8_t direction
Definition: direction.H:45
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