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-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 "surfaceFieldValue.H"
27 #include "processorFvPatch.H"
28 #include "processorCyclicFvPatch.H"
29 #include "sampledSurface.H"
30 #include "generatedFaceZone.H"
31 #include "mergePoints.H"
32 #include "indirectPrimitivePatch.H"
33 #include "PatchTools.H"
34 #include "fvMeshStitcher.H"
35 #include "polyTopoChangeMap.H"
36 #include "polyMeshMap.H"
37 #include "polyDistributionMap.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 namespace functionObjects
45 {
46 namespace fieldValues
47 {
51 }
52 }
53 }
54 
55 const Foam::NamedEnum
56 <
58  4
60 {
61  "faceZone",
62  "patch",
63  "patches",
64  "sampledSurface"
65 };
66 
67 const Foam::NamedEnum
68 <
70  15
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  "UI",
86  "areaNormalAverage",
87  "areaNormalIntegrate"
88 };
89 
90 
91 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
92 
93 Foam::labelList Foam::functionObjects::fieldValues::surfaceFieldValue::patchis
94 (
95  const wordReList& patchNames
96 ) const
97 {
98  labelList patchis;
99 
100  forAll(patchNames, i)
101  {
102  const labelList patchiis =
104 
105  if (patchiis.empty())
106  {
108  << type() << ' ' << this->name() << ": "
110  << "(" << patchNames[i] << "):" << nl
111  << " Unknown patch name: " << patchNames[i]
112  << ". Valid patch names are: "
113  << mesh_.boundaryMesh().names() << nl
114  << exit(FatalError);
115  }
116 
117  patchis.append(patchiis);
118  }
119 
120  return patchis;
121 }
122 
123 
124 void Foam::functionObjects::fieldValues::surfaceFieldValue::setFaceZoneFaces()
125 {
126  const faceZone& zone = faceZonePtr_->zone();
127 
128  // Ensure addressing is built on all processes
129  mesh_.polyBFacePatches();
130  mesh_.polyBFacePatchFaces();
131 
132  DynamicList<label> faceIds(zone.size());
133  DynamicList<label> facePatchIds(zone.size());
134  DynamicList<label> faceSigns(zone.size());
135 
136  forAll(zone, zoneFacei)
137  {
138  const label facei = zone[zoneFacei];
139  const label faceSign = zone.flipMap()[zoneFacei] ? -1 : 1;
140 
141  if (mesh_.isInternalFace(facei))
142  {
143  faceIds.append(facei);
144  facePatchIds.append(-1);
145  faceSigns.append(faceSign);
146  }
147  else
148  {
149  const label bFacei = facei - mesh_.nInternalFaces();
150 
151  const labelUList patches = mesh_.polyBFacePatches()[bFacei];
152  const labelUList patchFaces = mesh_.polyBFacePatchFaces()[bFacei];
153 
154  forAll(patches, i)
155  {
156  // Don't include processor patch faces twice
157  const fvPatch& fvp = mesh_.boundary()[patches[i]];
158  if
159  (
160  isType<processorFvPatch>(fvp)
161  && refCast<const processorFvPatch>(fvp).neighbour()
162  )
163  {
164  continue;
165  }
166 
167  faceIds.append(patchFaces[i]);
168  facePatchIds.append(patches[i]);
169  faceSigns.append(faceSign);
170  }
171  }
172  }
173 
174  faceId_.transfer(faceIds);
175  facePatchId_.transfer(facePatchIds);
176  faceSign_.transfer(faceSigns);
177 
178  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
179 }
180 
181 
182 void Foam::functionObjects::fieldValues::surfaceFieldValue::setPatchesFaces()
183 {
184  const labelList patchis = this->patchis(patchNames_);
185 
186  faceId_.clear();
187  facePatchId_.clear();
188  faceSign_.clear();
189 
190  forAll(patchis, i)
191  {
192  const label patchi = patchis[i];
193  const fvPatch& fvp = mesh_.boundary()[patchi];
194 
195  faceId_.append(identityMap(fvp.size()));
196  facePatchId_.append(labelList(fvp.size(), patchi));
197  faceSign_.append(labelList(fvp.size(), 1));
198 
199  // If we have selected a cyclic, also include any associated processor
200  // cyclic faces
201  forAll(mesh_.boundary(), pcPatchj)
202  {
203  const fvPatch& pcFvp = mesh_.boundary()[pcPatchj];
204 
205  if
206  (
207  isA<processorCyclicFvPatch>(pcFvp)
208  && refCast<const processorCyclicFvPatch>(pcFvp).referPatchIndex()
209  == patchi
210  )
211  {
212  faceId_.append(identityMap(pcFvp.size()));
213  facePatchId_.append(labelList(pcFvp.size(), pcPatchj));
214  faceSign_.append(labelList(pcFvp.size(), 1));
215  }
216  }
217  }
218 
219  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
220 }
221 
222 
223 void
224 Foam::functionObjects::fieldValues::surfaceFieldValue::setSampledSurfaceFaces()
225 {
226  surfacePtr_().update();
227 
228  nFaces_ = returnReduce(surfacePtr_().faces().size(), sumOp<label>());
229 }
230 
231 
232 void Foam::functionObjects::fieldValues::surfaceFieldValue::combineMeshGeometry
233 (
234  faceList& faces,
236 ) const
237 {
238  List<faceList> allFaces(Pstream::nProcs());
239  List<pointField> allPoints(Pstream::nProcs());
240 
241  labelList globalFacesIs(faceId_);
242  forAll(globalFacesIs, i)
243  {
244  if (facePatchId_[i] != -1)
245  {
246  label patchi = facePatchId_[i];
247  globalFacesIs[i] += mesh_.boundaryMesh()[patchi].start();
248  }
249  }
250 
251  // Add local faces and points to the all* lists
253  (
254  IndirectList<face>(mesh_.faces(), globalFacesIs),
255  mesh_.points()
256  );
257  allFaces[Pstream::myProcNo()] = pp.localFaces();
258  allPoints[Pstream::myProcNo()] = pp.localPoints();
259 
260  Pstream::gatherList(allFaces);
261  Pstream::gatherList(allPoints);
262 
263  // Renumber and flatten
264  label nFaces = 0;
265  label nPoints = 0;
266  forAll(allFaces, proci)
267  {
268  nFaces += allFaces[proci].size();
269  nPoints += allPoints[proci].size();
270  }
271 
272  faces.setSize(nFaces);
273  points.setSize(nPoints);
274 
275  nFaces = 0;
276  nPoints = 0;
277 
278  // My own data first
279  {
280  const faceList& fcs = allFaces[Pstream::myProcNo()];
281  forAll(fcs, i)
282  {
283  const face& f = fcs[i];
284  face& newF = faces[nFaces++];
285  newF.setSize(f.size());
286  forAll(f, fp)
287  {
288  newF[fp] = f[fp] + nPoints;
289  }
290  }
291 
292  const pointField& pts = allPoints[Pstream::myProcNo()];
293  forAll(pts, i)
294  {
295  points[nPoints++] = pts[i];
296  }
297  }
298 
299  // Other proc data follows
300  forAll(allFaces, proci)
301  {
302  if (proci != Pstream::myProcNo())
303  {
304  const faceList& fcs = allFaces[proci];
305  forAll(fcs, i)
306  {
307  const face& f = fcs[i];
308  face& newF = faces[nFaces++];
309  newF.setSize(f.size());
310  forAll(f, fp)
311  {
312  newF[fp] = f[fp] + nPoints;
313  }
314  }
315 
316  const pointField& pts = allPoints[proci];
317  forAll(pts, i)
318  {
319  points[nPoints++] = pts[i];
320  }
321  }
322  }
323 
324  // Merge
325  labelList oldToNew;
326  pointField newPoints;
327  const bool hasMerged = mergePoints
328  (
329  points,
330  small,
331  false,
332  oldToNew,
333  newPoints
334  );
335 
336  if (hasMerged)
337  {
338  points.transfer(newPoints);
339  forAll(faces, i)
340  {
341  inplaceRenumber(oldToNew, faces[i]);
342  }
343  }
344 }
345 
346 
347 void Foam::functionObjects::fieldValues::surfaceFieldValue::
348 combineSurfaceGeometry
349 (
350  faceList& faces,
352 ) const
353 {
354  if (selectionType_ == selectionTypes::sampledSurface)
355  {
356  const sampledSurface& s = surfacePtr_();
357 
358  if (Pstream::parRun())
359  {
360  // Dimension as fraction of mesh bounding box
361  scalar mergeDim = 1e-10*mesh_.bounds().mag();
362 
363  labelList pointsMap;
364 
366  (
367  mergeDim,
369  (
370  SubList<face>(s.faces(), s.faces().size()),
371  s.points()
372  ),
373  points,
374  faces,
375  pointsMap
376  );
377  }
378  else
379  {
380  faces = s.faces();
381  points = s.points();
382  }
383  }
384 }
385 
386 
387 Foam::scalar
388 Foam::functionObjects::fieldValues::surfaceFieldValue::area() const
389 {
390  if (selectionType_ == selectionTypes::sampledSurface)
391  {
392  return gSum(surfacePtr_().magSf());
393  }
394  else
395  {
396  return gSum(filterField(mesh_.magSf()));
397  }
398 }
399 
400 
401 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
402 
404 (
405  const label i
406 )
407 {
408  if (operation_ != operationType::none)
409  {
410  writeCommented(file(), "Selection");
411  file()
412  << setw(1) << ':' << setw(1) << ' '
413  << selectionTypeNames[selectionType_] << "("
414  << selectionName_.c_str() << ")" << endl;
415 
416  writeHeaderValue(file(), "Faces", nFaces_);
417 
418  writeHeaderValue(file(), "Area", area_);
419 
420  writeCommented(file(), "Time");
421 
422  if (writeNFaces_) file() << tab << "Faces";
423  if (writeArea_) file() << tab << "Area";
424 
425  forAll(fields_, fieldi)
426  {
427  file()
428  << tab << operationTypeNames_[operation_]
429  << "(" << fields_[fieldi] << ")";
430  }
431 
432  file() << endl;
433  }
434 }
435 
436 
438 (
439  const Field<scalar>& values,
440  const scalarField& signs,
441  const scalarField& weights,
442  const vectorField& Sf,
443  scalar& result
444 ) const
445 {
446  switch (operation_)
447  {
448  default:
449  {
450  // Fall through to same-type operations
451  return processValuesTypeType
452  (
453  values,
454  signs,
455  weights,
456  Sf,
457  result
458  );
459  }
460  }
461 }
462 
463 
465 (
466  const Field<vector>& values,
467  const scalarField& signs,
468  const scalarField& weights,
469  const vectorField& Sf,
470  scalar& result
471 ) const
472 {
473  switch (operation_)
474  {
475  case operationType::areaNormalAverage:
476  {
477  result = gSum(weights*values & Sf)/gSum(mag(weights*Sf));
478  return true;
479  }
480  case operationType::areaNormalIntegrate:
481  {
482  result = gSum(weights*values & Sf);
483  return true;
484  }
485  default:
486  {
487  // No fall through. Different types.
488  return false;
489  }
490  }
491 }
492 
493 
495 {
496  switch (selectionType_)
497  {
498  case selectionTypes::faceZone:
499  if (faceZonePtr_->zone().moveUpdate()) setFaceZoneFaces();
500  break;
501  case selectionTypes::patch:
503  break;
504  case selectionTypes::sampledSurface:
505  surfacePtr_->expire();
506  setSampledSurfaceFaces();
507  break;
508  }
509 
510  area_ = nFaces_ ? area() : NaN;
511 }
512 
513 
515 {
516  switch (selectionType_)
517  {
518  case selectionTypes::faceZone:
519  if (!faceZonePtr_->zone().moveUpdate()) setFaceZoneFaces();
520  break;
521  case selectionTypes::patch:
523  setPatchesFaces();
524  break;
525  case selectionTypes::sampledSurface:
526  break;
527  }
528 
529  moveMesh();
530 }
531 
532 
533 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
534 
536 (
537  const word& name,
538  const Time& runTime,
539  const dictionary& dict
540 )
541 :
542  fieldValue(name, runTime, dict, typeName),
543  surfaceWriterPtr_(nullptr),
544  selectionType_(selectionTypeNames.select(dict)),
545  selectionName_(string::null),
546  operation_(operationTypeNames_.read(dict.lookup("operation"))),
547  weightFieldNames_(),
548  nFaces_(0),
549  area_(NaN),
550  writeNFaces_(dict.lookupOrDefault("writeNumberOfFaces", false)),
551  writeArea_(dict.lookupOrDefault("writeArea", false)),
552  faceZonePtr_(nullptr),
553  patchNames_(),
554  faceId_(),
555  facePatchId_(),
556  faceSign_(),
557  surfacePtr_(nullptr)
558 {
559  read(dict);
560 }
561 
563 (
564  const word& name,
565  const objectRegistry& obr,
566  const dictionary& dict
567 )
568 :
569  fieldValue(name, obr, dict, typeName),
570  surfaceWriterPtr_(nullptr),
571  selectionType_(selectionTypeNames.select(dict)),
572  selectionName_(string::null),
573  operation_(operationTypeNames_.read(dict.lookup("operation"))),
574  weightFieldNames_(),
575  nFaces_(0),
576  area_(NaN),
577  writeNFaces_(dict.lookupOrDefault("writeNumberOfFaces", false)),
578  writeArea_(dict.lookupOrDefault("writeArea", false)),
579  faceZonePtr_(nullptr),
580  patchNames_(),
581  faceId_(),
582  facePatchId_(),
583  faceSign_(),
584  surfacePtr_(nullptr)
585 {
586  read(dict);
587 }
588 
589 
590 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
591 
593 {}
594 
595 
596 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
597 
599 (
600  const dictionary& dict
601 )
602 {
604 
605  const word selection(selectionTypeNames[selectionType_]);
606 
607  switch (selectionType_)
608  {
609  case selectionTypes::faceZone:
610  {
611  faceZonePtr_.reset(new generatedFaceZone(mesh_, dict));
612 
613  break;
614  }
615  case selectionTypes::patch:
616  {
617  const wordRe patchName = dict.lookup<wordRe>(selection);
618 
619  patchNames_ = wordReList(1, patchName);
620 
621  selectionName_ =
622  patchName.isPattern()
623  ? '"' + patchName + '"'
624  : string(patchName);
625 
626  break;
627  }
629  {
630  patchNames_ = dict.lookup<wordReList>(selection);
631 
632  const labelList patchis = this->patchis(patchNames_);
633 
634  selectionName_.clear();
635 
636  forAll(patchis, i)
637  {
638  const fvPatch& fvp = mesh_.boundary()[patchis[i]];
639  selectionName_.append((i ? " " : "") + fvp.name());
640  }
641 
642  break;
643  }
644  case selectionTypes::sampledSurface:
645  {
646  surfacePtr_ = sampledSurface::New(name(), mesh_, dict);
647 
648  selectionName_ = surfacePtr_().name();
649 
650  break;
651  }
652  }
653 
654  if (dict.found("weightFields"))
655  {
656  dict.lookup("weightFields") >> weightFieldNames_;
657  }
658  else if (dict.found("weightField"))
659  {
660  weightFieldNames_.setSize(1);
661 
662  dict.lookup("weightField") >> weightFieldNames_[0];
663  }
664 
665  if (writeFields_)
666  {
667  const word surfaceFormat(dict.lookup("surfaceFormat"));
668 
669  surfaceWriterPtr_.reset
670  (
671  surfaceWriter::New(surfaceFormat, dict).ptr()
672  );
673  }
674 
675  changeMesh();
676 
677  // Report configuration
678  Info<< type() << ' ' << name() << " read:" << nl;
679  Info<< " number of faces = " << nFaces_ << nl;
680  if (nFaces_)
681  {
682  Info<< " area = " << area_ << nl;
683  }
684  Info<< " operation = " << operationTypeNames_[operation_] << nl;
685  if (weightFieldNames_.size() == 1)
686  {
687  Info<< " weight field = " << weightFieldNames_[0] << nl;
688  }
689  if (weightFieldNames_.size() > 1)
690  {
691  Info<< " weight fields =";
692  forAll(weightFieldNames_, i) Info<< ' ' << weightFieldNames_[i];
693  Info<< nl;
694  }
695  Info<< endl;
696 
697  return true;
698 }
699 
700 
702 {
703  if (nFaces_ == 0)
704  {
705  return false;
706  }
707 
708  // Look to see if any fields exist. Use the flag to suppress output later.
709  bool anyFields = false;
710  forAll(fields_, i)
711  {
712  #define validFieldType(fieldType, none) \
713  anyFields = anyFields || validField<fieldType>(fields_[i]);
715  #undef validFieldType
716  }
717  if (!anyFields && fields_.size() > 1) // (error for 1 will happen below)
718  {
719  cannotFindObjects(fields_);
720  }
721 
722  // Initialise the file, write the header, etc...
723  if (anyFields && operation_ != operationType::none)
724  {
726  }
727 
728  // Update the surface
729  if (selectionType_ == selectionTypes::sampledSurface)
730  {
731  surfacePtr_().update();
732  }
733 
734  // Write the time
735  if (anyFields && operation_ != operationType::none && Pstream::master())
736  {
737  writeTime(file());
738  }
739 
740  // Write the number of faces and/or the area if necessary
741  if (anyFields && operation_ != operationType::none && Pstream::master())
742  {
743  if (writeNFaces_)
744  {
745  file() << tab << nFaces_;
746  }
747  if (writeArea_)
748  {
749  file() << tab << area_;
750  }
751  }
752  if (writeNFaces_)
753  {
754  Log << " number of faces = " << nFaces_ << endl;
755  }
756  if (writeArea_)
757  {
758  Log << " area = " << area_ << endl;
759  }
760 
761  // Construct the sign and weight fields and the surface normals
762  const scalarField signs
763  (
764  selectionType_ == selectionTypes::sampledSurface
765  ? scalarField(surfacePtr_().Sf().size(), 1)
766  : List<scalar>(faceSign_)
767  );
768  scalarField weights(signs.size(), 1);
769  forAll(weightFieldNames_, i)
770  {
771  weights *= getFieldValues<scalar>(weightFieldNames_[i]);
772  }
773  const vectorField Sf
774  (
775  selectionType_ == selectionTypes::sampledSurface
776  ? surfacePtr_().Sf()
777  : (signs*filterField(mesh_.Sf()))()
778  );
779 
780  // Create storage for the values
781  #define DeclareValues(fieldType, nullArg) \
782  PtrList<Field<fieldType>> fieldType##Values(fields_.size());
784  #undef DeclareValues
785 
786  // Process the fields in turn. Get the values and store if necessary.
787  // Compute the operation and write into the file and the log.
788  forAll(fields_, i)
789  {
790  const word& fieldName = fields_[i];
791  bool ok = false;
792 
793  #define writeValuesFieldType(fieldType, none) \
794  { \
795  const bool typeOk = validField<fieldType>(fieldName); \
796  \
797  if (typeOk) \
798  { \
799  tmp<Field<fieldType>> values = \
800  getFieldValues<fieldType>(fieldName); \
801  \
802  writeValues<fieldType> \
803  ( \
804  fieldName, \
805  values(), \
806  signs, \
807  weights, \
808  Sf \
809  ); \
810  \
811  if (writeFields_) \
812  { \
813  fieldType##Values.set \
814  ( \
815  i, \
816  getFieldValues<fieldType>(fieldName).ptr() \
817  ); \
818  } \
819  } \
820  \
821  ok = ok || typeOk; \
822  }
824  #undef writeValuesFieldType
825 
826  if (!ok)
827  {
828  cannotFindObject(fieldName);
829  }
830  }
831 
832  // Finalise the file and the log
833  if (anyFields && operation_ != operationType::none && Pstream::master())
834  {
835  file() << endl;
836  }
837  Log << endl;
838 
839  // Write a surface file with the values if specified
840  if (writeFields_)
841  {
842  faceList faces;
844 
845  if (selectionType_ == selectionTypes::sampledSurface)
846  {
847  combineSurfaceGeometry(faces, points);
848  }
849  else
850  {
851  combineMeshGeometry(faces, points);
852  }
853 
854  if (Pstream::master())
855  {
856  surfaceWriterPtr_->write
857  (
858  baseFileDir()/name()/time_.name(),
859  word(selectionTypeNames[selectionType_])
860  + "("
861  + (
862  selectionType_ == selectionTypes::patches
863  ? selectionName_.replace(" ", ",").c_str()
864  : selectionName_.c_str()
865  )
866  + ")",
867  points,
868  faces,
869  fields_,
870  false
871  #define ValuesParameter(fieldType, nullArg) \
872  , fieldType##Values
874  #undef ValuesParameter
875  );
876  }
877  }
878 
879  return true;
880 }
881 
882 
884 (
885  const polyMesh& mesh
886 )
887 {
888  if (&mesh == &this->mesh())
889  {
891  if (selectionType_ == selectionTypes::faceZone)
892  {
893  faceZonePtr_->movePoints();
894  }
895  moveMesh();
896  }
897 }
898 
899 
901 (
902  const polyTopoChangeMap& map
903 )
904 {
905  if (&map.mesh() == &mesh())
906  {
908  if (selectionType_ == selectionTypes::faceZone)
909  {
910  faceZonePtr_->topoChange(map);
911  }
912  changeMesh();
913  }
914 }
915 
916 
918 (
919  const polyMeshMap& map
920 )
921 {
922  if (&map.mesh() == &mesh())
923  {
924  fieldValue::mapMesh(map);
925  if (selectionType_ == selectionTypes::faceZone)
926  {
927  faceZonePtr_->mapMesh(map);
928  }
929  changeMesh();
930  }
931 }
932 
933 
935 (
936  const polyDistributionMap& map
937 )
938 {
939  if (&map.mesh() == &mesh())
940  {
942  if (selectionType_ == selectionTypes::faceZone)
943  {
944  faceZonePtr_->distribute(map);
945  }
946  changeMesh();
947  }
948 }
949 
950 
951 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
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:55
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.
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base-class for Time/database functionObjects.
virtual const word & type() const =0
Runtime type information.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual void movePoints(const polyMesh &mesh)
Update topology using the given map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
const word & name() const
Return the name of this functionObject.
virtual void topoChange(const polyTopoChangeMap &map)
Update topology using the given map.
Base class for field value -based function objects.
Definition: fieldValue.H:66
const dictionary & dict() const
Return the reference to the construction dictionary.
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: fieldValue.C:86
virtual bool write()
Write.
Definition: fieldValue.C:105
Surface, patch or faceZone selection class.
surfaceFieldValue(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
void changeMesh()
Update the surface following mesh change.
static const NamedEnum< selectionTypes, 4 > selectionTypeNames
Selection type names.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
void moveMesh()
Update the surface following mesh motion.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
static const NamedEnum< operationType, 15 > operationTypeNames_
Operation type names.
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.
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.
virtual bool read(const dictionary &)
Read from dictionary.
const fvMesh & mesh_
Reference to the fvMesh.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual const word & name() const
Return name.
Definition: fvPatch.H:126
faceZone selection or generation class
Registry of regIOobjects.
wordList names() const
Return the list of patch names.
labelList findIndices(const wordRe &, const bool usePatchGroups=true) const
Return patch indices for all matches. Optionally matches patchGroups.
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 polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const polyMesh & mesh() const
Return polyMesh.
static autoPtr< sampledSurface > New(const word &name, const polyMesh &, const dictionary &)
Return a reference to the selected surface.
A class for handling character strings derived from std::string.
Definition: string.H:79
static autoPtr< surfaceWriter > New(const word &writeType, const IOstream::streamFormat writeFormat, const IOstream::compressionType writeCompression)
Select given write options.
Definition: surfaceWriter.C:75
A wordRe is a word, but can also have a regular expression for matching words.
Definition: wordRe.H:77
static bool isPattern(const string &)
Test string for regular expression meta characters.
Definition: wordReI.H:34
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#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(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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:258
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
static const char tab
Definition: Ostream.H:266
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
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.
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
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:267
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
wordList patchNames(nPatches)
labelList f(nPoints)
dictionary dict
#define ValuesParameter(fieldType, nullArg)
#define DeclareValues(fieldType, nullArg)
#define writeValuesFieldType(fieldType, none)
#define validFieldType(fieldType, none)