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-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 "patchCutLayerAverage.H"
27 #include "cutPolyIntegral.H"
28 #include "OSspecific.H"
29 #include "volPointInterpolation.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 
55 Foam::fileName Foam::functionObjects::patchCutLayerAverage::outputPath() const
56 {
57  return
60  /(mesh_.name() != polyMesh::defaultRegion ? mesh_.name() : word())
61  /name()
62  /time_.name();
63 }
64 
65 
67 Foam::functionObjects::patchCutLayerAverage::calcNonInterpolatingWeights
68 (
69  const scalarField& pointXs,
70  const scalarField& faceMinXs,
71  const scalarField& faceMaxXs,
72  const labelList& faceMinOrder,
73  const scalarField& plotXs,
74  const bool normalise
75 ) const
76 {
77  const polyPatch& pp = mesh_.boundaryMesh()[patchName_];
78  const faceList& faces = pp.localFaces();
79  const vectorField::subField& faceAreas = pp.faceAreas();
80  const vectorField& faceNormals = pp.faceNormals();
81  const pointField& points = pp.localPoints();
82 
83  // Generate weights for each face in turn
84  DynamicList<weight> dynWeights(faces.size()*2);
85  label layeri = 0;
86  forAll(faceMinOrder, i)
87  {
88  const label facei = faceMinOrder[i];
89 
90  const scalar a = faceAreas[facei] & faceNormals[facei];
91 
92  // Find the next relevant layer
93  while (faceMinXs[facei] > plotXs[layeri + 1]) layeri ++;
94 
95  // Loop over all relevant layer intervals
96  label layerj = layeri;
97  while (faceMaxXs[facei] > plotXs[layerj])
98  {
99  // Add a new weight
100  dynWeights.append({facei, layerj, a});
101 
102  // Left interval
103  if (faceMinXs[facei] < plotXs[layerj])
104  {
105  dynWeights.last().value -=
107  (
108  faces[facei],
109  faceAreas[facei],
111  (
112  faces[facei],
113  pointXs,
114  plotXs[layerj]
115  ),
116  points,
117  pointXs,
118  plotXs[layerj],
119  true
120  ) & faceNormals[facei];
121  }
122 
123  // Right interval
124  if (faceMaxXs[facei] > plotXs[layerj + 1])
125  {
126  dynWeights.last().value -=
128  (
129  faces[facei],
130  faceAreas[facei],
132  (
133  faces[facei],
134  pointXs,
135  plotXs[layerj + 1]
136  ),
137  points,
138  pointXs,
139  plotXs[layerj + 1],
140  false
141  ) & faceNormals[facei];
142  }
143 
144  layerj ++;
145  }
146  }
147 
148  // Transfer to non-dynamic storage
149  List<weight> weights;
150  weights.transfer(dynWeights);
151 
152  // Normalise, if requested
153  if (normalise)
154  {
155  scalarField layerWeightSums(nLayers_, scalar(0));
156  forAll(weights, weighti)
157  {
158  layerWeightSums[weights[weighti].layeri] += weights[weighti].value;
159  }
160 
161  Pstream::listCombineGather(layerWeightSums, plusEqOp<scalar>());
162  Pstream::listCombineScatter(layerWeightSums);
163 
164  forAll(weights, weighti)
165  {
166  weights[weighti].value /=
167  (layerWeightSums[weights[weighti].layeri] + vSmall);
168  }
169  }
170 
171  return weights;
172 }
173 
174 
176 Foam::functionObjects::patchCutLayerAverage::calcInterpolatingWeights
177 (
178  const scalarField& pointXs,
179  const scalarField& faceMinXs,
180  const scalarField& faceMaxXs,
181  const labelList& faceMinOrder,
182  const scalarField& plotXs,
183  const bool normalise
184 ) const
185 {
186  const polyPatch& pp = mesh_.boundaryMesh()[patchName_];
187  const faceList& faces = pp.localFaces();
188  const vectorField::subField& faceAreas = pp.faceAreas();
189  const vectorField& faceNormals = pp.faceNormals();
190  const pointField& points = pp.localPoints();
191 
192  // Values of the (piecewise linear) interpolation basis function. Different
193  // for every interval. Is constantly being set on a sub-set of the
194  // currently relevant points in the loop below.
195  scalarField pointFs(points.size(), scalar(0));
196 
197  // Generate weights for each face in turn
198  DynamicList<weight> dynWeights(faces.size()*2);
199  label layeri = 0;
200  forAll(faceMinOrder, i)
201  {
202  const label facei = faceMinOrder[i];
203 
204  const scalar a = faceAreas[facei] & faceNormals[facei];
205 
206  // Find the next relevant layer
207  while (faceMinXs[facei] > plotXs[layeri + 1]) layeri ++;
208 
209  // Loop over all relevant layers
210  label layerj = layeri;
211  while (faceMaxXs[facei] > plotXs[max(layerj - 1, 0)])
212  {
213  // Add a new weight
214  dynWeights.append({facei, layerj, 0});
215 
216  // Left interval
217  if (layerj > 0 && faceMinXs[facei] < plotXs[layerj])
218  {
219  // Update the basis function on the relevant points and
220  // calculate the area-average value on the face
221  forAll(faces[facei], facePointi)
222  {
223  const label pointi = faces[facei][facePointi];
224  pointFs[pointi] =
225  (pointXs[pointi] - plotXs[layerj - 1])
226  /(plotXs[layerj] - plotXs[layerj - 1]);
227  }
228  const scalar faceF =
230  (
231  faces[facei],
232  points,
233  pointFs
234  ).second();
235 
236  // Add the whole face's contribution
237  dynWeights.last().value += faceF*a;
238 
239  // Cut off anything before the left point
240  if (faceMinXs[facei] < plotXs[layerj - 1])
241  {
242  dynWeights.last().value -=
244  (
245  faces[facei],
246  faceAreas[facei],
247  faceF,
249  (
250  faces[facei],
251  pointXs,
252  plotXs[layerj - 1]
253  ),
254  points,
255  pointFs,
256  pointXs,
257  plotXs[layerj - 1],
258  true
259  ).second() & faceNormals[facei];
260  }
261 
262  // Cut off anything after the middle point
263  if (faceMaxXs[facei] > plotXs[layerj])
264  {
265  dynWeights.last().value -=
267  (
268  faces[facei],
269  faceAreas[facei],
270  faceF,
272  (
273  faces[facei],
274  pointXs,
275  plotXs[layerj]
276  ),
277  points,
278  pointFs,
279  pointXs,
280  plotXs[layerj],
281  false
282  ).second() & faceNormals[facei];
283  }
284  }
285 
286  // Right interval
287  if (layerj < nLayers_ - 1 && faceMaxXs[facei] > plotXs[layerj])
288  {
289  // Update the basis function on the relevant points and
290  // calculate the area-average value on the face
291  forAll(faces[facei], facePointi)
292  {
293  const label pointi = faces[facei][facePointi];
294  pointFs[pointi] =
295  (plotXs[layerj + 1] - pointXs[pointi])
296  /(plotXs[layerj + 1] - plotXs[layerj]);
297  }
298  const scalar faceF =
300  (
301  faces[facei],
302  points,
303  pointFs
304  ).second();
305 
306  // Add the whole face's contribution
307  dynWeights.last().value += faceF*a;
308 
309  // Cut off anything before the middle point
310  if (faceMinXs[facei] < plotXs[layerj])
311  {
312  dynWeights.last().value -=
314  (
315  faces[facei],
316  faceAreas[facei],
317  faceF,
319  (
320  faces[facei],
321  pointXs,
322  plotXs[layerj]
323  ),
324  points,
325  pointFs,
326  pointXs,
327  plotXs[layerj],
328  true
329  ).second() & faceNormals[facei];
330  }
331 
332  // Cut off anything after the right point
333  if (faceMaxXs[facei] > plotXs[layerj + 1])
334  {
335  dynWeights.last().value -=
337  (
338  faces[facei],
339  faceAreas[facei],
340  faceF,
342  (
343  faces[facei],
344  pointXs,
345  plotXs[layerj + 1]
346  ),
347  points,
348  pointFs,
349  pointXs,
350  plotXs[layerj + 1],
351  false
352  ).second() & faceNormals[facei];
353  }
354  }
355 
356  layerj ++;
357  }
358  }
359 
360  // Transfer to non-dynamic storage
361  List<weight> weights;
362  weights.transfer(dynWeights);
363 
364  // Normalise, if requested. Otherwise, double the weight values on the ends
365  // to account for the fact that these points only have half a basis function
366  // contributing to their sums.
367  if (normalise)
368  {
369  scalarField layerWeightSums(nLayers_, scalar(0));
370  forAll(weights, weighti)
371  {
372  layerWeightSums[weights[weighti].layeri] += weights[weighti].value;
373  }
374 
375  Pstream::listCombineGather(layerWeightSums, plusEqOp<scalar>());
376  Pstream::listCombineScatter(layerWeightSums);
377 
378  forAll(weights, weighti)
379  {
380  weights[weighti].value /=
381  (layerWeightSums[weights[weighti].layeri] + vSmall);
382  }
383  }
384  else
385  {
386  forAll(weights, weighti)
387  {
388  if
389  (
390  weights[weighti].layeri == 0
391  || weights[weighti].layeri == nLayers_ - 1
392  )
393  {
394  weights[weighti].value *= 2;
395  }
396  }
397  }
398 
399  return weights;
400 }
401 
402 
404 Foam::functionObjects::patchCutLayerAverage::calcWeights
405 (
406  const scalarField& pointXs,
407  const scalarField& faceMinXs,
408  const scalarField& faceMaxXs,
409  const labelList& faceMinOrder,
410  const scalarField& plotXs,
411  const bool normalise
412 ) const
413 {
414  return
415  interpolate_
416  ? calcInterpolatingWeights
417  (
418  pointXs,
419  faceMinXs,
420  faceMaxXs,
421  faceMinOrder,
422  plotXs,
423  normalise
424  )
425  : calcNonInterpolatingWeights
426  (
427  pointXs,
428  faceMinXs,
429  faceMaxXs,
430  faceMinOrder,
431  plotXs,
432  normalise
433  );
434 }
435 
436 
437 void Foam::functionObjects::patchCutLayerAverage::calcWeights()
438 {
439  const polyPatch& pp = mesh_.boundaryMesh()[patchName_];
440  const faceList& faces = pp.localFaces();
441  const pointField& points = pp.localPoints();
442 
443  // If interpolating, then the layers and the plot points are coincident. If
444  // not interpolating, then the layers lie in between the plot points, so
445  // there is one more point than there are layers.
446  const label nPlot = interpolate_ ? nLayers_ : nLayers_ + 1;
447 
448  // Calculate or get the point coordinates
449  tmp<scalarField> tpointXs =
450  distanceName_ == word::null
451  ? points & direction_
452  : (
453  mesh_.foundObject<volScalarField>(distanceName_)
455  (
456  mesh_.lookupObject<volScalarField>(distanceName_)
457  )
458  : tmp<pointScalarField>
459  (
460  mesh_.lookupObject<pointScalarField>(distanceName_)
461  )
462  )().boundaryField()[pp.index()].patchInternalField();
463  const scalarField& pointXs = tpointXs();
464 
465  // Determine face min and max coordinates
466  scalarField faceMinXs(faces.size(), vGreat);
467  scalarField faceMaxXs(faces.size(), -vGreat);
468  forAll(faces, facei)
469  {
470  forAll(faces[facei], facePointi)
471  {
472  const label pointi = faces[facei][facePointi];
473  faceMinXs[facei] = min(faceMinXs[facei], pointXs[pointi]);
474  faceMaxXs[facei] = max(faceMaxXs[facei], pointXs[pointi]);
475  }
476  }
477 
478  // Create orderings of the faces based on their min and max coordinates
479  labelList faceMinOrder(faces.size());
480  sortedOrder(faceMinXs, faceMinOrder);
481 
482  // Assume equal spacing to begin with
483  const scalar xMin = gMin(pointXs), xMax = gMax(pointXs);
484  scalarField plotXs
485  (
486  (xMin + scalarList(identityMap(nPlot))/(nPlot - 1)*(xMax - xMin))
487  );
488 
489  // Names and fields for debug output of the counts, to observe the effect
490  // of iterative improvement of the spacing
492  #define DeclareTypeFieldValues(Type, nullArg) \
493  PtrList<Field<Type>> Type##FieldValues;
495  #undef DeclareTypeFieldValues
496 
497  // Iteratively optimise the spacing between the plot points to achieve an
498  // approximately equal number of data points in each interval
499  for (label iteri = 0; iteri < nOptimiseIter_ + debug; ++ iteri)
500  {
501  // Determine the count of faces that contribute to each layer
502  const List<weight> weights =
503  calcWeights
504  (
505  pointXs,
506  faceMinXs,
507  faceMaxXs,
508  faceMinOrder,
509  plotXs,
510  false
511  );
512  const scalarField layerCounts
513  (
514  applyWeights(weights, (1/pp.magFaceAreas())())
515  );
516 
517  if (debug)
518  {
519  const label nFields0 = (2 + !interpolate_)*iteri;
520  const label nFields = (2 + !interpolate_)*(iteri + 1);
521 
522  fieldNames.resize(nFields);
523  #define ResizeTypeFieldValues(Type, nullArg) \
524  Type##FieldValues.resize(nFields);
526  #undef ResizeTypeFieldValues
527 
528  if (!interpolate_)
529  {
530  const SubField<scalar> distance0s(plotXs, nLayers_);
531  const SubField<scalar> distance1s(plotXs, nLayers_, 1);
532 
533  fieldNames[nFields0] = "distance-" + Foam::name(iteri);
534  scalarFieldValues.set(nFields0, (distance0s + distance1s)/2);
535 
536  fieldNames[nFields0 + 1] = "thickness-" + Foam::name(iteri);
537  scalarFieldValues.set(nFields0 + 1, distance1s - distance0s);
538  }
539  else
540  {
541  fieldNames[nFields0] = "distance-" + Foam::name(iteri);
542  scalarFieldValues.set(nFields0, new scalarField(plotXs));
543  }
544 
545  fieldNames[nFields - 1] = "count-" + Foam::name(iteri);
546  scalarFieldValues.set(nFields - 1, new scalarField(layerCounts));
547 
548  if (iteri == nOptimiseIter_) break;
549  }
550 
551  // Do a cumulative sum of the layer counts across all plot points
552  scalarField plotSumCounts(nPlot, 0);
553  for (label ploti = 0; ploti < nPlot - 1; ++ ploti)
554  {
555  plotSumCounts[ploti + 1] =
556  plotSumCounts[ploti]
557  + (
558  interpolate_
559  ? (layerCounts[ploti + 1] + layerCounts[ploti])/2
560  : layerCounts[ploti]
561  );
562  }
563 
564  // Compute the desired count in each interval
565  const scalar plotDeltaCount = plotSumCounts.last()/(nPlot - 1);
566 
567  // Compute the new spacing between the points
568  scalarField plot0Xs(plotXs);
569  plotXs = -vGreat;
570  plotXs.first() = xMin;
571  label ploti = 1;
572  for (label ploti0 = 0; ploti0 < nPlot - 1; ++ ploti0)
573  {
574  while
575  (
576  ploti < nPlot
577  && plotSumCounts[ploti0 + 1] > ploti*plotDeltaCount
578  )
579  {
580  const scalar f =
581  (ploti*plotDeltaCount - plotSumCounts[ploti0])
582  /(plotSumCounts[ploti0 + 1] - plotSumCounts[ploti0]);
583 
584  plotXs[ploti] = (1 - f)*plot0Xs[ploti0] + f*plot0Xs[ploti0 + 1];
585 
586  ploti ++;
587  }
588  }
589  plotXs.last() = xMax;
590  }
591 
592  if (debug)
593  {
594  mkDir(outputPath());
595 
596  formatter_->write
597  (
598  outputPath(),
599  typeName + "_count",
600  coordSet(labelList(nLayers_, 1)),
601  fieldNames
602  #define TypeFieldValuesParameter(Type, nullArg) \
603  , Type##FieldValues
606  );
607  }
608 
609  // Finally, calculate the actual normalised interpolation weights
610  weights_.reset
611  (
612  new List<weight>
613  (
614  calcWeights
615  (
616  pointXs,
617  faceMinXs,
618  faceMaxXs,
619  faceMinOrder,
620  plotXs,
621  true
622  )
623  )
624  );
625 
626  // Calculate plot coordinates and widths
627  if (interpolate_)
628  {
629  layerDistances_.reset(new scalarField(plotXs));
630  }
631  else
632  {
633  const SubField<scalar> distance0s(plotXs, nLayers_);
634  const SubField<scalar> distance1s(plotXs, nLayers_, 1);
635  layerDistances_.reset(((distance0s + distance1s)/2).ptr());
636  layerThicknesses_.reset((distance1s - distance0s).ptr());
637  }
638 
639  // Calculate the plot positions
640  layerPositions_.reset
641  (
642  applyWeights(weights_, pointField(pp.faceCentres())).ptr()
643  );
644 
645  if (debug)
646  {
647  const List<weight> weights =
648  calcWeights
649  (
650  pointXs,
651  faceMinXs,
652  faceMaxXs,
653  faceMinOrder,
654  plotXs,
655  false
656  );
657 
658  volTensorField layers
659  (
660  IOobject
661  (
662  name() + ":layers",
663  mesh_.time().name(),
664  mesh()
665  ),
666  mesh(),
668  );
669 
670  tensorField& pLayers =
671  layers.boundaryFieldRef()[pp.index()];
672 
673  forAll(weights, weighti)
674  {
675  const weight& w = weights[weighti];
676 
677  pLayers[w.facei][w.layeri % tensor::nComponents] =
678  w.value/pp.magFaceAreas()[w.facei];
679  }
680 
681  Info<< name() << ": Writing " << layers.name() << endl;
682 
683  layers.write();
684  }
685 }
686 
687 
688 template<class Type>
690 Foam::functionObjects::patchCutLayerAverage::applyWeights
691 (
692  const List<weight>& weights,
693  const Field<Type>& faceValues
694 ) const
695 {
696  tmp<Field<Type>> tLayerValues(new Field<Type>(nLayers_, Zero));
697 
698  forAll(weights, weighti)
699  {
700  tLayerValues.ref()[weights[weighti].layeri] +=
701  weights[weighti].value*faceValues[weights[weighti].facei];
702  }
703 
704  Pstream::listCombineGather(tLayerValues.ref(), plusEqOp<Type>());
705  Pstream::listCombineScatter(tLayerValues.ref());
706 
707  return tLayerValues;
708 }
709 
710 
711 void Foam::functionObjects::patchCutLayerAverage::clear()
712 {
713  weights_.clear();
714  layerDistances_.clear();
715  layerThicknesses_.clear();
716  layerPositions_.clear();
717 }
718 
719 
720 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
721 
723 (
724  const word& name,
725  const Time& runTime,
726  const dictionary& dict
727 )
728 :
729  fvMeshFunctionObject(name, runTime, dict)
730 {
731  read(dict);
732 }
733 
734 
735 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
736 
738 {}
739 
740 
741 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
742 
744 {
745  patchName_ = dict.lookup<word>("patch");
746 
747  const bool haveDirection = dict.found("direction");
748  const bool haveDistance = dict.found("distance");
749  if (haveDirection == haveDistance)
750  {
752  << "keywords direction and distance both "
753  << (haveDirection ? "" : "un") << "defined in "
754  << "dictionary " << dict.name()
755  << exit(FatalIOError);
756  }
757  else if (haveDirection)
758  {
759  direction_ = normalised(dict.lookup<vector>("direction"));
760  distanceName_ = word::null;
761  }
762  else if (haveDistance)
763  {
764  direction_ = vector::nan;
765  distanceName_ = dict.lookup<word>("distance");
766  }
767 
768  nLayers_ = dict.lookup<label>("nPoints");
769 
770  interpolate_ = dict.lookupOrDefault<bool>("interpolate", false);
771 
772  fields_ = dict.lookup<wordList>("fields");
773 
774  axis_ =
776  [
777  dict.lookupOrDefault<word>
778  (
779  "axis",
781  )
782  ];
783 
784  formatter_ = setWriter::New(dict.lookup("setFormat"), dict);
785 
786  nOptimiseIter_ = dict.lookupOrDefault("nOptimiseIter", 2);
787 
788  clear();
789 
790  return true;
791 }
792 
793 
795 {
796  wordList result(fields_);
797 
798  if (distanceName_ != word::null)
799  {
800  result.append(distanceName_);
801  }
802 
803  return result;
804 }
805 
806 
808 {
809  return true;
810 }
811 
812 
814 {
815  if (!weights_.valid())
816  {
817  calcWeights();
818  }
819 
820  const polyPatch& pp = mesh_.boundaryMesh()[patchName_];
821 
822  const bool writeThickness =
823  !interpolate_
824  && (
826  || axis_ == coordSet::axisType::DISTANCE
827  );
828 
829  // Create list of field names
831  if (writeThickness)
832  {
833  fieldNames.append("thickness");
834  }
835  forAll(fields_, fieldi)
836  {
837  if
838  (
839  false
840  #define FoundTypeField(Type, nullArg) \
841  || foundObject<VolField<Type>>(fields_[fieldi])
843  #undef FoundTypeField
844  )
845  {
846  fieldNames.append(fields_[fieldi]);
847  }
848  else
849  {
850  cannotFindObject(fields_[fieldi]);
851  }
852  }
853 
854  // Calculate the field values
855  #define DeclareTypeFieldValues(Type, nullArg) \
856  PtrList<Field<Type>> Type##FieldValues(fieldNames.size());
858  #undef DeclareTypeFieldValues
859  if (writeThickness)
860  {
861  scalarFieldValues.set(0, new scalarField(layerThicknesses_));
862  }
863  for (label fieldi = writeThickness; fieldi < fieldNames.size(); ++ fieldi)
864  {
865  #define CollapseTypeFields(Type, nullArg) \
866  if (mesh_.foundObject<VolField<Type>>(fieldNames[fieldi])) \
867  { \
868  const VolField<Type>& field = \
869  mesh_.lookupObject<VolField<Type>>(fieldNames[fieldi]); \
870  \
871  Type##FieldValues.set \
872  ( \
873  fieldi, \
874  applyWeights(weights_, field.boundaryField()[pp.index()]) \
875  ); \
876  }
878  #undef CollapseTypeFields
879  }
880 
881  // Write
882  if (Pstream::master() && layerPositions_->size())
883  {
884  mkDir(outputPath());
885 
886  formatter_->write
887  (
888  outputPath(),
889  typeName,
890  coordSet
891  (
892  identityMap(layerPositions_->size()),
893  word::null,
894  layerPositions_,
896  layerDistances_,
898  ),
899  fieldNames
900  #define TypeFieldValuesParameter(Type, nullArg) , Type##FieldValues
903  );
904  }
905 
906  return true;
907 }
908 
909 
911 (
912  const polyMesh& mesh
913 )
914 {
915  if (&mesh == &mesh_)
916  {
917  clear();
918  }
919 }
920 
921 
923 (
924  const polyTopoChangeMap& map
925 )
926 {
927  if (&map.mesh() == &mesh_)
928  {
929  clear();
930  }
931 }
932 
933 
935 (
936  const polyMeshMap& map
937 )
938 {
939  if (&map.mesh() == &mesh_)
940  {
941  clear();
942  }
943 }
944 
945 
947 (
948  const polyDistributionMap& map
949 )
950 {
951  if (&map.mesh() == &mesh_)
952  {
953  clear();
954  }
955 }
956 
957 
958 // ************************************************************************* //
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:433
Macros for easy insertion into run-time selection tables.
static volPointInterpolation & New(const word &name, const fvMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Pre-declare SubField and related Field type.
Definition: Field.H:83
SubField< vector > subField
Declare type of subField.
Definition: Field.H:101
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
virtual Ostream & write(const char)=0
Write character.
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.
fileName globalPath() const
Return the global path.
Definition: TimePaths.H:132
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 direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:105
static const Form nan
Definition: VectorSpace.H:124
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const word & name() const
Return const reference to name.
A class for handling file names.
Definition: fileName.H:82
Abstract base-class for Time/database functionObjects.
const Time & time_
Reference to time.
const word & name() const
Return the name of this functionObject.
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
const word & name() const
Return reference to name.
Definition: fvMesh.H:434
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:270
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
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:197
tmp< PointField< Type > > interpolate(const VolField< Type > &) const
Interpolate volField using inverse distance weighting.
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
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
static List< word > fieldNames
Definition: globalFoam.H:46
const pointField & points
tUEqn clear()
Tuple2< vector, std::tuple< AreaIntegralType< Types > ... > > faceCutAreaIntegral(const face &f, const vector &fArea, const std::tuple< Types ... > &fPsis, const List< labelPair > &fCuts, const pointField &ps, const std::tuple< const Field< Types > &... > &pPsis, const scalarField &pAlphas, const scalar isoAlpha, const bool below)
Compute the face-cut-area and face-cut-area-integral of the given properties.
vector faceCutArea(const face &f, const vector &fArea, const List< labelPair > &fCuts, const pointField &ps, const scalarField &pAlphas, const scalar isoAlpha, const bool below)
Compute the face-cut-area.
Tuple2< vector, std::tuple< Types ... > > faceAreaAverage(const FaceValues< point > &fPs, const point &fPAvg, const std::tuple< FaceValues< Types > ... > &fPsis, const std::tuple< Types ... > &fPsiAvgs)
Compute the face-area and face-area-averages of the given properties.
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)
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
List< word > wordList
A List of words.
Definition: fileName.H:54
dimensioned< tensor > dimensionedTensor
Dimensioned tensor 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 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:258
const dimensionSet dimless
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
VolField< tensor > volTensorField
Definition: volFieldsFwd.H:66
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
PointField< scalar > pointScalarField
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
dimensionSet normalised(const dimensionSet &)
Definition: dimensionSet.C:507
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
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
List< face > faceList
Definition: faceListFwd.H:41
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:29
const scalar xMax
Definition: createFields.H:30