surfaceFieldValue.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 "surfaceFieldValue.H"
27 #include "processorFvPatch.H"
28 #include "processorCyclicFvPatch.H"
29 #include "sampledSurface.H"
30 #include "mergePoints.H"
31 #include "indirectPrimitivePatch.H"
32 #include "PatchTools.H"
33 #include "fvMeshStitcher.H"
34 #include "polyTopoChangeMap.H"
35 #include "polyMeshMap.H"
36 #include "polyDistributionMap.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 namespace functionObjects
44 {
45 namespace fieldValues
46 {
50 }
51 }
52 }
53 
54 template<>
55 const char* Foam::NamedEnum
56 <
58  3
59 >::names[] =
60 {
61  "faceZone",
62  "patch",
63  "sampledSurface"
64 };
65 
66 template<>
67 const char* Foam::NamedEnum
68 <
70  14
71 >::names[] =
72 {
73  "none",
74  "sum",
75  "sumMag",
76  "orientedSum",
77  "average",
78  "areaAverage",
79  "areaIntegrate",
80  "min",
81  "max",
82  "minMag",
83  "maxMag",
84  "CoV",
85  "areaNormalAverage",
86  "areaNormalIntegrate"
87 };
88 
89 const Foam::NamedEnum
90 <
92  3
94 
95 const Foam::NamedEnum
96 <
98  14
100 
101 
102 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
103 
104 void Foam::functionObjects::fieldValues::surfaceFieldValue::setFaceZoneFaces()
105 {
107 
108  if (zoneId < 0)
109  {
111  << type() << " " << name() << ": "
113  << "(" << selectionName_ << "):" << nl
114  << " Unknown face zone name: " << selectionName_
115  << ". Valid face zones are: " << mesh_.faceZones().names()
116  << nl << exit(FatalError);
117  }
118 
119  // Ensure addressing is built on all processes
122 
123  const faceZone& zone = mesh_.faceZones()[zoneId];
124 
125  DynamicList<label> faceIds(zone.size());
126  DynamicList<label> facePatchIds(zone.size());
127  DynamicList<label> faceSigns(zone.size());
128 
129  forAll(zone, zoneFacei)
130  {
131  const label facei = zone[zoneFacei];
132  const label faceSign = zone.flipMap()[zoneFacei] ? -1 : 1;
133 
134  if (mesh_.isInternalFace(facei))
135  {
136  faceIds.append(facei);
137  facePatchIds.append(-1);
138  faceSigns.append(faceSign);
139  }
140  else
141  {
142  const label bFacei = facei - mesh_.nInternalFaces();
143 
144  const labelUList patches = mesh_.polyBFacePatches()[bFacei];
145  const labelUList patchFaces = mesh_.polyBFacePatchFaces()[bFacei];
146 
147  forAll(patches, i)
148  {
149  // Don't include processor patch faces twice
150  const fvPatch& fvp = mesh_.boundary()[patches[i]];
151  if
152  (
153  isType<processorFvPatch>(fvp)
154  && refCast<const processorFvPatch>(fvp).neighbour()
155  )
156  {
157  continue;
158  }
159 
160  faceIds.append(patchFaces[i]);
161  facePatchIds.append(patches[i]);
162  faceSigns.append(faceSign);
163  }
164  }
165  }
166 
167  faceId_.transfer(faceIds);
168  facePatchId_.transfer(facePatchIds);
169  faceSign_.transfer(faceSigns);
170 
171  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
172 }
173 
174 
175 void Foam::functionObjects::fieldValues::surfaceFieldValue::setPatchFaces()
176 {
177  const label patchId = mesh_.boundaryMesh().findPatchID(selectionName_);
178 
179  if (patchId < 0)
180  {
182  << type() << " " << name() << ": "
183  << selectionTypeNames[selectionType_]
184  << "(" << selectionName_ << "):" << nl
185  << " Unknown patch name: " << selectionName_
186  << ". Valid patch names are: "
187  << mesh_.boundaryMesh().names() << nl
188  << exit(FatalError);
189  }
190 
191  const fvPatch& fvp = mesh_.boundary()[patchId];
192 
193  faceId_ = identityMap(fvp.size());
194  facePatchId_ = labelList(fvp.size(), patchId);
195  faceSign_ = labelList(fvp.size(), 1);
196 
197  // If we have selected a cyclic, also include any associated processor
198  // cyclic faces
199  forAll(mesh_.boundary(), patchi)
200  {
201  const fvPatch& fvp = mesh_.boundary()[patchi];
202 
203  if
204  (
205  isA<processorCyclicFvPatch>(fvp)
206  && refCast<const processorCyclicFvPatch>(fvp).referPatchID() == patchId
207  )
208  {
209  faceId_.append(identityMap(fvp.size()));
210  facePatchId_.append(labelList(fvp.size(), patchi));
211  faceSign_.append(labelList(fvp.size(), 1));
212  }
213  }
214 
215  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
216 }
217 
218 
219 void Foam::functionObjects::fieldValues::surfaceFieldValue::sampledSurfaceFaces
220 (
221  const dictionary& dict
222 )
223 {
224  surfacePtr_ = sampledSurface::New
225  (
226  name(),
227  mesh_,
228  dict.subDict("sampledSurfaceDict")
229  );
230  surfacePtr_().update();
231  nFaces_ = returnReduce(surfacePtr_().faces().size(), sumOp<label>());
232 }
233 
234 
235 void Foam::functionObjects::fieldValues::surfaceFieldValue::combineMeshGeometry
236 (
237  faceList& faces,
239 ) const
240 {
241  List<faceList> allFaces(Pstream::nProcs());
242  List<pointField> allPoints(Pstream::nProcs());
243 
244  labelList globalFacesIs(faceId_);
245  forAll(globalFacesIs, i)
246  {
247  if (facePatchId_[i] != -1)
248  {
249  label patchi = facePatchId_[i];
250  globalFacesIs[i] += mesh_.boundaryMesh()[patchi].start();
251  }
252  }
253 
254  // Add local faces and points to the all* lists
256  (
257  IndirectList<face>(mesh_.faces(), globalFacesIs),
258  mesh_.points()
259  );
260  allFaces[Pstream::myProcNo()] = pp.localFaces();
261  allPoints[Pstream::myProcNo()] = pp.localPoints();
262 
263  Pstream::gatherList(allFaces);
264  Pstream::gatherList(allPoints);
265 
266  // Renumber and flatten
267  label nFaces = 0;
268  label nPoints = 0;
269  forAll(allFaces, proci)
270  {
271  nFaces += allFaces[proci].size();
272  nPoints += allPoints[proci].size();
273  }
274 
275  faces.setSize(nFaces);
276  points.setSize(nPoints);
277 
278  nFaces = 0;
279  nPoints = 0;
280 
281  // My own data first
282  {
283  const faceList& fcs = allFaces[Pstream::myProcNo()];
284  forAll(fcs, i)
285  {
286  const face& f = fcs[i];
287  face& newF = faces[nFaces++];
288  newF.setSize(f.size());
289  forAll(f, fp)
290  {
291  newF[fp] = f[fp] + nPoints;
292  }
293  }
294 
295  const pointField& pts = allPoints[Pstream::myProcNo()];
296  forAll(pts, i)
297  {
298  points[nPoints++] = pts[i];
299  }
300  }
301 
302  // Other proc data follows
303  forAll(allFaces, proci)
304  {
305  if (proci != Pstream::myProcNo())
306  {
307  const faceList& fcs = allFaces[proci];
308  forAll(fcs, i)
309  {
310  const face& f = fcs[i];
311  face& newF = faces[nFaces++];
312  newF.setSize(f.size());
313  forAll(f, fp)
314  {
315  newF[fp] = f[fp] + nPoints;
316  }
317  }
318 
319  const pointField& pts = allPoints[proci];
320  forAll(pts, i)
321  {
322  points[nPoints++] = pts[i];
323  }
324  }
325  }
326 
327  // Merge
328  labelList oldToNew;
329  pointField newPoints;
330  bool hasMerged = mergePoints
331  (
332  points,
333  small,
334  false,
335  oldToNew,
336  newPoints
337  );
338 
339  if (hasMerged)
340  {
341  if (debug)
342  {
343  Pout<< "Merged from " << points.size()
344  << " down to " << newPoints.size() << " points" << endl;
345  }
346 
347  points.transfer(newPoints);
348  forAll(faces, i)
349  {
350  inplaceRenumber(oldToNew, faces[i]);
351  }
352  }
353 }
354 
355 
356 void Foam::functionObjects::fieldValues::surfaceFieldValue::
357 combineSurfaceGeometry
358 (
359  faceList& faces,
361 ) const
362 {
363  if (selectionType_ == selectionTypes::sampledSurface)
364  {
365  const sampledSurface& s = surfacePtr_();
366 
367  if (Pstream::parRun())
368  {
369  // Dimension as fraction of mesh bounding box
370  scalar mergeDim = 1e-10*mesh_.bounds().mag();
371 
372  labelList pointsMap;
373 
375  (
376  mergeDim,
378  (
379  SubList<face>(s.faces(), s.faces().size()),
380  s.points()
381  ),
382  points,
383  faces,
384  pointsMap
385  );
386  }
387  else
388  {
389  faces = s.faces();
390  points = s.points();
391  }
392  }
393 }
394 
395 
396 Foam::scalar
397 Foam::functionObjects::fieldValues::surfaceFieldValue::totalArea() const
398 {
399  if (selectionType_ == selectionTypes::sampledSurface)
400  {
401  return gSum(surfacePtr_().magSf());
402  }
403  else
404  {
405  return gSum(filterField(mesh_.magSf()));
406  }
407 }
408 
409 
410 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
411 
413 (
414  const dictionary& dict
415 )
416 {
417  switch (selectionType_)
418  {
419  case selectionTypes::faceZone:
420  {
421  dict.lookupBackwardsCompatible({"faceZone", "name"})
422  >> selectionName_;
423  setFaceZoneFaces();
424  break;
425  }
426  case selectionTypes::patch:
427  {
428  dict.lookupBackwardsCompatible({"patch", "name"}) >> selectionName_;
429  setPatchFaces();
430  break;
431  }
432  case selectionTypes::sampledSurface:
433  {
434  sampledSurfaceFaces(dict);
435  selectionName_ = surfacePtr_().name();
436  break;
437  }
438  default:
439  {
441  << type() << " " << name() << ": "
442  << selectionTypeNames[selectionType_]
443  << "(" << selectionName_ << "):" << nl
444  << " Unknown selection type. Valid selection types are:"
445  << selectionTypeNames.sortedToc() << nl << exit(FatalError);
446  }
447  }
448 
449  if (nFaces_ == 0 && (!mesh_.stitcher().stitches() || !mesh_.conformal()))
450  {
452  << type() << " " << name() << ": "
453  << selectionTypeNames[selectionType_]
454  << "(" << selectionName_ << "):" << nl
455  << " selection has no faces" << exit(FatalError);
456  }
457 
458  if (selectionType_ == selectionTypes::sampledSurface)
459  {
460  surfacePtr_().update();
461  }
462 
463  totalArea_ = totalArea();
464 
465  Info<< type() << " " << name() << ":" << nl
466  << " total faces = " << nFaces_
467  << nl
468  << " total area = " << totalArea_
469  << nl;
470 
471  if (dict.readIfPresent("weightFields", weightFieldNames_))
472  {
473  Info<< name() << " " << operationTypeNames_[operation_]
474  << " weight fields " << weightFieldNames_;
475  }
476  else if (dict.found("weightField"))
477  {
478  weightFieldNames_.setSize(1);
479  dict.lookup("weightField") >> weightFieldNames_[0];
480 
481  Info<< name() << " " << operationTypeNames_[operation_]
482  << " weight field " << weightFieldNames_[0];
483  }
484 
485  if (dict.readIfPresent("scaleFactor", scaleFactor_))
486  {
487  Info<< " scale factor = " << scaleFactor_ << nl;
488  }
489 
490  Info<< nl << endl;
491 
492  if (writeFields_)
493  {
494  const word surfaceFormat(dict.lookup("surfaceFormat"));
495 
496  surfaceWriterPtr_.reset
497  (
498  surfaceWriter::New(surfaceFormat, dict).ptr()
499  );
500  }
501 }
502 
503 
505 (
506  const label i
507 )
508 {
509  if (operation_ != operationType::none)
510  {
511  writeCommented(file(), "Selection type : ");
512  file()
513  << selectionTypeNames[selectionType_] << " "
514  << selectionName_ << endl;
515  writeCommented(file(), "Faces : ");
516  file() << nFaces_ << endl;
517  writeCommented(file(), "Area : ");
518  file() << totalArea_ << endl;
519 
520  writeCommented(file(), "Time");
521  if (writeArea_)
522  {
523  file() << tab << "Area";
524  }
525 
526  forAll(fields_, fieldi)
527  {
528  file()
529  << tab << operationTypeNames_[operation_]
530  << "(" << fields_[fieldi] << ")";
531  }
532 
533  file() << endl;
534  }
535 }
536 
537 
539 (
540  const Field<scalar>& values,
541  const scalarField& signs,
542  const scalarField& weights,
543  const vectorField& Sf,
544  scalar& result
545 ) const
546 {
547  switch (operation_)
548  {
549  default:
550  {
551  // Fall through to same-type operations
552  return processValuesTypeType
553  (
554  values,
555  signs,
556  weights,
557  Sf,
558  result
559  );
560  }
561  }
562 }
563 
564 
566 (
567  const Field<vector>& values,
568  const scalarField& signs,
569  const scalarField& weights,
570  const vectorField& Sf,
571  scalar& result
572 ) const
573 {
574  switch (operation_)
575  {
576  case operationType::areaNormalAverage:
577  {
578  result = gSum(weights*values & Sf)/gSum(mag(weights*Sf));
579  return true;
580  }
581  case operationType::areaNormalIntegrate:
582  {
583  result = gSum(weights*values & Sf);
584  return true;
585  }
586  default:
587  {
588  // No fall through. Different types.
589  return false;
590  }
591  }
592 }
593 
594 
595 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
596 
598 (
599  const word& name,
600  const Time& runTime,
601  const dictionary& dict
602 )
603 :
604  fieldValue(name, runTime, dict, typeName),
605  dict_(dict),
606  surfaceWriterPtr_(nullptr),
607  selectionType_
608  (
609  selectionTypeNames.read
610  (
611  dict.lookupBackwardsCompatible({"select", "regionType"})
612  )
613  ),
614  selectionName_(word::null),
615  operation_(operationTypeNames_.read(dict.lookup("operation"))),
616  weightFieldNames_(),
617  scaleFactor_(1),
618  writeArea_(dict.lookupOrDefault("writeArea", false)),
619  nFaces_(0),
620  faceId_(),
621  facePatchId_(),
622  faceSign_()
623 {
624  read(dict_);
625 }
626 
628 (
629  const word& name,
630  const objectRegistry& obr,
631  const dictionary& dict
632 )
633 :
634  fieldValue(name, obr, dict, typeName),
635  dict_(dict),
636  surfaceWriterPtr_(nullptr),
637  selectionType_
638  (
639  selectionTypeNames.read
640  (
641  dict.lookupBackwardsCompatible({"select", "regionType"})
642  )
643  ),
644  selectionName_(word::null),
645  operation_(operationTypeNames_.read(dict.lookup("operation"))),
646  weightFieldNames_(),
647  scaleFactor_(1),
648  writeArea_(dict.lookupOrDefault("writeArea", false)),
649  nFaces_(0),
650  faceId_(),
651  facePatchId_(),
652  faceSign_()
653 {
654  read(dict_);
655 }
656 
657 
658 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
659 
661 {}
662 
663 
664 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
665 
667 (
668  const dictionary& dict
669 )
670 {
672  initialise(dict);
673 
674  return true;
675 }
676 
677 
679 {
680  if (operation_ != operationType::none)
681  {
683  }
684 
685  if (selectionType_ == selectionTypes::sampledSurface)
686  {
687  surfacePtr_().update();
688  }
689 
690  if (operation_ != operationType::none && Pstream::master())
691  {
692  writeTime(file());
693  }
694 
695  if (writeArea_)
696  {
697  totalArea_ = totalArea();
698  if (operation_ != operationType::none && Pstream::master())
699  {
700  file() << tab << totalArea_;
701  }
702  Log << " total area = " << totalArea_ << endl;
703  }
704 
705  // Write the surface geometry
706  if (writeFields_)
707  {
708  faceList faces;
710 
711  if (selectionType_ == selectionTypes::sampledSurface)
712  {
713  combineSurfaceGeometry(faces, points);
714  }
715  else
716  {
717  combineMeshGeometry(faces, points);
718  }
719 
720  if (Pstream::master())
721  {
722  surfaceWriterPtr_->write
723  (
724  outputDir(),
725  selectionTypeNames[selectionType_] + ("_" + selectionName_),
726  points,
727  faces
728  );
729  }
730  }
731 
732  // Construct the sign and weight fields and the surface normals
733  const scalarField signs
734  (
735  selectionType_ == selectionTypes::sampledSurface
736  ? scalarField(surfacePtr_().Sf().size(), 1)
737  : List<scalar>(faceSign_)
738  );
739  scalarField weights(signs.size(), 1);
740  forAll(weightFieldNames_, i)
741  {
742  weights *= getFieldValues<scalar>(weightFieldNames_[i]);
743  }
744  const vectorField Sf
745  (
746  selectionType_ == selectionTypes::sampledSurface
747  ? surfacePtr_().Sf()
748  : (signs*filterField(mesh_.Sf()))()
749  );
750 
751  // Process the fields
752  forAll(fields_, i)
753  {
754  const word& fieldName = fields_[i];
755  bool ok = false;
756 
757  #define writeValuesFieldType(fieldType, none) \
758  ok = \
759  ok \
760  || writeValues<fieldType> \
761  ( \
762  fieldName, \
763  signs, \
764  weights, \
765  Sf \
766  );
768  #undef writeValuesFieldType
769 
770  if (!ok)
771  {
773  << "Requested field " << fieldName
774  << " not found in database and not processed"
775  << endl;
776  }
777  }
778 
779  if (operation_ != operationType::none && Pstream::master())
780  {
781  file() << endl;
782  }
783 
784  Log << endl;
785 
786  return true;
787 }
788 
789 
791 (
792  const polyMesh& mesh
793 )
794 {
795  if (&mesh == &mesh_)
796  {
797  // It may be necessary to reset if the mesh moves. The total area might
798  // change, as might non-conformal faces.
799  initialise(dict_);
800  }
801 }
802 
803 
805 (
806  const polyTopoChangeMap& map
807 )
808 {
809  if (&map.mesh() == &mesh_)
810  {
811  initialise(dict_);
812  }
813 }
814 
815 
817 (
818  const polyMeshMap& map
819 )
820 {
821  if (&map.mesh() == &mesh_)
822  {
823  initialise(dict_);
824  }
825 }
826 
827 
829 (
830  const polyDistributionMap& map
831 )
832 {
833  if (&map.mesh() == &mesh_)
834  {
835  initialise(dict_);
836  }
837 }
838 
839 
840 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: MeshZones.C:341
wordList names() const
Return a list of zone names.
Definition: MeshZones.C:256
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:54
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< FaceList, PointField > &p, Field< typename PrimitivePatch< FaceList, PointField >::PointType > &mergedPoints, List< typename PrimitivePatch< FaceList, PointField >::FaceType > &mergedFaces, labelList &pointMergeMap)
Gather points and faces onto master and merge into single patch.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
void append(T *)
Append an element at the end of the list.
Definition: PtrListI.H:39
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 label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Abstract base-class for Time/database functionObjects.
virtual const word & type() const =0
Runtime type information.
const word & name() const
Return the name of this functionObject.
Base class for field value -based function objects.
Definition: fieldValue.H:66
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: fieldValue.C:88
virtual bool write()
Write.
Definition: fieldValue.C:112
surfaceFieldValue(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
const labelList & faceSign() const
Return the list of +1/-1 representing face flip map.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
void initialise(const dictionary &dict)
Initialise, e.g. face addressing.
static const NamedEnum< operationType, 14 > operationTypeNames_
Operation type names.
word selectionName_
Name of face selection (patch, faceZone, etc.)
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual void writeFileHeader(const label i)
Output file header information.
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
static const NamedEnum< selectionTypes, 3 > selectionTypeNames
Selection type names.
bool processValues(const Field< Type > &values, const scalarField &signs, const scalarField &weights, const vectorField &Sf, ResultType &result) const
Apply the operation to the values, and return true if successful.
labelList facePatchId_
Local list of patch ID per face.
labelList faceSign_
List of +1/-1 representing face flip map.
virtual bool read(const dictionary &)
Read from dictionary.
const fvMesh & mesh_
Reference to the fvMesh.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:895
const UCompactListList< label > & polyBFacePatches() const
Return poly-bFace-patch addressing.
Definition: fvMesh.C:953
const UCompactListList< label > & polyBFacePatchFaces() const
Return poly-bFace-patch-face addressing.
Definition: fvMesh.C:1029
Registry of regIOobjects.
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:447
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const polyMesh & mesh() const
Return polyMesh.
label nInternalFaces() const
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
static autoPtr< sampledSurface > New(const word &name, const polyMesh &, const dictionary &)
Return a reference to the selected surface.
static autoPtr< surfaceWriter > New(const word &writeType, const IOstream::streamFormat writeFormat, const IOstream::compressionType writeCompression)
Select given write options.
Definition: surfaceWriter.C:75
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
label patchi
const pointField & points
label nPoints
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
label patchId(-1)
const fvPatchList & patches
Merge points. See below.
#define WarningInFunction
Report a warning using Foam::Warning.
#define Log
Report write to Foam::Info if the local log switch is true.
defineTypeNameAndDebug(fieldValueDelta, 0)
addToRunTimeSelectionTable(functionObject, fieldValueDelta, dictionary)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Type gSum(const FieldField< Field, Type > &f)
const doubleScalar e
Definition: doubleScalar.H:105
bool read(const char *, int32_t &)
Definition: int32IO.C:85
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
static const char tab
Definition: Ostream.H:259
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
FOR_ALL_FIELD_TYPES(DefineContiguousFvWallLocationDataType)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
Addressing for a faceList slice.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
label mergePoints(const UList< Type > &points, const scalar mergeTol, const bool verbose, labelList &pointMap, const Type &origin=Type::zero)
Sorts and merges points. All points closer than/equal mergeTol get merged.
dimensioned< scalar > mag(const dimensioned< Type > &)
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
List< face > faceList
Definition: faceListFwd.H:43
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
UList< label > labelUList
Definition: UList.H:65
static const char nl
Definition: Ostream.H:260
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
labelList f(nPoints)
dictionary dict
#define writeValuesFieldType(fieldType, none)