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-2024 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().toc()
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().findIndex(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).referPatchIndex()
207  == patchId
208  )
209  {
210  faceId_.append(identityMap(fvp.size()));
211  facePatchId_.append(labelList(fvp.size(), patchi));
212  faceSign_.append(labelList(fvp.size(), 1));
213  }
214  }
215 
216  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
217 }
218 
219 
220 void Foam::functionObjects::fieldValues::surfaceFieldValue::sampledSurfaceFaces
221 (
222  const dictionary& dict
223 )
224 {
225  surfacePtr_ = sampledSurface::New
226  (
227  name(),
228  mesh_,
229  dict.subDict("sampledSurfaceDict")
230  );
231  surfacePtr_().update();
232  nFaces_ = returnReduce(surfacePtr_().faces().size(), sumOp<label>());
233 }
234 
235 
236 void Foam::functionObjects::fieldValues::surfaceFieldValue::combineMeshGeometry
237 (
238  faceList& faces,
240 ) const
241 {
242  List<faceList> allFaces(Pstream::nProcs());
243  List<pointField> allPoints(Pstream::nProcs());
244 
245  labelList globalFacesIs(faceId_);
246  forAll(globalFacesIs, i)
247  {
248  if (facePatchId_[i] != -1)
249  {
250  label patchi = facePatchId_[i];
251  globalFacesIs[i] += mesh_.boundaryMesh()[patchi].start();
252  }
253  }
254 
255  // Add local faces and points to the all* lists
257  (
258  IndirectList<face>(mesh_.faces(), globalFacesIs),
259  mesh_.points()
260  );
261  allFaces[Pstream::myProcNo()] = pp.localFaces();
262  allPoints[Pstream::myProcNo()] = pp.localPoints();
263 
264  Pstream::gatherList(allFaces);
265  Pstream::gatherList(allPoints);
266 
267  // Renumber and flatten
268  label nFaces = 0;
269  label nPoints = 0;
270  forAll(allFaces, proci)
271  {
272  nFaces += allFaces[proci].size();
273  nPoints += allPoints[proci].size();
274  }
275 
276  faces.setSize(nFaces);
277  points.setSize(nPoints);
278 
279  nFaces = 0;
280  nPoints = 0;
281 
282  // My own data first
283  {
284  const faceList& fcs = allFaces[Pstream::myProcNo()];
285  forAll(fcs, i)
286  {
287  const face& f = fcs[i];
288  face& newF = faces[nFaces++];
289  newF.setSize(f.size());
290  forAll(f, fp)
291  {
292  newF[fp] = f[fp] + nPoints;
293  }
294  }
295 
296  const pointField& pts = allPoints[Pstream::myProcNo()];
297  forAll(pts, i)
298  {
299  points[nPoints++] = pts[i];
300  }
301  }
302 
303  // Other proc data follows
304  forAll(allFaces, proci)
305  {
306  if (proci != Pstream::myProcNo())
307  {
308  const faceList& fcs = allFaces[proci];
309  forAll(fcs, i)
310  {
311  const face& f = fcs[i];
312  face& newF = faces[nFaces++];
313  newF.setSize(f.size());
314  forAll(f, fp)
315  {
316  newF[fp] = f[fp] + nPoints;
317  }
318  }
319 
320  const pointField& pts = allPoints[proci];
321  forAll(pts, i)
322  {
323  points[nPoints++] = pts[i];
324  }
325  }
326  }
327 
328  // Merge
329  labelList oldToNew;
330  pointField newPoints;
331  bool hasMerged = mergePoints
332  (
333  points,
334  small,
335  false,
336  oldToNew,
337  newPoints
338  );
339 
340  if (hasMerged)
341  {
342  if (debug)
343  {
344  Pout<< "Merged from " << points.size()
345  << " down to " << newPoints.size() << " points" << endl;
346  }
347 
348  points.transfer(newPoints);
349  forAll(faces, i)
350  {
351  inplaceRenumber(oldToNew, faces[i]);
352  }
353  }
354 }
355 
356 
357 void Foam::functionObjects::fieldValues::surfaceFieldValue::
358 combineSurfaceGeometry
359 (
360  faceList& faces,
362 ) const
363 {
364  if (selectionType_ == selectionTypes::sampledSurface)
365  {
366  const sampledSurface& s = surfacePtr_();
367 
368  if (Pstream::parRun())
369  {
370  // Dimension as fraction of mesh bounding box
371  scalar mergeDim = 1e-10*mesh_.bounds().mag();
372 
373  labelList pointsMap;
374 
376  (
377  mergeDim,
379  (
380  SubList<face>(s.faces(), s.faces().size()),
381  s.points()
382  ),
383  points,
384  faces,
385  pointsMap
386  );
387  }
388  else
389  {
390  faces = s.faces();
391  points = s.points();
392  }
393  }
394 }
395 
396 
397 Foam::scalar
398 Foam::functionObjects::fieldValues::surfaceFieldValue::totalArea() const
399 {
400  if (selectionType_ == selectionTypes::sampledSurface)
401  {
402  return gSum(surfacePtr_().magSf());
403  }
404  else
405  {
406  return gSum(filterField(mesh_.magSf()));
407  }
408 }
409 
410 
411 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
412 
414 (
415  const dictionary& dict
416 )
417 {
418  switch (selectionType_)
419  {
420  case selectionTypes::faceZone:
421  {
422  dict.lookupBackwardsCompatible({"faceZone", "name"})
423  >> selectionName_;
424  setFaceZoneFaces();
425  break;
426  }
427  case selectionTypes::patch:
428  {
429  dict.lookupBackwardsCompatible({"patch", "name"}) >> selectionName_;
430  setPatchFaces();
431  break;
432  }
433  case selectionTypes::sampledSurface:
434  {
435  sampledSurfaceFaces(dict);
436  selectionName_ = surfacePtr_().name();
437  break;
438  }
439  default:
440  {
442  << type() << " " << name() << ": "
443  << selectionTypeNames[selectionType_]
444  << "(" << selectionName_ << "):" << nl
445  << " Unknown selection type. Valid selection types are:"
446  << selectionTypeNames.sortedToc() << nl << exit(FatalError);
447  }
448  }
449 
450  if (nFaces_ == 0 && (!mesh_.stitcher().stitches() || !mesh_.conformal()))
451  {
453  << type() << " " << name() << ": "
454  << selectionTypeNames[selectionType_]
455  << "(" << selectionName_ << "):" << nl
456  << " selection has no faces" << exit(FatalError);
457  }
458 
459  if (selectionType_ == selectionTypes::sampledSurface)
460  {
461  surfacePtr_().update();
462  }
463 
464  totalArea_ = totalArea();
465 
466  Info<< type() << " " << name() << ":" << nl
467  << " total faces = " << nFaces_
468  << nl
469  << " total area = " << totalArea_
470  << nl;
471 
472  if (dict.readIfPresent("weightFields", weightFieldNames_))
473  {
474  Info<< name() << " " << operationTypeNames_[operation_]
475  << " weight fields " << weightFieldNames_;
476  }
477  else if (dict.found("weightField"))
478  {
479  weightFieldNames_.setSize(1);
480  dict.lookup("weightField") >> weightFieldNames_[0];
481 
482  Info<< name() << " " << operationTypeNames_[operation_]
483  << " weight field " << weightFieldNames_[0];
484  }
485 
486  if (dict.readIfPresent("scaleFactor", scaleFactor_))
487  {
488  Info<< " scale factor = " << scaleFactor_ << nl;
489  }
490 
491  Info<< nl << endl;
492 
493  if (writeFields_)
494  {
495  const word surfaceFormat(dict.lookup("surfaceFormat"));
496 
497  surfaceWriterPtr_.reset
498  (
499  surfaceWriter::New(surfaceFormat, dict).ptr()
500  );
501  }
502 }
503 
504 
506 (
507  const label i
508 )
509 {
510  if (operation_ != operationType::none)
511  {
512  writeCommented(file(), "Selection type : ");
513  file()
514  << selectionTypeNames[selectionType_] << " "
515  << selectionName_ << endl;
516  writeCommented(file(), "Faces : ");
517  file() << nFaces_ << endl;
518  writeCommented(file(), "Area : ");
519  file() << totalArea_ << endl;
520 
521  writeCommented(file(), "Time");
522  if (writeArea_)
523  {
524  file() << tab << "Area";
525  }
526 
527  forAll(fields_, fieldi)
528  {
529  file()
530  << tab << operationTypeNames_[operation_]
531  << "(" << fields_[fieldi] << ")";
532  }
533 
534  file() << endl;
535  }
536 }
537 
538 
540 (
541  const Field<scalar>& values,
542  const scalarField& signs,
543  const scalarField& weights,
544  const vectorField& Sf,
545  scalar& result
546 ) const
547 {
548  switch (operation_)
549  {
550  default:
551  {
552  // Fall through to same-type operations
553  return processValuesTypeType
554  (
555  values,
556  signs,
557  weights,
558  Sf,
559  result
560  );
561  }
562  }
563 }
564 
565 
567 (
568  const Field<vector>& values,
569  const scalarField& signs,
570  const scalarField& weights,
571  const vectorField& Sf,
572  scalar& result
573 ) const
574 {
575  switch (operation_)
576  {
577  case operationType::areaNormalAverage:
578  {
579  result = gSum(weights*values & Sf)/gSum(mag(weights*Sf));
580  return true;
581  }
582  case operationType::areaNormalIntegrate:
583  {
584  result = gSum(weights*values & Sf);
585  return true;
586  }
587  default:
588  {
589  // No fall through. Different types.
590  return false;
591  }
592  }
593 }
594 
595 
596 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
597 
599 (
600  const word& name,
601  const Time& runTime,
602  const dictionary& dict
603 )
604 :
605  fieldValue(name, runTime, dict, typeName),
606  dict_(dict),
607  surfaceWriterPtr_(nullptr),
608  selectionType_
609  (
610  selectionTypeNames.read
611  (
612  dict.lookupBackwardsCompatible({"select", "regionType"})
613  )
614  ),
615  selectionName_(word::null),
616  operation_(operationTypeNames_.read(dict.lookup("operation"))),
617  weightFieldNames_(),
618  scaleFactor_(1),
619  writeArea_(dict.lookupOrDefault("writeArea", false)),
620  nFaces_(0),
621  faceId_(),
622  facePatchId_(),
623  faceSign_()
624 {
625  read(dict_);
626 }
627 
629 (
630  const word& name,
631  const objectRegistry& obr,
632  const dictionary& dict
633 )
634 :
635  fieldValue(name, obr, dict, typeName),
636  dict_(dict),
637  surfaceWriterPtr_(nullptr),
638  selectionType_
639  (
640  selectionTypeNames.read
641  (
642  dict.lookupBackwardsCompatible({"select", "regionType"})
643  )
644  ),
645  selectionName_(word::null),
646  operation_(operationTypeNames_.read(dict.lookup("operation"))),
647  weightFieldNames_(),
648  scaleFactor_(1),
649  writeArea_(dict.lookupOrDefault("writeArea", false)),
650  nFaces_(0),
651  faceId_(),
652  facePatchId_(),
653  faceSign_()
654 {
655  read(dict_);
656 }
657 
658 
659 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
660 
662 {}
663 
664 
665 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
666 
668 (
669  const dictionary& dict
670 )
671 {
673  initialise(dict);
674 
675  return true;
676 }
677 
678 
680 {
681  // Look to see if any fields exist. Use the flag to suppress output later.
682  bool anyFields = false;
683  forAll(fields_, i)
684  {
685  #define validFieldType(fieldType, none) \
686  anyFields = anyFields || validField<fieldType>(fields_[i]);
688  #undef validFieldType
689  }
690  if (!anyFields && fields_.size() > 1) // (error for 1 will happen below)
691  {
692  cannotFindObjects(fields_);
693  }
694 
695  // Initialise the file, write the header, etc...
696  if (anyFields && operation_ != operationType::none)
697  {
699  }
700 
701  // Update the surface
702  if (selectionType_ == selectionTypes::sampledSurface)
703  {
704  surfacePtr_().update();
705  }
706 
707  // Write the time
708  if (anyFields && operation_ != operationType::none && Pstream::master())
709  {
710  writeTime(file());
711  }
712 
713  // Write the total area if necessary
714  if (writeArea_)
715  {
716  totalArea_ = totalArea();
717  if (anyFields && operation_ != operationType::none && Pstream::master())
718  {
719  file() << tab << totalArea_;
720  }
721  Log << " total area = " << totalArea_ << endl;
722  }
723 
724  // Write the surface geometry
725  if (writeFields_)
726  {
727  faceList faces;
729 
730  if (selectionType_ == selectionTypes::sampledSurface)
731  {
732  combineSurfaceGeometry(faces, points);
733  }
734  else
735  {
736  combineMeshGeometry(faces, points);
737  }
738 
739  if (Pstream::master())
740  {
741  surfaceWriterPtr_->write
742  (
743  outputDir(),
744  selectionTypeNames[selectionType_] + ("_" + selectionName_),
745  points,
746  faces
747  );
748  }
749  }
750 
751  // Construct the sign and weight fields and the surface normals
752  const scalarField signs
753  (
754  selectionType_ == selectionTypes::sampledSurface
755  ? scalarField(surfacePtr_().Sf().size(), 1)
756  : List<scalar>(faceSign_)
757  );
758  scalarField weights(signs.size(), 1);
759  forAll(weightFieldNames_, i)
760  {
761  weights *= getFieldValues<scalar>(weightFieldNames_[i]);
762  }
763  const vectorField Sf
764  (
765  selectionType_ == selectionTypes::sampledSurface
766  ? surfacePtr_().Sf()
767  : (signs*filterField(mesh_.Sf()))()
768  );
769 
770  // Process the fields
771  forAll(fields_, i)
772  {
773  const word& fieldName = fields_[i];
774  bool ok = false;
775 
776  #define writeValuesFieldType(fieldType, none) \
777  ok = ok || writeValues<fieldType>(fieldName, signs, weights, Sf);
779  #undef writeValuesFieldType
780 
781  if (!ok)
782  {
783  cannotFindObject(fieldName);
784  }
785  }
786 
787  // Finalise
788  if (anyFields && operation_ != operationType::none && Pstream::master())
789  {
790  file() << endl;
791  }
792  Log << endl;
793 
794  return true;
795 }
796 
797 
799 (
800  const polyMesh& mesh
801 )
802 {
803  if (&mesh == &mesh_)
804  {
805  // It may be necessary to reset if the mesh moves. The total area might
806  // change, as might non-conformal faces.
807  initialise(dict_);
808  }
809 }
810 
811 
813 (
814  const polyTopoChangeMap& map
815 )
816 {
817  if (&map.mesh() == &mesh_)
818  {
819  initialise(dict_);
820  }
821 }
822 
823 
825 (
826  const polyMeshMap& map
827 )
828 {
829  if (&map.mesh() == &mesh_)
830  {
831  initialise(dict_);
832  }
833 }
834 
835 
837 (
838  const polyDistributionMap& map
839 )
840 {
841  if (&map.mesh() == &mesh_)
842  {
843  initialise(dict_);
844  }
845 }
846 
847 
848 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
wordList toc() const
Return the table of contents.
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
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.
label findIndex(const word &key) const
Return the index of the given the key or -1 if not found.
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:162
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:857
const UCompactListList< label > & polyBFacePatches() const
Return poly-bFace-patch addressing.
Definition: fvMesh.C:911
const UCompactListList< label > & polyBFacePatchFaces() const
Return poly-bFace-patch-face addressing.
Definition: fvMesh.C:987
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 faceZoneList & faceZones() const
Return face zones.
Definition: polyMesh.H:446
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:334
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 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:106
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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
static const char tab
Definition: Ostream.H:265
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
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
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
List< face > faceList
Definition: faceListFwd.H:41
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
UList< label > labelUList
Definition: UList.H:65
static const char nl
Definition: Ostream.H:266
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)
#define validFieldType(fieldType, none)