patchCutLayerAverage.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) 2022-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 "patchCutLayerAverage.H"
27 #include "cutPolyIntegral.H"
28 #include "OSspecific.H"
29 #include "volFields.H"
30 #include "writeFile.H"
31 #include "polyTopoChangeMap.H"
32 #include "polyMeshMap.H"
33 #include "polyDistributionMap.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace functionObjects
41 {
44  (
48  );
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
56 Foam::functionObjects::patchCutLayerAverage::calcNonInterpolatingWeights
57 (
58  const scalarField& pointXs,
59  const scalarField& faceMinXs,
60  const scalarField& faceMaxXs,
61  const labelList& faceMinOrder,
62  const scalarField& plotXs,
63  const bool normalise
64 ) const
65 {
66  const polyPatch& pp = mesh_.boundaryMesh()[patchName_];
67  const faceList& faces = pp.localFaces();
68  const vectorField::subField& faceAreas = pp.faceAreas();
69  const vectorField& faceNormals = pp.faceNormals();
70  const pointField& points = pp.localPoints();
71 
72  // One-value for every point. Needed for the cutPoly routines.
73  const scalarField pointOnes(points.size(), scalar(1));
74 
75  // Generate weights for each face in turn
76  DynamicList<weight> dynWeights(faces.size()*2);
77  label layeri = 0;
78  forAll(faceMinOrder, i)
79  {
80  const label facei = faceMinOrder[i];
81 
82  const scalar a = faceAreas[facei] & faceNormals[facei];
83 
84  // Find the next relevant layer
85  while (faceMinXs[facei] > plotXs[layeri + 1])
86  {
87  layeri ++;
88  }
89 
90  // Loop over all relevant layer intervals
91  label layerj = layeri;
92  while (faceMaxXs[facei] > plotXs[layerj])
93  {
94  const scalar x0 = plotXs[layerj], x1 = plotXs[layerj + 1];
95 
96  // Add a new weight
97  dynWeights.append({facei, layerj, a});
98 
99  // Left interval
100  if (faceMinXs[facei] < x0)
101  {
102  dynWeights.last().value -=
104  (
105  faces[facei],
106  faceAreas[facei],
107  scalar(1),
108  cutPoly::faceCuts(faces[facei], pointXs, x0),
109  points,
110  pointOnes,
111  pointXs,
112  x0,
113  true
114  ).first() & faceNormals[facei];
115  }
116 
117  // Right interval
118  if (faceMaxXs[facei] > x1)
119  {
120  dynWeights.last().value -=
122  (
123  faces[facei],
124  faceAreas[facei],
125  scalar(1),
126  cutPoly::faceCuts(faces[facei], pointXs, x1),
127  points,
128  pointOnes,
129  pointXs,
130  x1,
131  false
132  ).first() & faceNormals[facei];
133  }
134 
135  layerj ++;
136  }
137  }
138 
139  // Transfer to non-dynamic storage
140  List<weight> weights;
141  weights.transfer(dynWeights);
142 
143  // Normalise, if requested
144  if (normalise)
145  {
146  scalarField layerWeightSums(nLayers_, scalar(0));
147  forAll(weights, weighti)
148  {
149  layerWeightSums[weights[weighti].layeri] += weights[weighti].value;
150  }
151 
152  Pstream::listCombineGather(layerWeightSums, plusEqOp<scalar>());
153  Pstream::listCombineScatter(layerWeightSums);
154 
155  forAll(weights, weighti)
156  {
157  weights[weighti].value /=
158  (layerWeightSums[weights[weighti].layeri] + vSmall);
159  }
160  }
161 
162  return weights;
163 }
164 
165 
167 Foam::functionObjects::patchCutLayerAverage::calcInterpolatingWeights
168 (
169  const scalarField& pointXs,
170  const scalarField& faceMinXs,
171  const scalarField& faceMaxXs,
172  const labelList& faceMinOrder,
173  const scalarField& plotXs,
174  const bool normalise
175 ) const
176 {
177  const polyPatch& pp = mesh_.boundaryMesh()[patchName_];
178  const faceList& faces = pp.localFaces();
179  const vectorField::subField& faceAreas = pp.faceAreas();
180  const vectorField& faceNormals = pp.faceNormals();
181  const pointField& points = pp.localPoints();
182 
183  // Values of the (piecewise linear) interpolation basis function. Different
184  // for every interval. Is constantly being set on a sub-set of the
185  // currently relevant points in the loop below.
186  scalarField pointFs(points.size(), scalar(0));
187 
188  // Generate weights for each face in turn
189  DynamicList<weight> dynWeights(faces.size()*2);
190  label layeri = 0;
191  forAll(faceMinOrder, i)
192  {
193  const label facei = faceMinOrder[i];
194 
195  const scalar a = faceAreas[facei] & faceNormals[facei];
196 
197  // Find the next relevant plot point
198  while (faceMinXs[facei] > plotXs[layeri + 1])
199  {
200  layeri ++;
201  }
202 
203  // Loop over all relevant plot points
204  label layerj = layeri;
205  while (layerj == 0 || faceMaxXs[facei] > plotXs[layerj - 1])
206  {
207  const scalar xLeft = layerj > 0 ? plotXs[layerj - 1] : NaN;
208  const scalar xMid = plotXs[layerj];
209  const scalar xRight =
210  layerj < nLayers_ - 1 ? plotXs[layerj + 1] : NaN;
211 
212  // Add a new weight
213  dynWeights.append({facei, layerj, 0});
214 
215  // Left interval
216  if (layerj > 0 && faceMinXs[facei] < xMid)
217  {
218  forAll(faces[facei], facePointi)
219  {
220  const label pointi = faces[facei][facePointi];
221  pointFs[pointi] =
222  (pointXs[pointi] - xLeft)/(xMid - xLeft);
223  }
224 
225  const scalar f = faces[facei].average(points, pointFs);
226 
227  // Add the whole face's contribution
228  dynWeights.last().value += f*a;
229 
230  // Cut off anything before the left point
231  if (faceMinXs[facei] < xLeft)
232  {
233  dynWeights.last().value -=
235  (
236  faces[facei],
237  faceAreas[facei],
238  f,
239  cutPoly::faceCuts(faces[facei], pointXs, xLeft),
240  points,
241  pointFs,
242  pointXs,
243  xLeft,
244  true
245  ).second() & faceNormals[facei];
246  }
247 
248  // Cut off anything after the middle point
249  if (faceMaxXs[facei] > xMid)
250  {
251  dynWeights.last().value -=
253  (
254  faces[facei],
255  faceAreas[facei],
256  f,
257  cutPoly::faceCuts(faces[facei], pointXs, xMid),
258  points,
259  pointFs,
260  pointXs,
261  xMid,
262  false
263  ).second() & faceNormals[facei];
264  }
265  }
266 
267  // Right interval
268  if (layerj < nLayers_ - 1 && faceMaxXs[facei] > xMid)
269  {
270  forAll(faces[facei], facePointi)
271  {
272  const label pointi = faces[facei][facePointi];
273  pointFs[pointi] =
274  (xRight - pointXs[pointi])/(xRight - xMid);
275  }
276 
277  const scalar f = faces[facei].average(points, pointFs);
278 
279  // Add the whole face's contribution
280  dynWeights.last().value += f*a;
281 
282  // Cut off anything before the middle point
283  if (faceMinXs[facei] < xMid)
284  {
285  dynWeights.last().value -=
287  (
288  faces[facei],
289  faceAreas[facei],
290  f,
291  cutPoly::faceCuts(faces[facei], pointXs, xMid),
292  points,
293  pointFs,
294  pointXs,
295  xMid,
296  true
297  ).second() & faceNormals[facei];
298  }
299 
300  // Cut off anything after the right point
301  if (faceMaxXs[facei] > xRight)
302  {
303  dynWeights.last().value -=
305  (
306  faces[facei],
307  faceAreas[facei],
308  f,
309  cutPoly::faceCuts(faces[facei], pointXs, xRight),
310  points,
311  pointFs,
312  pointXs,
313  xRight,
314  false
315  ).second() & faceNormals[facei];
316  }
317  }
318 
319  layerj ++;
320  }
321  }
322 
323  // Transfer to non-dynamic storage
324  List<weight> weights;
325  weights.transfer(dynWeights);
326 
327  // Normalise, if requested. Otherwise, double the weight values on the ends
328  // to account for the fact that these points only have half a basis function
329  // contributing to their sums.
330  if (normalise)
331  {
332  scalarField layerWeightSums(nLayers_, scalar(0));
333  forAll(weights, weighti)
334  {
335  layerWeightSums[weights[weighti].layeri] += weights[weighti].value;
336  }
337 
338  Pstream::listCombineGather(layerWeightSums, plusEqOp<scalar>());
339  Pstream::listCombineScatter(layerWeightSums);
340 
341  forAll(weights, weighti)
342  {
343  weights[weighti].value /=
344  (layerWeightSums[weights[weighti].layeri] + vSmall);
345  }
346  }
347  else
348  {
349  forAll(weights, weighti)
350  {
351  if
352  (
353  weights[weighti].layeri == 0
354  || weights[weighti].layeri == nLayers_ - 1
355  )
356  {
357  weights[weighti].value *= 2;
358  }
359  }
360  }
361 
362  return weights;
363 }
364 
365 
367 Foam::functionObjects::patchCutLayerAverage::calcWeights
368 (
369  const scalarField& pointXs,
370  const scalarField& faceMinXs,
371  const scalarField& faceMaxXs,
372  const labelList& faceMinOrder,
373  const scalarField& plotXs,
374  const bool normalise
375 ) const
376 {
377  return
378  interpolate_
379  ? calcInterpolatingWeights
380  (
381  pointXs,
382  faceMinXs,
383  faceMaxXs,
384  faceMinOrder,
385  plotXs,
386  normalise
387  )
388  : calcNonInterpolatingWeights
389  (
390  pointXs,
391  faceMinXs,
392  faceMaxXs,
393  faceMinOrder,
394  plotXs,
395  normalise
396  );
397 }
398 
399 
400 template<class Type>
402 Foam::functionObjects::patchCutLayerAverage::applyWeights
403 (
404  const List<weight>& weights,
405  const Field<Type>& faceValues
406 ) const
407 {
408  tmp<Field<Type>> tLayerValues(new Field<Type>(nLayers_, Zero));
409 
410  forAll(weights, weighti)
411  {
412  tLayerValues.ref()[weights[weighti].layeri] +=
413  weights[weighti].value*faceValues[weights[weighti].facei];
414  }
415 
416  Pstream::listCombineGather(tLayerValues.ref(), plusEqOp<Type>());
417  Pstream::listCombineScatter(tLayerValues.ref());
418 
419  return tLayerValues;
420 }
421 
422 
423 void Foam::functionObjects::patchCutLayerAverage::initialise()
424 {
425  const polyPatch& pp = mesh_.boundaryMesh()[patchName_];
426  const faceList& faces = pp.localFaces();
427  const pointField& points = pp.localPoints();
428 
429  // If interpolating, then the layers and the plot points are coincident. If
430  // not interpolating, then the layers lie in between the plot points, so
431  // there is one more point than there are layers.
432  const label nPlot = interpolate_ ? nLayers_ : nLayers_ + 1;
433 
434  // Calculate point coordinates
435  const scalarField pointXs(points & direction_);
436 
437  // Determine face min and max coordinates
438  scalarField faceMinXs(faces.size(), vGreat);
439  scalarField faceMaxXs(faces.size(), -vGreat);
440  forAll(faces, facei)
441  {
442  forAll(faces[facei], facePointi)
443  {
444  const label pointi = faces[facei][facePointi];
445  faceMinXs[facei] = min(faceMinXs[facei], pointXs[pointi]);
446  faceMaxXs[facei] = max(faceMaxXs[facei], pointXs[pointi]);
447  }
448  }
449 
450  // Create orderings of the faces based on their min and max coordinates
451  labelList faceMinOrder(faces.size());
452  sortedOrder(faceMinXs, faceMinOrder);
453 
454  // Assume equal spacing to begin with
455  const scalar xMin = gMin(pointXs), xMax = gMax(pointXs);
456  scalarField plotXs
457  (
458  (xMin + scalarList(identityMap(nPlot))/(nPlot - 1)*(xMax - xMin))
459  );
460 
461  // Names and fields for debug output of the counts, to observe the effect
462  // of iterative improvement of the spacing
464  #define DeclareTypeFieldValues(Type, nullArg) \
465  PtrList<Field<Type>> Type##FieldValues;
467  #undef DeclareTypeFieldValues
468 
469  // Iteratively optimise the spacing between the plot points to achieve an
470  // approximately equal number of data points in each interval
471  for (label iteri = 0; iteri < nOptimiseIter_ + debug; ++ iteri)
472  {
473  // Determine the count of faces that contribute to each layer
474  const List<weight> weights =
475  calcWeights
476  (
477  pointXs,
478  faceMinXs,
479  faceMaxXs,
480  faceMinOrder,
481  plotXs,
482  false
483  );
484  const scalarField layerCounts
485  (
486  applyWeights(weights, (1/pp.magFaceAreas())())
487  );
488 
489  if (debug)
490  {
491  const label nFields0 = (2 + !interpolate_)*iteri;
492  const label nFields = (2 + !interpolate_)*(iteri + 1);
493 
494  fieldNames.resize(nFields);
495  #define ResizeTypeFieldValues(Type, nullArg) \
496  Type##FieldValues.resize(nFields);
498  #undef ResizeTypeFieldValues
499 
500  if (!interpolate_)
501  {
502  const SubField<scalar> distance0s(plotXs, nLayers_);
503  const SubField<scalar> distance1s(plotXs, nLayers_, 1);
504 
505  fieldNames[nFields0] = "distance-" + Foam::name(iteri);
506  scalarFieldValues.set(nFields0, (distance0s + distance1s)/2);
507 
508  fieldNames[nFields0 + 1] = "thickness-" + Foam::name(iteri);
509  scalarFieldValues.set(nFields0 + 1, distance1s - distance0s);
510  }
511  else
512  {
513  fieldNames[nFields0] = "distance-" + Foam::name(iteri);
514  scalarFieldValues.set(nFields0, new scalarField(plotXs));
515  }
516 
517  fieldNames[nFields - 1] = "count-" + Foam::name(iteri);
518  scalarFieldValues.set(nFields - 1, new scalarField(layerCounts));
519 
520  if (iteri == nOptimiseIter_) break;
521  }
522 
523  // Do a cumulative sum of the layer counts across all plot points
524  scalarField plotSumCounts(nPlot, 0);
525  for (label ploti = 0; ploti < nPlot - 1; ++ ploti)
526  {
527  plotSumCounts[ploti + 1] =
528  plotSumCounts[ploti]
529  + (
530  interpolate_
531  ? (layerCounts[ploti + 1] + layerCounts[ploti])/2
532  : layerCounts[ploti]
533  );
534  }
535 
536  // Compute the desired count in each interval
537  const scalar plotDeltaCount = plotSumCounts.last()/(nPlot - 1);
538 
539  // Compute the new spacing between the points
540  scalarField plot0Xs(plotXs);
541  plotXs = -vGreat;
542  plotXs.first() = xMin;
543  label ploti = 1;
544  for (label ploti0 = 0; ploti0 < nPlot - 1; ++ ploti0)
545  {
546  while (plotSumCounts[ploti0 + 1] > ploti*plotDeltaCount)
547  {
548  const scalar f =
549  (ploti*plotDeltaCount - plotSumCounts[ploti0])
550  /(plotSumCounts[ploti0 + 1] - plotSumCounts[ploti0]);
551 
552  plotXs[ploti] = (1 - f)*plot0Xs[ploti0] + f*plot0Xs[ploti0 + 1];
553 
554  ploti ++;
555  }
556  }
557  plotXs.last() = xMax;
558  }
559 
560  if (debug)
561  {
562  mkDir(outputPath());
563 
564  formatter_->write
565  (
566  outputPath(),
567  typeName + "_count",
568  coordSet(labelList(nLayers_, 1)),
569  fieldNames
570  #define TypeFieldValuesParameter(Type, nullArg) \
571  , Type##FieldValues
574  );
575  }
576 
577  // Finally, calculate the actual normalised interpolation weights
578  weights_ =
579  calcWeights
580  (
581  pointXs,
582  faceMinXs,
583  faceMaxXs,
584  faceMinOrder,
585  plotXs,
586  true
587  );
588 
589  // Calculate plot coordinates and widths
590  if (interpolate_)
591  {
592  layerDistances_ = plotXs;
593  }
594  else
595  {
596  const SubField<scalar> distance0s(plotXs, nLayers_);
597  const SubField<scalar> distance1s(plotXs, nLayers_, 1);
598  layerDistances_ = (distance0s + distance1s)/2;
599  layerThicknesses_.reset((distance1s - distance0s).ptr());
600  }
601 
602  // Calculate the plot positions
603  layerPositions_ = applyWeights(weights_, pointField(pp.faceCentres()));
604 }
605 
606 
607 Foam::fileName Foam::functionObjects::patchCutLayerAverage::outputPath() const
608 {
609  return
610  time_.globalPath()
612  /(mesh_.name() != polyMesh::defaultRegion ? mesh_.name() : word())
613  /name()
614  /time_.name();
615 }
616 
617 
618 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
619 
621 (
622  const word& name,
623  const Time& runTime,
624  const dictionary& dict
625 )
626 :
627  fvMeshFunctionObject(name, runTime, dict)
628 {
629  read(dict);
630 }
631 
632 
633 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
634 
636 {}
637 
638 
639 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
640 
642 {
643  patchName_ = dict.lookup<word>("patch");
644  direction_ = normalised(dict.lookup<vector>("direction"));
645  nLayers_ = dict.lookup<label>("nPoints");
646 
647  interpolate_ = dict.lookupOrDefault<bool>("interpolate", false);
648 
649  fields_ = dict.lookup<wordList>("fields");
650 
651  axis_ =
653  [
654  dict.lookupOrDefault<word>
655  (
656  "axis",
658  )
659  ];
660 
661  formatter_ = setWriter::New(dict.lookup("setFormat"), dict);
662 
663  nOptimiseIter_ = dict.lookupOrDefault("nOptimiseIter", 2);
664 
665  initialise();
666 
667  return true;
668 }
669 
670 
672 {
673  return fields_;
674 }
675 
676 
678 {
679  return true;
680 }
681 
682 
684 {
685  const polyPatch& pp = mesh_.boundaryMesh()[patchName_];
686 
687  const bool writeThickness =
688  !interpolate_
689  && (
691  || axis_ == coordSet::axisType::DISTANCE
692  );
693 
694  // Create list of field names
696  if (writeThickness)
697  {
698  fieldNames.append("thickness");
699  }
700  forAll(fields_, fieldi)
701  {
702  if
703  (
704  false
705  #define FoundTypeField(Type, nullArg) \
706  || foundObject<VolField<Type>>(fields_[fieldi])
708  #undef FoundTypeField
709  )
710  {
711  fieldNames.append(fields_[fieldi]);
712  }
713  else
714  {
715  cannotFindObject(fields_[fieldi]);
716  }
717  }
718 
719  // Calculate the field values
720  #define DeclareTypeFieldValues(Type, nullArg) \
721  PtrList<Field<Type>> Type##FieldValues(fieldNames.size());
723  #undef DeclareTypeFieldValues
724  if (writeThickness)
725  {
726  scalarFieldValues.set(0, new scalarField(layerThicknesses_));
727  }
728  for (label fieldi = writeThickness; fieldi < fieldNames.size(); ++ fieldi)
729  {
730  #define CollapseTypeFields(Type, nullArg) \
731  if (mesh_.foundObject<VolField<Type>>(fieldNames[fieldi])) \
732  { \
733  const VolField<Type>& field = \
734  mesh_.lookupObject<VolField<Type>>(fieldNames[fieldi]); \
735  \
736  Type##FieldValues.set \
737  ( \
738  fieldi, \
739  applyWeights(weights_, field.boundaryField()[pp.index()]) \
740  ); \
741  }
743  #undef CollapseTypeFields
744  }
745 
746  // Write
747  if (Pstream::master() && layerPositions_.size())
748  {
749  mkDir(outputPath());
750 
751  formatter_->write
752  (
753  outputPath(),
754  typeName,
755  coordSet
756  (
757  identityMap(layerPositions_.size()),
758  word::null,
759  layerPositions_,
761  layerDistances_,
763  ),
764  fieldNames
765  #define TypeFieldValuesParameter(Type, nullArg) , Type##FieldValues
768  );
769  }
770 
771  return true;
772 }
773 
774 
776 (
777  const polyMesh& mesh
778 )
779 {
780  if (&mesh == &mesh_)
781  {
782  initialise();
783  }
784 }
785 
786 
788 (
789  const polyTopoChangeMap& map
790 )
791 {
792  if (&map.mesh() == &mesh_)
793  {
794  initialise();
795  }
796 }
797 
798 
800 (
801  const polyMeshMap& map
802 )
803 {
804  if (&map.mesh() == &mesh_)
805  {
806  initialise();
807  }
808 }
809 
810 
812 (
813  const polyDistributionMap& map
814 )
815 {
816  if (&map.mesh() == &mesh_)
817  {
818  initialise();
819  }
820 }
821 
822 
823 // ************************************************************************* //
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.
Pre-declare SubField and related Field type.
Definition: Field.H:82
SubField< vector > subField
Declare type of subField.
Definition: Field.H:100
Generic GeometricField class.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
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 size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
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
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.
This function object writes graphs of patch face values, area-averaged in planes perpendicular to a g...
patchCutLayerAverage(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
virtual wordList fields() const
Return the list of fields required.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
virtual bool execute()
Execute, currently does nothing.
virtual bool write()
Write the patchCutLayerAverage.
virtual bool read(const dictionary &)
Read the patchCutLayerAverage data.
static const word outputPrefix
Directory prefix.
Definition: writeFile.H:72
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:266
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:403
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:270
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:276
const scalarField::subField magFaceAreas() const
Return face area magnitudes.
Definition: polyPatch.C:288
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const polyMesh & mesh() const
Return polyMesh.
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
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
volScalarField scalarField(fieldObject, mesh)
static List< word > fieldNames
Definition: globalFoam.H:46
const pointField & points
Tuple2< vector, typename outerProduct< vector, Type >::type > faceCutAreaIntegral(const face &f, const vector &fArea, const Type &fPsi, const List< labelPair > &fCuts, const pointField &ps, const Field< Type > &pPsis, const scalarField &pAlphas, const scalar isoAlpha, const bool below)
Compute the face-area and face-area-integral of the given property over.
List< labelPair > faceCuts(const face &f, const scalarField &pAlphas, const scalar isoAlpha)
Return the cuts for a given face. This returns a list of pairs of.
Definition: cutPoly.C:61
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
static Type NaN()
Return a primitive with all components set to NaN.
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
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 mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
FOR_ALL_FIELD_TYPES(DefineContiguousFvWallLocationDataType)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Form normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:413
Field< vector > vectorField
Specialisation of Field<T> for vector.
Type gMin(const FieldField< Field, Type > &f)
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
void sortedOrder(const UList< T > &, labelList &order)
Generate the (stable) sort order for the list.
quaternion normalise(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:603
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
List< face > faceList
Definition: faceListFwd.H:43
Type gMax(const FieldField< Field, Type > &f)
#define TypeFieldValuesParameter(Type, nullArg)
#define FoundTypeField(Type, nullArg)
#define DeclareTypeFieldValues(Type, nullArg)
#define CollapseTypeFields(Type, nullArg)
#define ResizeTypeFieldValues(Type, nullArg)
labelList f(nPoints)
dictionary dict
const scalar xMin
Definition: createFields.H:39
const scalar xMax
Definition: createFields.H:40