vtkPVFoam.H
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 Class
25  Foam::vtkPVFoam
26 
27 Description
28  Provides a reader interface for OpenFOAM to VTK interaction.
29 
30 SourceFiles
31  vtkPVFoam.C
32  vtkPVFoam.H
33  vtkPVFoamFields.C
34  vtkPVFoamMesh.C
35  vtkPVFoamMeshLagrangian.C
36  vtkPVFoamTemplates.C
37  vtkPVFoamMeshSet.C
38  vtkPVFoamMeshVolume.C
39  vtkPVFoamMeshZone.C
40  vtkPVFoamSurfaceField.H
41  vtkPVFoamLagrangianFields.H
42  vtkPVFoamPatchField.H
43  vtkPVFoamPointFields.H
44  vtkPVFoamPoints.H
45  vtkPVFoamUpdateInfo.C
46  vtkPVFoamUtils.C
47  vtkPVFoamVolFields.H
48  vtkPVFoamAddToSelection.H
49 
50 \*---------------------------------------------------------------------------*/
51 
52 #ifndef vtkPVFoam_H
53 #define vtkPVFoam_H
54 
55 // do not include legacy strstream headers
56 #ifndef VTK_EXCLUDE_STRSTREAM_HEADERS
57  #define VTK_EXCLUDE_STRSTREAM_HEADERS
58 #endif
59 
60 #include "className.H"
61 #include "fileName.H"
62 #include "stringList.H"
63 #include "wordList.H"
64 #include "primitivePatch.H"
66 #include "dictionary.H"
67 #include "PtrList.H"
68 #include "HashSet.H"
69 #include "IOField.H"
70 #include "volFieldsFwd.H"
71 #include "surfaceFieldsFwd.H"
72 #include "pointFieldsFwd.H"
73 
74 // * * * * * * * * * * * * * Forward Declarations * * * * * * * * * * * * * //
75 
76 class vtkDataArraySelection;
77 class vtkDataSet;
78 class vtkPoints;
79 class vtkPVFoamReader;
80 class vtkRenderer;
81 class vtkTextActor;
82 class vtkMultiBlockDataSet;
83 class vtkPolyData;
84 class vtkUnstructuredGrid;
85 class vtkIndent;
86 
87 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
88 
89 namespace Foam
90 {
91 
92 // OpenFOAM class forward declarations
93 class argList;
94 class Time;
95 class processorRunTimes;
96 class domainDecomposition;
97 class fvMesh;
98 class IOobjectList;
99 class polyPatch;
100 class faceSet;
101 class pointSet;
102 class fvFieldReconstructor;
103 class pointFieldReconstructor;
104 class lagrangianFieldReconstructor;
105 class LagrangianMesh;
106 class LagrangianFieldReconstructor;
107 
108 /*---------------------------------------------------------------------------*\
109  Class vtkPVFoam Declaration
110 \*---------------------------------------------------------------------------*/
111 
112 class vtkPVFoam
113 {
114  // Private classes
115 
116  //- Bookkeeping for GUI checklists and the multi-block organisation
117  class arrayRange
118  {
119  const char *name_;
120  int block_;
121  int start_;
122  int size_;
123 
124  public:
125 
126  arrayRange(const char *name, const int blockNo=0)
127  :
128  name_(name),
129  block_(blockNo),
130  start_(0),
131  size_(0)
132  {}
133 
134  //- Return the block holding these datasets
135  int block() const
136  {
137  return block_;
138  }
139 
140  //- Assign block number, return previous value
141  int block(int blockNo)
142  {
143  int prev = block_;
144  block_ = blockNo;
145  return prev;
146  }
147 
148  //- Return block name
149  const char* name() const
150  {
151  return name_;
152  }
153 
154  //- Return array start index
155  int start() const
156  {
157  return start_;
158  }
159 
160  //- Return array end index
161  int end() const
162  {
163  return start_ + size_;
164  }
165 
166  //- Return sublist size
167  int size() const
168  {
169  return size_;
170  }
171 
172  bool empty() const
173  {
174  return !size_;
175  }
176 
177  //- Reset the size to zero and optionally assign a new start
178  void reset(const int startAt = 0)
179  {
180  start_ = startAt;
181  size_ = 0;
182  }
183 
184  //- Increment the size
185  void operator+=(const int n)
186  {
187  size_ += n;
188  }
189  };
190 
191  //- Bookkeeping for polyhedral cell decomposition
192  // hide in extra pointMap (cellSet/cellZone) for now
193  class polyDecomp
194  {
195  labelList superCells_;
196  labelList addPointCellLabels_;
197  labelList pointMap_;
198 
199  public:
200 
201  polyDecomp()
202  {}
203 
204  //- Label of original cell for decomposed cells
205  labelList& superCells()
206  {
207  return superCells_;
208  }
209 
210  //- Label of original cell for decomposed cells
211  const labelList& superCells() const
212  {
213  return superCells_;
214  }
215 
216  //- Cell-centre labels for additional points of decomposed cells
217  labelList& addPointCellLabels()
218  {
219  return addPointCellLabels_;
220  }
221 
222  //- Cell-centre labels for additional points of decomposed cells
223  const labelList& addPointCellLabels() const
224  {
225  return addPointCellLabels_;
226  }
227 
228  //- Point labels for subsetted meshes
229  labelList& pointMap()
230  {
231  return pointMap_;
232  }
233 
234  //- Point labels for subsetted meshes
235  const labelList& pointMap() const
236  {
237  return pointMap_;
238  }
239 
240 
241  //- Clear
242  void clear()
243  {
244  superCells_.clear();
245  addPointCellLabels_.clear();
246  pointMap_.clear();
247  }
248  };
249 
250 
251  // Private Data
252 
253  //- Access to the controlling vtkPVFoamReader
254  vtkPVFoamReader* reader_;
255 
256  //- Configuration dictionary
257  dictionary configDict_;
258 
259  //- OpenFOAM time controls
260  autoPtr<processorRunTimes> procDbsPtr_;
261 
262  //- OpenFOAM meshes
263  autoPtr<domainDecomposition> procMeshesPtr_;
264 
265  //- The mesh in the optional meshes sub-directory
266  word meshMesh_;
267 
268  //- Path to the mesh in the optional meshes sub-directory
269  word meshPath_;
270 
271  //- The mesh region
272  word meshRegion_;
273 
274  //- The mesh directory
275  fileName meshDir_;
276 
277  //- The time index
278  int timeIndex_;
279 
280  //- Selected geometrical parts (internalMesh, patches, ...)
281  boolList partStatus_;
282 
283  //- Datasets corresponding to selected geometrical pieces
284  // a negative number indicates that no vtkmesh exists for this piece
285  labelList partDataset_;
286 
287  //- First instance and size of various mesh parts
288  // used to index into partStatus_ and partDataset_
289  arrayRange arrayRangeVolume_;
290  arrayRange arrayRangePatches_;
291  arrayRange arrayRangelagrangian_;
292  arrayRange arrayRangeLagrangian_;
293  arrayRange arrayRangeCellZones_;
294  arrayRange arrayRangeFaceZones_;
295  arrayRange arrayRangePointZones_;
296  arrayRange arrayRangeCellSets_;
297  arrayRange arrayRangeFaceSets_;
298  arrayRange arrayRangePointSets_;
299 
300  //- Decomposed cells information (mesh regions)
301  List<polyDecomp> regionPolyDecomp_;
302 
303  //- Decomposed cells information (cellZones)
304  List<polyDecomp> zonePolyDecomp_;
305 
306  //- Decomposed cells information (cellSets)
307  List<polyDecomp> setPolyDecomp_;
308 
309  //- List of patch names for rendering to window
310  List<vtkTextActor*> patchTextActorsPtrs_;
311 
312  //- Finite-volume field reconstructor
313  autoPtr<fvFieldReconstructor> fvReconstructorPtr_;
314 
315  //- Point field reconstructor
316  autoPtr<pointFieldReconstructor> pointReconstructorPtr_;
317 
318  //- lagrangian field reconstructors
319  PtrList<lagrangianFieldReconstructor> lagrangianReconstructors_;
320 
321  //- Lagrangian meshes
322  PtrList<LagrangianMesh> LagrangianMeshes_;
323 
324  //- Lagrangian field reconstructors
325  PtrList<LagrangianFieldReconstructor> LagrangianReconstructors_;
326 
327 
328  // Private Member Functions
329 
330  // Convenience method use to convert the readers from VTK 5
331  // multiblock API to the current composite data infrastructure
332  static void AddToBlock
333  (
334  vtkMultiBlockDataSet* output,
335  vtkDataSet* dataset,
336  const arrayRange&,
337  const label datasetNo,
338  const std::string& datasetName
339  );
340 
341  // Convenience method use to convert the readers from VTK 5
342  // multiblock API to the current composite data infrastructure
343  static vtkDataSet* GetDataSetFromBlock
344  (
345  vtkMultiBlockDataSet* output,
346  const arrayRange&,
347  const label datasetNo
348  );
349 
350  // Convenience method use to convert the readers from VTK 5
351  // multiblock API to the current composite data infrastructure
352  static label GetNumberOfDataSets
353  (
354  vtkMultiBlockDataSet* output,
355  const arrayRange&
356  );
357 
358 
359  // Update information helper functions
360 
361  //- Update the mesh parts selected in the GUI
362  void topoChangePartsStatus();
363 
364  //- Internal mesh info
365  void updateInfoInternalMesh(vtkDataArraySelection*);
366 
367  //- Lagrangian info
368  void updateInfolagrangian(vtkDataArraySelection*);
369 
370  //- Lagrangian info
371  void updateInfoLagrangian(vtkDataArraySelection*);
372 
373  //- Patch info
374  void updateInfoPatches
375  (
376  vtkDataArraySelection*,
377  stringList&,
378  const bool
379  );
380 
381  //- Set info
382  void updateInfoSets(vtkDataArraySelection*);
383 
384  //- Zone info
385  void updateInfoZones(vtkDataArraySelection*);
386 
387  //- Get non-empty zone names for zoneType from file
388  template<class ZonesType>
389  wordList getZoneNames(const word& zonesName) const;
390 
391  //- Get non-empty zone names from mesh info
392  template<class ZonesType>
393  wordList getZoneNames(const ZonesType&) const;
394 
395  //- Add objects of Type to paraview array selection
396  template<class Type>
397  label addToSelection
398  (
399  vtkDataArraySelection*,
400  const IOobjectList&,
401  const string& suffix=string::null
402  );
403 
404  //- Add objects of Type to paraview array selection
405  template<class meshType>
406  void addFieldsToSelection
407  (
408  vtkDataArraySelection *select,
409  const IOobjectList& objects,
410  const string& suffix=string::null
411  );
412 
413  //- Add objects of Type to paraview array selection
414  template<class meshType>
415  void addInternalFieldsToSelection
416  (
417  vtkDataArraySelection *select,
418  const IOobjectList& objects,
419  const string& suffix=string::null
420  );
421 
422  //- Field info
423  void updateInfoFields();
424 
425  //- Lagrangian field info
426  void updateInfolagrangianFields();
427 
428  //- Lagrangian field info
429  void updateInfoLagrangianFields();
430 
431 
432  // Update helper functions
433 
434  //- Clear the reconstructor objects
435  void clearReconstructors();
436 
437  //- Clear the OpenFOAM mesh if "Cache Mesh" is false and/or if we
438  // have switched from reading between the reconstructed and the
439  // decomposed case.
440  void clearFoamMesh();
441 
442  //- OpenFOAM mesh
443  void updateFoamMesh();
444 
445  //- Reset data counters
446  void resetCounters();
447 
448  //- Reduce memory footprint after conversion
449  void reduceMemory();
450 
451  //- Volume fields
452  void updateVolFields(vtkMultiBlockDataSet*);
453 
454  //- Point fields
455  void updatePointFields(vtkMultiBlockDataSet*);
456 
457  //- Lagrangian fields
458  void updatelagrangianFields(vtkMultiBlockDataSet*);
459 
460  //- Lagrangian fields
461  void updateLagrangianFields(vtkMultiBlockDataSet*);
462 
463 
464  // Mesh conversion functions
465 
466  //- Volume mesh
467  void convertMeshVolume(vtkMultiBlockDataSet*, int& blockNo);
468 
469  //- Lagrangian mesh
470  void convertMeshlagrangian(vtkMultiBlockDataSet*, int& blockNo);
471 
472  //- Lagrangian mesh
473  void convertMeshLagrangian(vtkMultiBlockDataSet*, int& blockNo);
474 
475  //- Patches
476  void convertMeshPatches(vtkMultiBlockDataSet*, int& blockNo);
477 
478  //- Cell zones
479  void convertMeshCellZones(vtkMultiBlockDataSet*, int& blockNo);
480 
481  //- Face zones
482  void convertMeshFaceZones(vtkMultiBlockDataSet*, int& blockNo);
483 
484  //- Point zones
485  void convertMeshPointZones(vtkMultiBlockDataSet*, int& blockNo);
486 
487  //- Cell sets
488  void convertMeshCellSets(vtkMultiBlockDataSet*, int& blockNo);
489 
490  //- Face sets
491  void convertMeshFaceSets(vtkMultiBlockDataSet*, int& blockNo);
492 
493  //- Point sets
494  void convertMeshPointSets(vtkMultiBlockDataSet*, int& blockNo);
495 
496 
497  // Add mesh functions
498 
499  //- Add internal mesh/cell set
500  vtkUnstructuredGrid* volumeVTKMesh(const fvMesh&, polyDecomp&);
501 
502  //- Add Lagrangian mesh
503  vtkPolyData* lagrangianVTKMesh
504  (
505  const word&,
507  );
508 
509  //- Add Lagrangian mesh
510  vtkPolyData* LagrangianVTKMesh
511  (
512  const word&,
515  );
516 
517  //- Add patch mesh
518  template<class PatchType>
519  vtkPolyData* patchVTKMesh(const PatchType&);
520 
521  //- Add point zone
522  vtkPolyData* pointZoneVTKMesh(const fvMesh&, const labelList&);
523 
524  //- Add face set
525  vtkPolyData* faceSetVTKMesh(const fvMesh&, const faceSet&);
526 
527  //- Add point mesh
528  vtkPolyData* pointSetVTKMesh(const fvMesh&, const pointSet&);
529 
530 
531  // Field conversion functions
532 
533  //- Convert fields
534  void convertFields(vtkMultiBlockDataSet*);
535 
536  //- Convert Lagrangian fields
537  void convertlagrangianFields(vtkMultiBlockDataSet*);
538 
539  //- Convert Lagrangian fields
540  void convertLagrangianFields(vtkMultiBlockDataSet*);
541 
542 
543  //- Add the fields in the selected time directory to the selection
544  // lists
545  template<class GeoField>
546  label addObjectsToSelection
547  (
548  vtkDataArraySelection*,
549  const IOobjectList&,
550  const string& suffix=string::null
551  );
552 
553 
554  // Convert OpenFOAM fields
555 
556  //- Volume fields - all types
557  template<class Type>
558  void convertVolFields
559  (
561  const IOobjectList&,
562  const bool interpFields,
563  vtkMultiBlockDataSet* output
564  );
565 
566  //- Volume internal fields - all types
567  template<class Type>
568  void convertVolInternalFields
569  (
570  const IOobjectList&,
571  vtkMultiBlockDataSet* output
572  );
573 
574  //- Volume field - all selected parts
575  template<class Type>
576  void convertVolFieldBlock
577  (
578  const VolField<Type>&,
580  vtkMultiBlockDataSet* output,
581  const arrayRange&,
582  const List<polyDecomp>& decompLst
583  );
584 
585  //- Volume internal field - all selected parts
586  template<class Type>
587  void convertVolInternalFieldBlock
588  (
589  const typename VolField<Type>::Internal&,
590  vtkMultiBlockDataSet* output,
591  const arrayRange&,
592  const List<polyDecomp>& decompLst
593  );
594 
595  //- Volume field
596  template<class Type>
597  void convertVolInternalField
598  (
599  const typename VolField<Type>::Internal&,
600  vtkMultiBlockDataSet* output,
601  const arrayRange&,
602  const label datasetNo,
603  const polyDecomp&
604  );
605 
606  //- Patch field
607  template<class Type>
608  void convertPatchField
609  (
610  const word& name,
611  const Field<Type>&,
612  vtkMultiBlockDataSet* output,
613  const arrayRange&,
614  const label datasetNo
615  );
616 
617  //- Face set/zone field from volField
618  template<class Type>
619  void convertSurfaceField
620  (
621  const VolField<Type>&,
622  vtkMultiBlockDataSet* output,
623  const arrayRange&,
624  const label datasetNo,
625  const fvMesh&,
626  const labelList& faceLabels
627  );
628 
629  //- Face set/zone field from surfaceField
630  template<class Type>
631  void convertSurfaceField
632  (
633  const SurfaceField<Type>& tf,
634  vtkMultiBlockDataSet* output,
635  const arrayRange& range,
636  const label datasetNo,
637  const fvMesh& mesh,
638  const labelList& faceLabels
639  );
640 
641  //- Surface fields - all types
642  template<class Type>
643  void convertSurfaceFields
644  (
645  const IOobjectList& objects,
646  vtkMultiBlockDataSet* output
647  );
648 
649  //- Lagrangian fields - all types
650  template<class Type>
651  void convertlagrangianFields
652  (
653  const IOobjectList&,
654  vtkMultiBlockDataSet* output,
655  const label datasetNo
656  );
657 
658  //- Lagrangian fields - all types
659  template<class Type, template<class> class GeoField>
660  void convertLagrangianFields
661  (
662  const IOobjectList&,
663  vtkMultiBlockDataSet* output,
664  const label datasetNo
665  );
666 
667  //- Lagrangian field
668  template<class Type>
669  void convertlagrangianField
670  (
671  const IOField<Type>&,
672  vtkMultiBlockDataSet* output,
673  const arrayRange&,
674  const label datasetNo
675  );
676 
677  //- Lagrangian field
678  template<class Type, template<class> class GeoField>
679  void convertLagrangianField
680  (
681  const GeoField<Type>&,
682  vtkMultiBlockDataSet* output,
683  const arrayRange&,
684  const label datasetNo
685  );
686 
687  //- Point fields - all types
688  template<class Type>
689  void convertPointFields
690  (
691  const IOobjectList&,
692  vtkMultiBlockDataSet* output
693  );
694 
695  //- Point field - all selected parts
696  template<class Type>
697  void convertPointFieldBlock
698  (
699  const PointField<Type>&,
700  vtkMultiBlockDataSet* output,
701  const arrayRange&,
702  const List<polyDecomp>&
703  );
704 
705  //- Point fields
706  template<class Type>
707  void convertPointField
708  (
709  const PointField<Type>&,
710  const VolField<Type>&,
711  vtkMultiBlockDataSet* output,
712  const arrayRange&,
713  const label datasetNo,
714  const polyDecomp&
715  );
716 
717  //- Patch point field
718  template<class Type>
719  void convertPatchPointField
720  (
721  const word& name,
722  const Field<Type>&,
723  vtkMultiBlockDataSet* output,
724  const arrayRange&,
725  const label datasetNo
726  );
727 
728 
729  // GUI selection helper functions
730 
731  //- Only keep what is listed in hashSet
732  static void pruneObjectList
733  (
734  IOobjectList&,
735  const wordHashSet&
736  );
737 
738  //- Get the list of selected objects
739  IOobjectList getObjects
740  (
741  const wordHashSet& selected,
742  const fvMesh& mesh,
743  const fileName& instance,
744  const fileName& local = ""
745  );
746 
747  //- Retrieve the current selections
748  static wordHashSet getSelected(vtkDataArraySelection*);
749 
750  //- Retrieve a sub-list of the current selections
751  static wordHashSet getSelected
752  (
753  vtkDataArraySelection*,
754  const arrayRange&
755  );
756 
757  //- Retrieve the current selections
758  static stringList getSelectedArrayEntries
759  (
760  vtkDataArraySelection*,
761  const bool debug = vtkPVFoam::debug
762  );
763 
764  //- Retrieve a sub-list of the current selections
765  static stringList getSelectedArrayEntries
766  (
767  vtkDataArraySelection*,
768  const arrayRange&,
769  const bool debug = vtkPVFoam::debug
770  );
771 
772  //- Set selection(s)
773  static void setSelectedArrayEntries
774  (
775  vtkDataArraySelection*,
776  const stringList&
777  );
778 
779  //- Get the first word from the mesh parts selection
780  word getPartName(const int);
781 
782 
783 public:
784 
785  //- Static data members
786 
787  ClassName("vtkPVFoam");
788 
789 
790  // Constructors
791 
792  //- Construct from components
793  vtkPVFoam
794  (
795  const char* const FileName,
796  vtkPVFoamReader* reader
797  );
798 
799  //- Disallow default bitwise copy construction
800  vtkPVFoam(const vtkPVFoam&) = delete;
801 
802 
803  //- Destructor
804  ~vtkPVFoam();
805 
806 
807  // Member Functions
808 
809  //- Update information for selection dialogs
810  void updateInfo();
811 
812  //- Update
813  void Update
814  (
815  vtkMultiBlockDataSet* output,
816  vtkMultiBlockDataSet* lagrangianOutput,
817  vtkMultiBlockDataSet* LagrangianOutput
818  );
819 
820  //- Clean any storage
821  void CleanUp();
822 
823  //- Allocate and return a list of selected times
824  // returns the count via the parameter
825  double* findTimes(const bool first, int& nTimeSteps);
826 
827  //- Add/remove patch names to/from the view
828  void renderPatchNames(vtkRenderer*, const bool show);
829 
830  //- Set the runTime to the first plausible request time,
831  // returns the timeIndex
832  // sets to "constant" on error
833  int setTime(int count, const double requestTimes[]);
834 
835  //- The current time index
836  int timeIndex() const
837  {
838  return timeIndex_;
839  }
840 
841 
842  // Access
843 
844  //- Debug information
845  void PrintSelf(ostream&, vtkIndent) const;
846 
847 
848  // Member Operators
849 
850  //- Disallow default bitwise assignment
851  void operator=(const vtkPVFoam&) = delete;
852 };
853 
854 
855 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
856 
857 } // End namespace Foam
858 
859 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
860 
861 #ifdef NoRepository
862  #include "vtkPVFoamTemplates.C"
863 #endif
864 
865 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
866 
867 #endif
868 
869 // ************************************************************************* //
scalar range
label n
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Pre-declare SubField and related Field type.
Definition: Field.H:83
Generic GeometricField class.
A HashTable with keys but without contents.
Definition: HashSet.H:62
A primitive field of type <Type> with automated input and output.
Definition: IOField.H:53
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
Interpolation class within a primitive patch. Allows interpolation from points to faces and vice vers...
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:66
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A list of face labels.
Definition: faceSet.H:51
A class for handling file names.
Definition: fileName.H:82
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
A set of point labels.
Definition: pointSet.H:51
static const string null
An empty string.
Definition: string.H:88
Provides a reader interface for OpenFOAM to VTK interaction.
Definition: vtkPVFoam.H:112
ClassName("vtkPVFoam")
Static data members.
double * findTimes(const bool first, int &nTimeSteps)
Allocate and return a list of selected times.
int setTime(int count, const double requestTimes[])
Set the runTime to the first plausible request time,.
void operator=(const vtkPVFoam &)=delete
Disallow default bitwise assignment.
void PrintSelf(ostream &, vtkIndent) const
Debug information.
int timeIndex() const
The current time index.
Definition: vtkPVFoam.H:835
void Update(vtkMultiBlockDataSet *output, vtkMultiBlockDataSet *lagrangianOutput, vtkMultiBlockDataSet *LagrangianOutput)
Update.
void renderPatchNames(vtkRenderer *, const bool show)
Add/remove patch names to/from the view.
vtkPVFoam(const char *const FileName, vtkPVFoamReader *reader)
Construct from components.
void CleanUp()
Clean any storage.
~vtkPVFoam()
Destructor.
void updateInfo()
Update information for selection dialogs.
A class for handling words, derived from string.
Definition: word.H:62
Macro definitions for declaring ClassName(), NamespaceName(), etc.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const tensorField & tf
tUEqn clear()
Namespace for OpenFOAM.
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
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
void operator+=(fvMatrix< Type > &fvEqn, const CarrierEqn< Type > &cEqn)
Add to a finite-volume equation.
Definition: CarrierEqn.C:71
objects