vtkPVFoamUpdateInfo.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 "vtkPVFoam.H"
27 #include "vtkPVFoamReader.h"
29 
30 // OpenFOAM includes
31 #include "processorRunTimes.H"
32 #include "domainDecomposition.H"
33 #include "cellSet.H"
34 #include "faceSet.H"
35 #include "pointSet.H"
36 #include "IOobjectList.H"
38 #include "entry.H"
39 #include "cloud.H"
40 #include "pointMesh.H"
41 #include "LagrangianMesh.H"
42 
43 // VTK includes
44 #include "vtkDataArraySelection.h"
45 
46 // * * * * * * * * * * * * * * * Private Classes * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 //- A class for reading zone information without requiring a mesh
52 template<class ZonesType>
53 class zonesEntries
54 :
55  public regIOobject,
56  public PtrList<entry>
57 {
58 
59 public:
60 
61  // Constructors
62 
63  explicit zonesEntries(const IOobject& io)
64  :
65  regIOobject(io),
66  PtrList<entry>(readStream(ZonesType::typeName))
67  {
68  close();
69  }
70 
71  // Member functions
72 
73  bool writeData(Ostream&) const
74  {
76  return false;
77  }
78 };
79 
80 }
81 
82 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
83 
84 template<class ZonesType>
85 Foam::wordList Foam::vtkPVFoam::getZoneNames(const ZonesType& zones) const
86 {
87  wordList names(zones.size());
88  label nZone = 0;
89 
90  forAll(zones, zoneI)
91  {
92  if (zones[zoneI].size())
93  {
94  names[nZone++] = zones[zoneI].name();
95  }
96  }
97  names.setSize(nZone);
98 
99  return names;
100 }
101 
102 
103 template<class ZonesType>
104 Foam::wordList Foam::vtkPVFoam::getZoneNames(const word& zonesName) const
105 {
106  const Time& runTime =
107  reader_->GetDecomposedCase()
108  ? procDbsPtr_->proc0Time()
109  : procDbsPtr_->completeTime();
110 
111  // Mesh not loaded - read from file
112  typeIOobject<ZonesType> ioObj
113  (
114  zonesName,
115  runTime.findInstance
116  (
117  meshDir_,
118  zonesName,
120  ),
121  meshDir_,
122  runTime,
125  false
126  );
127 
128  if (ioObj.headerOk())
129  {
130  const zonesEntries<ZonesType> zones(ioObj);
131 
132  wordList names(zones.size());
133 
134  forAll(zones, zoneI)
135  {
136  names[zoneI] = zones[zoneI].keyword();
137  }
138 
139  return names;
140  }
141  else
142  {
143  return wordList();
144  }
145 }
146 
147 
148 void Foam::vtkPVFoam::updateInfoInternalMesh
149 (
150  vtkDataArraySelection* arraySelection
151 )
152 {
153  // Determine mesh parts (internalMesh, patches...)
154  // - Add internal mesh as first entry
155  arrayRangeVolume_.reset(arraySelection->GetNumberOfArrays());
156  arraySelection->AddArray
157  (
158  "internalMesh"
159  );
160  arrayRangeVolume_ += 1;
161 }
162 
163 
164 void Foam::vtkPVFoam::updateInfolagrangian
165 (
166  vtkDataArraySelection* arraySelection
167 )
168 {
169  UPtrList<const Time> runTimes;
170  if (reader_->GetDecomposedCase())
171  {
172  runTimes.setSize(procDbsPtr_->nProcs());
173  forAll(procDbsPtr_->procTimes(), proci)
174  {
175  runTimes.set(proci, &procDbsPtr_->procTimes()[proci]);
176  }
177  }
178  else
179  {
180  runTimes.setSize(1);
181  runTimes.set(0, &procDbsPtr_->completeTime());
182  }
183 
184  // Use the db directly since this might be called without a mesh,
185  // but the region must get added back in
186  const fileName lagrangianPrefix =
187  meshRegion_ == polyMesh::defaultRegion
188  ? fileName(lagrangian::cloud::prefix)
189  : meshRegion_/lagrangian::cloud::prefix;
190 
191  arrayRangelagrangian_.reset(arraySelection->GetNumberOfArrays());
192 
193  // Generate a list of lagrangian clouds across all times
194  HashSet<fileName> lagrangianDirs;
195 
196  // Get times list. Flush first to force refresh.
197  fileHandler().flush();
198  forAll(runTimes, runTimei)
199  {
200  const instantList times = runTimes[runTimei].times();
201 
202  forAll(times, timei)
203  {
204  lagrangianDirs +=
206  (
207  fileHandler().filePath
208  (
209  runTimes[runTimei].path()
210  /times[timei].name()
211  /lagrangianPrefix
212  ),
214  );
215  }
216  }
217 
218  forAllConstIter(HashSet<fileName>, lagrangianDirs, cloudIter)
219  {
220  arraySelection->AddArray
221  (
222  (cloudIter.key() + " - lagrangian").c_str()
223  );
224  }
225 
226  arrayRangelagrangian_ += lagrangianDirs.size();
227 }
228 
229 
230 void Foam::vtkPVFoam::updateInfoLagrangian
231 (
232  vtkDataArraySelection* arraySelection
233 )
234 {
235  const Time& runTime =
236  reader_->GetDecomposedCase()
237  ? procDbsPtr_->proc0Time()
238  : procDbsPtr_->completeTime();
239 
240  // Use the db directly since this might be called without a mesh,
241  // but the region must get added back in
242  const fileName LagrangianPrefix =
243  meshRegion_ == polyMesh::defaultRegion
244  ? fileName(LagrangianMesh::prefix)
245  : meshRegion_/LagrangianMesh::prefix;
246 
247  arrayRangeLagrangian_.reset(arraySelection->GetNumberOfArrays());
248 
249  // Generate a list of Lagrangian meshes across all times
250  HashSet<fileName> LagrangianDirs;
251 
252  // Get times list. Flush first to force refresh.
253  fileHandler().flush();
254  instantList times = runTime.times();
255  forAll(times, timei)
256  {
257  LagrangianDirs +=
259  (
260  fileHandler().filePath
261  (
262  runTime.path()
263  /times[timei].name()
264  /LagrangianPrefix
265  ),
267  );
268  }
269 
270  forAllConstIter(HashSet<fileName>, LagrangianDirs, LagrangianIter)
271  {
272  arraySelection->AddArray
273  (
274  (LagrangianIter.key() + " - Lagrangian").c_str()
275  );
276  }
277 
278  arrayRangeLagrangian_ += LagrangianDirs.size();
279 }
280 
281 
282 void Foam::vtkPVFoam::updateInfoPatches
283 (
284  vtkDataArraySelection* arraySelection,
285  stringList& enabledEntries,
286  const bool first
287 )
288 {
289  HashSet<string> enabledEntriesSet(enabledEntries);
290 
291  arrayRangePatches_.reset(arraySelection->GetNumberOfArrays());
292 
293  int nPatches = 0;
294  if (procMeshesPtr_.valid())
295  {
296  const fvMesh& mesh = procMeshesPtr_->completeMesh();
297 
298  const polyBoundaryMesh& patches = mesh.boundaryMesh();
299  const HashTable<labelList, word>& groups = patches.groupPatchIndices();
300  const wordList allPatchNames = patches.names();
301 
302  // Add patch groups
303  for
304  (
305  HashTable<labelList, word>::const_iterator iter = groups.begin();
306  iter != groups.end();
307  ++iter
308  )
309  {
310  const word& groupName = iter.key();
311  const labelList& patchIDs = iter();
312 
313  // Count the number of faces in this group
314  label nFaces = 0;
315  forAll(patchIDs, i)
316  {
317  nFaces += patches[patchIDs[i]].size();
318  }
319 
320  // Valid patch if nFace > 0 - add patch to GUI list
321  if (nFaces)
322  {
323  string vtkGrpName = groupName + " - group";
324  arraySelection->AddArray(vtkGrpName.c_str());
325 
326  ++nPatches;
327 
328  if (enabledEntriesSet.found(vtkGrpName))
329  {
330  if (!reader_->GetShowGroupsOnly())
331  {
332  enabledEntriesSet.erase(vtkGrpName);
333  forAll(patchIDs, i)
334  {
335  const polyPatch& pp = patches[patchIDs[i]];
336 
337  if (pp.size())
338  {
339  enabledEntriesSet.insert
340  (
341  pp.name() + " - " + pp.type()
342  );
343  }
344  }
345  }
346  }
347  }
348  }
349 
350  // Add patches
351  if (!reader_->GetShowGroupsOnly())
352  {
354  {
355  const polyPatch& pp = patches[patchi];
356 
357  if (pp.size())
358  {
359  arraySelection->AddArray
360  (
361  (pp.name() + " - " + pp.type()).c_str()
362  );
363 
364  ++nPatches;
365  }
366  }
367  }
368  }
369  else
370  {
371  const bool decomposed = reader_->GetDecomposedCase();
372 
373  const Time& runTime =
374  decomposed
375  ? procDbsPtr_->proc0Time()
376  : procDbsPtr_->completeTime();
377 
378  // Mesh not loaded - read from file
379  // but this could fail if we've supplied a bad region name
380  typeIOobject<polyBoundaryMesh> ioObj
381  (
382  "boundary",
383  runTime.findInstance
384  (
385  meshDir_,
386  "boundary",
388  ),
389  meshDir_,
390  runTime,
393  false
394  );
395 
396  // This should only ever fail if the mesh region doesn't exist
397  if (ioObj.headerOk())
398  {
399  polyBoundaryMeshEntries patchEntries(ioObj);
400 
401  // Read patches and determine sizes
402  wordList names(patchEntries.size());
403  labelList sizes(patchEntries.size());
404  boolList isProc(patchEntries.size());
405  forAll(patchEntries, patchi)
406  {
407  const dictionary& patchDict = patchEntries[patchi].dict();
408 
409  names[patchi] = patchEntries[patchi].keyword();
410  sizes[patchi] = patchDict.lookup<label>("nFaces");
411  isProc[patchi] = patchDict.found("myProcNo");
412  }
413 
414  // Build a map from the group name to a list of the indices of the
415  // patches in the group
416  HashTable<labelList> groups(patchEntries.size());
417  forAll(patchEntries, patchi)
418  {
419  const dictionary& patchDict = patchEntries[patchi].dict();
420 
421  const wordList groupNames =
422  patchDict.lookupOrDefault("inGroups", wordList());
423 
424  forAll(groupNames, groupI)
425  {
427  groups.find(groupNames[groupI]);
428 
429  if (iter != groups.end())
430  {
431  iter().append(patchi);
432  }
433  else
434  {
435  groups.insert(groupNames[groupI], labelList(1, patchi));
436  }
437  }
438  }
439 
440  // Add (non-zero) patch groups to the list of mesh parts
441  forAllConstIter(HashTable<labelList>, groups, iter)
442  {
443  const word& groupName = iter.key();
444  const labelList& patchIDs = iter();
445 
446  label nFaces = 0;
447  forAll(patchIDs, i)
448  {
449  if (!isProc[patchIDs[i]])
450  {
451  nFaces += sizes[patchIDs[i]];
452  }
453  }
454 
455  // Valid patch if nFace > 0 - add patch to GUI list
456  if (nFaces)
457  {
458  string vtkGrpName = groupName + " - group";
459  arraySelection->AddArray(vtkGrpName.c_str());
460 
461  ++nPatches;
462 
463  if (enabledEntriesSet.found(vtkGrpName))
464  {
465  if (!reader_->GetShowGroupsOnly())
466  {
467  enabledEntriesSet.erase(vtkGrpName);
468 
469  forAll(patchIDs, i)
470  {
471  if (sizes[patchIDs[i]])
472  {
473  enabledEntriesSet.insert
474  (
475  names[patchIDs[i]]
476  + " - "
477  + patchEntries[patchIDs[i]]
478  .dict()
479  .lookup<word>("type")
480  );
481  }
482  }
483  }
484  }
485  }
486  }
487 
488  // Add (non-zero) patches to the list of mesh parts
489  if (!reader_->GetShowGroupsOnly())
490  {
491  const wordReList defaultPatchTypes
492  (
493  configDict_.lookupOrDefault
494  (
495  "defaultPatchTypes",
496  wordReList(wordList({"patch", "wall"}))
497  )
498  );
499 
500  forAll(names, patchi)
501  {
502  // Valid patch if nFace > 0 - add patch to GUI list
503  if (sizes[patchi])
504  {
505  const word patchType
506  (
507  patchEntries[patchi].dict().lookup("type")
508  );
509 
510  const string vtkPatchName
511  (
512  names[patchi] + " - " + patchType
513  );
514 
515  arraySelection->AddArray(vtkPatchName.c_str());
516 
517  if (first && findStrings(defaultPatchTypes, patchType))
518  {
519  enabledEntriesSet.insert(vtkPatchName);
520  }
521 
522  ++nPatches;
523  }
524  }
525  }
526  }
527  }
528 
529  arrayRangePatches_ += nPatches;
530 
531  // Update enabled entries in case of group selection
532  enabledEntries = enabledEntriesSet.toc();
533 }
534 
535 
536 void Foam::vtkPVFoam::updateInfoZones
537 (
538  vtkDataArraySelection* arraySelection
539 )
540 {
541  if (!reader_->GetIncludeZones()) return;
542 
543  const fvMesh& mesh =
544  !procMeshesPtr_.valid()
545  ? NullObjectRef<fvMesh>()
546  : reader_->GetDecomposedCase()
547  ? procMeshesPtr_->haveProcs()
548  ? procMeshesPtr_->procMeshes().first()
549  : NullObjectRef<fvMesh>()
550  : procMeshesPtr_->completeMesh();
551 
552  // cellZones information
553  const wordList cellZoneNames =
554  notNull(mesh)
555  ? getZoneNames(mesh.cellZones())
556  : getZoneNames<cellZoneList>("cellZones");
557 
558  arrayRangeCellZones_.reset(arraySelection->GetNumberOfArrays());
559  forAll(cellZoneNames, i)
560  {
561  arraySelection->AddArray
562  (
563  (cellZoneNames[i] + " - cellZone").c_str()
564  );
565  }
566  arrayRangeCellZones_ += cellZoneNames.size();
567 
568  // faceZones information
569  const wordList faceZoneNames =
570  notNull(mesh)
571  ? getZoneNames(mesh.faceZones())
572  : getZoneNames<faceZoneList>("faceZones");
573 
574  arrayRangeFaceZones_.reset(arraySelection->GetNumberOfArrays());
575  forAll(faceZoneNames, i)
576  {
577  arraySelection->AddArray
578  (
579  (faceZoneNames[i] + " - faceZone").c_str()
580  );
581  }
582  arrayRangeFaceZones_ += faceZoneNames.size();
583 
584  // pointZones information
585  const wordList pointZoneNames =
586  notNull(mesh)
587  ? getZoneNames(mesh.pointZones())
588  : getZoneNames<pointZoneList>("pointZones");
589 
590  arrayRangePointZones_.reset(arraySelection->GetNumberOfArrays());
591  forAll(pointZoneNames, i)
592  {
593  arraySelection->AddArray
594  (
595  (pointZoneNames[i] + " - pointZone").c_str()
596  );
597  }
598  arrayRangePointZones_ += pointZoneNames.size();
599 }
600 
601 
602 void Foam::vtkPVFoam::updateInfoSets
603 (
604  vtkDataArraySelection* arraySelection
605 )
606 {
607  if (!reader_->GetIncludeSets()) return;
608 
609  const Time& runTime =
610  reader_->GetDecomposedCase()
611  ? procDbsPtr_->proc0Time()
612  : procDbsPtr_->completeTime();
613 
614  // Add names of sets. Search for last time directory with a sets
615  // subdirectory. Take care not to search beyond the last mesh.
616 
617  const word facesInstance = runTime.findInstance
618  (
619  meshDir_,
620  "faces",
622  );
623 
624  const word setsInstance = runTime.findInstance
625  (
626  meshDir_/"sets",
627  word::null,
629  facesInstance
630  );
631 
632  const IOobjectList objects(runTime, setsInstance, meshDir_/"sets");
633 
634  arrayRangeCellSets_.reset(arraySelection->GetNumberOfArrays());
635  arrayRangeCellSets_ +=
636  addToSelection<cellSet>(arraySelection, objects, " - cellSet");
637 
638  arrayRangeFaceSets_.reset(arraySelection->GetNumberOfArrays());
639  arrayRangeFaceSets_ +=
640  addToSelection<faceSet>(arraySelection, objects, " - faceSet");
641 
642  arrayRangePointSets_.reset(arraySelection->GetNumberOfArrays());
643  arrayRangePointSets_ +=
644  addToSelection<pointSet>(arraySelection, objects, " - pointSet");
645 }
646 
647 
648 void Foam::vtkPVFoam::updateInfoFields()
649 {
651 
652  vtkDataArraySelection* fieldSelection = reader_->GetFieldSelection();
653 
654  // Use the db directly since this might be called without a mesh,
655  // but the region must get added back in
656  word regionPrefix;
657  if (meshRegion_ != polyMesh::defaultRegion)
658  {
659  regionPrefix = meshRegion_;
660  }
661 
662  const Time& runTime =
663  reader_->GetDecomposedCase()
664  ? procDbsPtr_->proc0Time()
665  : procDbsPtr_->completeTime();
666 
667  const instantList times = runTime.times();
668 
669  stringList enabledEntries;
670  if (fieldSelection->GetNumberOfArrays() == 0 && !procMeshesPtr_.valid())
671  {
672  const wordReList defaultFieldRes
673  (
674  configDict_.lookupOrDefault
675  (
676  "defaultFields",
677  wordReList(wordList{"p", "U"})
678  )
679  );
680 
681  wordHashSet objectNameSet;
682  forAll(times, timei)
683  {
684  // Search for list of objects for this time and mesh region
685  IOobjectList objects(runTime, times[timei].name(), regionPrefix);
686 
687  forAllConstIter(IOobjectList, objects, iter)
688  {
689  objectNameSet.insert(iter.key());
690  }
691  }
692 
693  const wordList objectNames(objectNameSet.toc());
694  const labelList defaultFields
695  (
696  findStrings(defaultFieldRes, objectNames)
697  );
698 
699  enabledEntries.setSize(defaultFields.size());
700  forAll(defaultFields, i)
701  {
702  enabledEntries[i] = objectNames[defaultFields[i]];
703  }
704  }
705  else
706  {
707  // Preserve the enabled selections
708  enabledEntries = getSelectedArrayEntries(fieldSelection, false);
709  }
710 
711  fieldSelection->RemoveAllArrays();
712 
713  forAll(times, timei)
714  {
715  // Search for list of objects for this time and mesh region
716  IOobjectList objects(runTime, times[timei].name(), regionPrefix);
717 
718  addFieldsToSelection<volMesh>(fieldSelection, objects);
719  addInternalFieldsToSelection<volMesh>(fieldSelection, objects);
720  addFieldsToSelection<surfaceMesh>(fieldSelection, objects);
721  addFieldsToSelection<pointMesh>(fieldSelection, objects);
722  }
723 
724  // Restore the enabled selections
725  setSelectedArrayEntries(fieldSelection, enabledEntries);
726 
727  if (debug) getSelectedArrayEntries(fieldSelection);
728 }
729 
730 
731 void Foam::vtkPVFoam::updateInfolagrangianFields()
732 {
734 
735  UPtrList<const Time> runTimes;
736  if (reader_->GetDecomposedCase())
737  {
738  runTimes.setSize(procDbsPtr_->nProcs());
739  forAll(procDbsPtr_->procTimes(), proci)
740  {
741  runTimes.set(proci, &procDbsPtr_->procTimes()[proci]);
742  }
743  }
744  else
745  {
746  runTimes.setSize(1);
747  runTimes.set(0, &procDbsPtr_->completeTime());
748  }
749 
750  vtkDataArraySelection* fieldSelection =
751  reader_->GetlagrangianFieldSelection();
752 
753  // Preserve the enabled selections
754  stringList enabledEntries = getSelectedArrayEntries(fieldSelection, false);
755  fieldSelection->RemoveAllArrays();
756 
757  // Use the db directly since this might be called without a mesh,
758  // but the region must get added back in
759  fileName lagrangianPrefix(lagrangian::cloud::prefix);
760  if (meshRegion_ != polyMesh::defaultRegion)
761  {
762  lagrangianPrefix = meshRegion_/lagrangian::cloud::prefix;
763  }
764 
765  // Add the available fields from all clouds and all time directories.
766  // Differing sets of fields from multiple clouds get combined into a single
767  // set. ParaView will display "(partial)" after field names that only apply
768  // to some of the clouds.
769  const arrayRange& range = arrayRangelagrangian_;
770  fileHandler().flush();
771  for (label partId = range.start(); partId < range.end(); ++ partId)
772  {
773  forAll(runTimes, runTimei)
774  {
775  const instantList times = runTimes[runTimei].times();
776 
777  forAll(times, timei)
778  {
779  IOobjectList objects
780  (
781  runTimes[runTimei],
782  times[timei].name(),
783  lagrangianPrefix/getPartName(partId)
784  );
785 
786  #define ADD_TO_SELECTION(Type, nullArg) \
787  addToSelection<IOField<Type>>(fieldSelection, objects);
790  #undef ADD_TO_SELECTION
791  }
792  }
793  }
794 
795  // Restore the enabled selections
796  setSelectedArrayEntries(fieldSelection, enabledEntries);
797 
798  if (debug) getSelectedArrayEntries(fieldSelection);
799 }
800 
801 
802 void Foam::vtkPVFoam::updateInfoLagrangianFields()
803 {
805 
806  const Time& runTime =
807  reader_->GetDecomposedCase()
808  ? procDbsPtr_->proc0Time()
809  : procDbsPtr_->completeTime();
810 
811  const instantList times = runTime.times();
812 
813  vtkDataArraySelection* fieldSelection =
814  reader_->GetLagrangianFieldSelection();
815 
816  // Preserve the enabled selections
817  stringList enabledEntries = getSelectedArrayEntries(fieldSelection, false);
818  fieldSelection->RemoveAllArrays();
819 
820  // Use the db directly since this might be called without a mesh,
821  // but the region must get added back in
822  const fileName LagrangianPrefix =
823  meshRegion_ == polyMesh::defaultRegion
824  ? fileName(LagrangianMesh::prefix)
825  : meshRegion_/LagrangianMesh::prefix;
826 
827  // Add the available fields from all clouds and all time directories.
828  // Differing sets of fields from multiple clouds get combined into a single
829  // set. ParaView will display "(partial)" after field names that only apply
830  // to some of the clouds.
831  const arrayRange& range = arrayRangeLagrangian_;
832 
833  fileHandler().flush();
834 
835  for (label partId = range.start(); partId < range.end(); ++ partId)
836  {
837  forAll(times, timei)
838  {
839  IOobjectList objects
840  (
841  runTime,
842  times[timei].name(),
843  LagrangianPrefix/getPartName(partId)
844  );
845 
846  #define ADD_TO_SELECTION(Type, GeoField) \
847  addToSelection<GeoField<Type>>(fieldSelection, objects);
852  #undef ADD_TO_SELECTION
853  }
854  }
855 
856  // Restore the enabled selections
857  setSelectedArrayEntries(fieldSelection, enabledEntries);
858 
859  if (debug) getSelectedArrayEntries(fieldSelection);
860 }
861 
862 
863 // ************************************************************************* //
scalar range
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
friend class iterator
Declare friendship with the iterator.
Definition: HashTable.H:194
friend class const_iterator
Declare friendship with the const_iterator.
Definition: HashTable.H:197
static const word prefix
Instance prefix.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
virtual void flush() const
Forcibly wait until all output done. Flush any cached data.
virtual fileNameList readDir(const fileName &, const fileType=fileType::file, const bool filterVariants=true, const bool followLink=true) const =0
Read a directory and return the entries as a string list.
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:69
const pointZoneList & pointZones() const
Return point zones.
Definition: polyMesh.H:437
const cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:449
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:270
const faceZoneList & faceZones() const
Return face zones.
Definition: polyMesh.H:443
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
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 NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
label patchi
const fvPatchList & patches
#define DebugInFunction
Report an information message using Foam::Info.
Namespace for OpenFOAM.
const fileOperation & fileHandler()
Get current file handler.
List< word > wordList
A List of words.
Definition: fileName.H:54
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
List< instant > instantList
List of instants.
Definition: instantList.H:42
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:64
DimensionedField< Type, LagrangianMesh > LagrangianInternalField
const T & NullObjectRef()
Return const reference to the nullObject of type T.
Definition: nullObjectI.H:27
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.
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:210
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
List< string > stringList
A List of strings.
Definition: stringList.H:50
GeometricField< Type, LagrangianMesh > LagrangianField
label nPatches
Definition: readKivaGrid.H:396
objects
dictionary dict
const bool writeData(readBool(pdfDictionary.lookup("writeData")))
#define ADD_TO_SELECTION(Type, nullArg)