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-2022 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"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  namespace functionObjects
42  {
43  defineTypeNameAndDebug(layerAverage, 0);
44  addToRunTimeSelectionTable(functionObject, layerAverage, dictionary);
45  }
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::functionObjects::layerAverage::calcLayers()
52 {
53  // Initialise the faces from which the layers extrude
54  DynamicList<label> startFaces;
55  DynamicList<layerInfo> startFacesInfo;
56  forAll(patchIDs_, i)
57  {
58  const polyPatch& pp = mesh_.boundaryMesh()[patchIDs_[i]];
59  forAll(pp, j)
60  {
61  startFaces.append(pp.start() + j);
62  startFacesInfo.append(layerInfo(0, -1));
63  }
64  }
65  forAll(zoneIDs_, i)
66  {
67  const faceZone& fz = mesh_.faceZones()[zoneIDs_[i]];
68  forAll(fz, j)
69  {
70  startFaces.append(fz[j]);
71  startFacesInfo.append(layerInfo(0, fz.flipMap()[j] ? -1 : +1));
72  }
73  }
74 
75  // Wave to generate layer indices
76  List<layerInfo> allFaceLayerInfo(mesh_.nFaces());
77  List<layerInfo> allCellLayerInfo(mesh_.nCells());
78  FaceCellWave<layerInfo> wave
79  (
80  mesh_,
81  startFaces,
82  startFacesInfo,
83  allFaceLayerInfo,
84  allCellLayerInfo,
85  mesh_.globalData().nTotalCells() + 1
86  );
87 
88  // Copy indices out of the wave and determine the total number of layers
89  nLayers_ = 0;
90  cellLayer_ = labelList(mesh_.nCells(), -1);
91  forAll(cellLayer_, celli)
92  {
93  if (allCellLayerInfo[celli].valid(wave.data()))
94  {
95  const label layeri = allCellLayerInfo[celli].cellLayer();
96  nLayers_ = max(nLayers_, layeri + 1);
97  cellLayer_[celli] = layeri;
98  }
99  }
100  reduce(nLayers_, maxOp<label>());
101 
102  // Report
103  if (nLayers_ != 0)
104  {
105  Info<< " Detected " << nLayers_ << " layers" << nl << endl;
106  }
107  else
108  {
109  WarningInFunction<< "No layers detected" << endl;
110  }
111 
112  // Write the indices for debugging
113  if (debug)
114  {
115  tmp<volScalarField> cellLayer =
117  (
118  "cellLayer",
119  mesh_,
121  );
122  cellLayer.ref().primitiveFieldRef() = List<scalar>(cellLayer_);
123  cellLayer.ref().write();
124  }
125 
126  // Sum number of entries per layer
127  layerCount_ = sum(scalarField(mesh_.nCells(), 1));
128 
129  // Average the cell centres
130  layerCentre_ = sum(mesh_.cellCentres())/layerCount_;
131 
132  // If symmetric, keep only half of the coordinates
133  if (symmetric_)
134  {
135  layerCentre_.setSize(nLayers_/2);
136  }
137 }
138 
139 
140 template<>
142 Foam::functionObjects::layerAverage::symmetricCoeff<Foam::vector>() const
143 {
144  direction i = -1;
145 
146  switch (axis_)
147  {
151  i = label(axis_) - label(coordSet::axisType::X);
152  break;
157  << "Symmetric layer average requested with "
158  << coordSet::axisTypeNames_[axis_] << " axis. Symmetric "
159  << "averaging is only possible along coordinate axes."
160  << exit(FatalError);
161  break;
162  }
163 
164  vector result = vector::one;
165  result[i] = -1;
166  return result;
167 }
168 
169 
170 template<>
172 Foam::functionObjects::layerAverage::symmetricCoeff<Foam::symmTensor>() const
173 {
174  return sqr(symmetricCoeff<vector>());
175 }
176 
177 
178 template<>
180 Foam::functionObjects::layerAverage::symmetricCoeff<Foam::tensor>() const
181 {
182  return symmetricCoeff<symmTensor>();
183 }
184 
185 
186 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
187 
189 (
190  const word& name,
191  const Time& runTime,
192  const dictionary& dict
193 )
194 :
195  fvMeshFunctionObject(name, runTime, dict)
196 {
197  read(dict);
198 }
199 
200 
201 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
202 
204 {}
205 
206 
207 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
208 
210 {
211  Info<< type() << " " << name() << ":" << nl;
212 
213  patchIDs_ =
214  mesh_.boundaryMesh().patchSet
215  (
216  dict.lookupOrDefault<wordReList>("patches", wordReList())
217  ).toc();
218 
219  zoneIDs_ =
221  (
222  dict.lookupOrDefault<wordReList>("zones", wordReList()),
223  mesh_.faceZones().names()
224  );
225 
226  if (patchIDs_.empty() && zoneIDs_.empty())
227  {
229  << "No patches or zones specified" << endl;
230  }
231 
232  symmetric_ = dict.lookupOrDefault<bool>("symmetric", false);
233 
234  axis_ =
236  [
237  dict.lookupOrDefault<word>
238  (
239  "axis",
241  )
242  ];
243 
244  fields_ = dict.lookup<wordList>("fields");
245 
246  formatter_ = setWriter::New(dict.lookup("setFormat"), dict);
247 
248  calcLayers();
249 
250  return true;
251 }
252 
253 
255 {
256  return fields_;
257 }
258 
259 
261 {
262  return true;
263 }
264 
265 
267 {
268  // Create list of available fields
270  forAll(fields_, fieldi)
271  {
272  if
273  (
274  false
275  #define FoundTypeField(Type, nullArg) \
276  || foundObject<VolField<Type>>(fields_[fieldi])
278  #undef FoundTypeField
279  )
280  {
281  fieldNames.append(fields_[fieldi]);
282  }
283  else
284  {
285  cannotFindObject(fields_[fieldi]);
286  }
287  }
288 
289  // Collapse the fields
290  #define DeclareTypeValueSets(Type, nullArg) \
291  PtrList<Field<Type>> Type##ValueSets(fieldNames.size());
293  #undef DeclareTypeValueSets
294  forAll(fieldNames, fieldi)
295  {
296  #define CollapseTypeFields(Type, nullArg) \
297  if (mesh_.foundObject<VolField<Type>>(fieldNames[fieldi])) \
298  { \
299  const VolField<Type>& field = \
300  mesh_.lookupObject<VolField<Type>>(fieldNames[fieldi]); \
301  \
302  Type##ValueSets.set \
303  ( \
304  fieldi, \
305  average(field.primitiveField()) \
306  ); \
307  }
309  #undef CollapseTypeFields
310  }
311 
312  // Output directory
313  fileName outputPath =
314  mesh_.time().globalPath()/writeFile::outputPrefix/name();
315  if (mesh_.name() != fvMesh::defaultRegion)
316  {
317  outputPath = outputPath/mesh_.name();
318  }
319 
320  // Write
321  if (Pstream::master() && layerCentre_.size())
322  {
323  scalarField layerDistance(layerCentre_.size(), 0);
324  for (label i = 1; i < layerCentre_.size(); ++ i)
325  {
326  layerDistance[i] =
327  layerDistance[i-1] + mag(layerCentre_[i] - layerCentre_[i-1]);
328  }
329 
330  formatter_->write
331  (
332  outputPath/mesh_.time().timeName(),
333  typeName,
334  coordSet
335  (
336  identity(layerCentre_.size()),
337  word::null,
338  layerCentre_,
340  layerDistance,
342  ),
343  fieldNames
344  #define TypeValueSetsParameter(Type, nullArg) , Type##ValueSets
346  #undef TypeValueSetsParameter
347  );
348  }
349 
350  return true;
351 }
352 
353 
355 {
356  if (&mesh == &mesh_)
357  {
358  Info<< type() << " " << name() << ":" << nl;
359  calcLayers();
360  }
361 }
362 
363 
365 (
366  const polyTopoChangeMap& map
367 )
368 {
369  if (&map.mesh() == &mesh_)
370  {
371  Info<< type() << " " << name() << ":" << nl;
372  calcLayers();
373  }
374 }
375 
376 
378 {
379  if (&map.mesh() == &mesh_)
380  {
381  Info<< type() << " " << name() << ":" << nl;
382  calcLayers();
383  }
384 }
385 
386 
387 // ************************************************************************* //
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:279
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: layerAverage.C:365
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A class for handling file names.
Definition: fileName.H:79
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: layerAverage.C:377
dimensionedSymmTensor sqr(const dimensionedVector &dv)
uint8_t direction
Definition: direction.H:45
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:325
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Generic GeometricField class.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
const dimensionSet dimless
const Type & data() const
#define DeclareTypeValueSets(Type, nullArg)
fvMesh & mesh
static const NamedEnum< axisType, 6 > axisTypeNames_
String representation of axis enums.
Definition: coordSet.H:69
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
#define TypeValueSetsParameter(Type, nullArg)
virtual bool read(const dictionary &)
Read the field average data.
Definition: layerAverage.C:209
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
static List< word > fieldNames
Definition: globalFoam.H:46
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
bool read(const char *, int32_t &)
Definition: int32IO.C:85
Holds list of sampling positions.
Definition: coordSet.H:50
layerAverage(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: layerAverage.C:189
static const word outputPrefix
Directory prefix.
Definition: writeFile.H:72
A class for handling words, derived from string.
Definition: word.H:59
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
#define FoundTypeField(Type, nullArg)
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
static const word null
An empty word.
Definition: word.H:77
List< label > labelList
A List of labels.
Definition: labelList.H:56
const polyMesh & mesh() const
Return polyMesh.
virtual bool execute()
Do nothing.
Definition: layerAverage.C:260
static const char nl
Definition: Ostream.H:260
#define CollapseTypeFields(Type, nullArg)
virtual wordList fields() const
Return the list of fields required.
Definition: layerAverage.C:254
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
defineTypeNameAndDebug(Qdot, 0)
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: layerAverage.C:354
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
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
virtual bool write()
Calculate and write the graphs.
Definition: layerAverage.C:266
messageStream Info
FOR_ALL_FIELD_TYPES(DefineFvWallInfoType)
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
virtual ~layerAverage()
Destructor.
Definition: layerAverage.C:203
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Templated 3D tensor derived from MatrixSpace adding construction from 9 components, element access using xx(), xy() etc. member functions and the inner-product (dot-product) and outer-product of two Vectors (tensor-product) operators.
Definition: complexI.H:215
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
const polyMesh & mesh() const
Return polyMesh.
Definition: polyMeshMap.H:72
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.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864