1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
29 Provides a reader interface for OpenFOAM to VTK interaction.
37 vtkPV3FoamMeshLagrangian.C
40 vtkPV3FoamMeshVolume.C
43 vtkPV3FoamLagrangianFields.H
44 vtkPV3FoamPatchField.H
45 vtkPV3FoamPointFields.H
47 vtkPV3FoamUpdateInfo.C
48 vtkPV3FoamUpdateInfoFields.H
51 vtkPV3FoamAddToSelection.H
54 vtkDataArrayTemplateImplicit.txx
56 \*---------------------------------------------------------------------------*/
61 // do not include legacy strstream headers
62 #ifndef VTK_EXCLUDE_STRSTREAM_HEADERS
63 # define VTK_EXCLUDE_STRSTREAM_HEADERS
66 #include "className.H"
68 #include "stringList.H"
70 #include "primitivePatch.H"
71 #include "PrimitivePatchInterpolation.H"
72 #include "volPointInterpolation.H"
74 #undef VTKPV3FOAM_DUALPORT
76 // * * * * * * * * * * * * * Forward Declarations * * * * * * * * * * * * * //
78 class vtkDataArraySelection;
81 class vtkPV3FoamReader;
84 class vtkMultiBlockDataSet;
86 class vtkUnstructuredGrid;
89 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
94 // Foam class forward declarations
103 template<class Type> class IOField;
104 template<class Type> class List;
106 /*---------------------------------------------------------------------------*\
107 Class vtkPV3Foam Declaration
108 \*---------------------------------------------------------------------------*/
114 //- Bookkeeping for GUI checklists and the multi-block organization
124 partInfo(const char *name, const int blockNo=0)
132 //- Return the block holding these datasets
138 //- Assign block number, return previous value
139 int block(int blockNo)
146 const char* name() const
158 return start_ + size_;
177 //- Assign new start and reset the size
178 void operator=(const int i)
184 //- Increment the size
185 void operator+=(const int n)
191 //- bookkeeping for polyhedral cell decomposition
192 // hide in extra pointMap (cellSet/cellZone) for now
195 labelList superCells_;
196 labelList addPointCellLabels_;
204 //- Label of original cell for decomposed cells
205 labelList& superCells()
210 //- Label of original cell for decomposed cells
211 const labelList& superCells() const
216 //- Cell-centre labels for additional points of decomposed cells
217 labelList& addPointCellLabels()
219 return addPointCellLabels_;
222 //- Cell-centre labels for additional points of decomposed cells
223 const labelList& addPointCellLabels() const
225 return addPointCellLabels_;
228 //- Point labels for subsetted meshes
229 labelList& pointMap()
234 //- Point labels for subsetted meshes
235 const labelList& pointMap() const
245 addPointCellLabels_.clear();
253 //- Access to the controlling vtkPV3FoamReader
254 vtkPV3FoamReader* reader_;
256 //- Foam time control
257 autoPtr<Time> dbPtr_;
265 //- The mesh directory for the region
271 //- Track changes in mesh geometry
274 //- Track changes in fields
277 //- Selected geometrical parts (internalMesh, patches, ...)
278 boolList partStatus_;
280 //- Datasets corresponding to selected geometrical pieces
281 // a negative number indicates that no vtkmesh exists for this piece
282 labelList partDataset_;
284 //- First instance and size of various mesh parts
285 // used to index into partStatus_ and partDataset_
286 partInfo partInfoVolume_;
287 partInfo partInfoPatches_;
288 partInfo partInfoLagrangian_;
289 partInfo partInfoCellZones_;
290 partInfo partInfoFaceZones_;
291 partInfo partInfoPointZones_;
292 partInfo partInfoCellSets_;
293 partInfo partInfoFaceSets_;
294 partInfo partInfoPointSets_;
296 //- Decomposed cells information (mesh regions)
298 List<polyDecomp> regionPolyDecomp_;
300 //- Decomposed cells information (cellZone meshes)
301 List<polyDecomp> zonePolyDecomp_;
303 //- Decomposed cells information (cellSet meshes)
304 List<polyDecomp> csetPolyDecomp_;
306 //- List of patch names for rendering to window
307 List<vtkTextActor*> patchTextActorsPtrs_;
309 // Private Member Functions
311 // Convenience method use to convert the readers from VTK 5
312 // multiblock API to the current composite data infrastructure
313 static void AddToBlock
315 vtkMultiBlockDataSet* output,
318 const label datasetNo,
319 const string& datasetName
322 // Convenience method use to convert the readers from VTK 5
323 // multiblock API to the current composite data infrastructure
324 static vtkDataSet* GetDataSetFromBlock
326 vtkMultiBlockDataSet* output,
328 const label datasetNo
331 // Convenience method use to convert the readers from VTK 5
332 // multiblock API to the current composite data infrastructure
333 static label GetNumberOfDataSets
335 vtkMultiBlockDataSet* output,
339 //- Reset data counters
340 void resetCounters();
342 // Update information helper functions
344 //- Update the mesh parts selected in the GUI
345 void updateMeshPartsStatus();
347 //- Internal mesh info
348 void updateInfoInternalMesh();
351 void updateInfoLagrangian();
354 void updateInfoPatches();
357 void updateInfoSets();
360 void updateInfoZones();
362 //- Read zone names for zoneType from file
363 wordList readZoneNames(const word& zoneType);
365 //- Add objects of Type to paraview array selection
369 vtkDataArraySelection*,
371 const string& suffix=string::null
375 template<template<class> class patchType, class meshType>
376 void updateInfoFields(vtkDataArraySelection*);
378 //- Lagrangian field info
379 void updateInfoLagrangianFields();
382 // Update helper functions
385 void updateFoamMesh();
387 //- Reduce memory footprint after conversion
391 void updateVolFields(vtkMultiBlockDataSet*);
394 void updatePointFields(vtkMultiBlockDataSet*);
396 //- Lagrangian fields
397 void updateLagrangianFields(vtkMultiBlockDataSet*);
400 // Mesh conversion functions
403 void convertMeshVolume(vtkMultiBlockDataSet*, int& blockNo);
406 void convertMeshLagrangian(vtkMultiBlockDataSet*, int& blockNo);
409 void convertMeshPatches(vtkMultiBlockDataSet*, int& blockNo);
412 void convertMeshCellZones(vtkMultiBlockDataSet*, int& blockNo);
415 void convertMeshFaceZones(vtkMultiBlockDataSet*, int& blockNo);
417 //- Point zone meshes
418 void convertMeshPointZones(vtkMultiBlockDataSet*, int& blockNo);
421 void convertMeshCellSets(vtkMultiBlockDataSet*, int& blockNo);
424 void convertMeshFaceSets(vtkMultiBlockDataSet*, int& blockNo);
427 void convertMeshPointSets(vtkMultiBlockDataSet*, int& blockNo);
430 // Add mesh functions
432 //- Add internal mesh/cell set meshes
433 vtkUnstructuredGrid* volumeVTKMesh(const fvMesh&, polyDecomp&);
435 //- Add Lagrangian mesh
436 vtkPolyData* lagrangianVTKMesh
439 const word& cloudName
443 vtkPolyData* patchVTKMesh(const polyPatch&);
445 //- Add face zone mesh
446 vtkPolyData* faceZoneVTKMesh
449 const labelList& faceLabels
453 vtkPolyData* pointZoneVTKMesh
456 const labelList& pointLabels
459 //- Add face set mesh
460 vtkPolyData* faceSetVTKMesh
467 vtkPolyData* pointSetVTKMesh
473 // Field conversion functions
475 //- Convert volume fields
476 void convertVolFields(vtkMultiBlockDataSet*);
478 //- Convert point fields
479 void convertPointFields(vtkMultiBlockDataSet*);
481 //- Convert Lagrangian fields
482 void convertLagrangianFields(vtkMultiBlockDataSet*);
485 //- Add the fields in the selected time directory to the selection
487 template<class GeoField>
488 label addObjectsToSelection
490 vtkDataArraySelection*,
492 const string& suffix=string::null
496 // Convert Foam fields
498 //- Volume fields - all types
500 void convertVolFields
503 const PtrList<PrimitivePatchInterpolation<primitivePatch> >&,
505 vtkMultiBlockDataSet* output
508 //- Volume field - all selected parts
510 void convertVolFieldBlock
512 const GeometricField<Type, fvPatchField, volMesh>&,
513 autoPtr<GeometricField<Type, pointPatchField, pointMesh> >&,
514 vtkMultiBlockDataSet* output,
515 const partInfo& selector,
516 const List<polyDecomp>& decompLst
523 const GeometricField<Type, fvPatchField, volMesh>&,
524 vtkMultiBlockDataSet* output,
526 const label datasetNo,
532 void convertPatchField
536 vtkMultiBlockDataSet* output,
538 const label datasetNo
541 //- face set/zone field
543 void convertFaceField
545 const GeometricField<Type, fvPatchField, volMesh>&,
546 vtkMultiBlockDataSet* output,
548 const label datasetNo,
550 const labelList& faceLabels
553 //- face set/zone field
555 void convertFaceField
557 const GeometricField<Type, fvPatchField, volMesh>&,
558 vtkMultiBlockDataSet* output,
560 const label datasetNo,
565 //- Lagrangian fields - all types
567 void convertLagrangianFields
570 vtkMultiBlockDataSet* output,
571 const label datasetNo
576 void convertLagrangianField
578 const IOField<Type>&,
579 vtkMultiBlockDataSet* output,
581 const label datasetNo
584 //- Point fields - all types
586 void convertPointFields
591 vtkMultiBlockDataSet* output
594 //- Point field - all selected parts
596 void convertPointFieldBlock
598 const GeometricField<Type, pointPatchField, pointMesh>&,
599 vtkMultiBlockDataSet* output,
600 const partInfo& selector,
601 const List<polyDecomp>&
606 void convertPointField
608 const GeometricField<Type, pointPatchField, pointMesh>&,
609 const GeometricField<Type, fvPatchField, volMesh>&,
610 vtkMultiBlockDataSet* output,
612 const label datasetNo,
616 //- Patch point field
618 void convertPatchPointField
622 vtkMultiBlockDataSet* output,
624 const label datasetNo
628 // GUI selection helper functions
630 //- Extract up to the first non-word characters
631 inline static word getFirstWord(const char*);
633 //- Only keep what is listed in hashSet
634 static void pruneObjectList
640 //- Retrieve the current selections
641 static wordHashSet getSelected(vtkDataArraySelection*);
643 //- Retrieve a sub-list of the current selections
644 static wordHashSet getSelected
646 vtkDataArraySelection*,
650 //- Retrieve the current selections
651 static stringList getSelectedArrayEntries(vtkDataArraySelection*);
653 //- Retrieve a sub-list of the current selections
654 static stringList getSelectedArrayEntries
656 vtkDataArraySelection*,
661 static void setSelectedArrayEntries
663 vtkDataArraySelection*,
667 //- Get the first word from the mesh parts selection
668 word getPartName(int);
671 //- Disallow default bitwise copy construct
672 vtkPV3Foam(const vtkPV3Foam&);
674 //- Disallow default bitwise assignment
675 void operator=(const vtkPV3Foam&);
680 //- Static data members
682 ClassName("vtkPV3Foam");
687 //- Construct from components
690 const char* const FileName,
691 vtkPV3FoamReader* reader
707 vtkMultiBlockDataSet* output,
708 vtkMultiBlockDataSet* lagrangianOutput
711 //- Clean any storage
714 //- Allocate and return a list of selected times
715 // returns the count via the parameter
716 double* findTimes(int& nTimeSteps);
718 //- Add patch names to the display
719 void addPatchNames(vtkRenderer* renderer);
721 //- Remove patch names from the display
722 void removePatchNames(vtkRenderer* renderer);
724 //- set the runTime to the first plausible request time,
725 // returns the timeIndex
726 // sets to "constant" on error
727 int setTime(int count, const double requestTimes[]);
730 //- The current time index
731 int timeIndex() const
739 //- Debug information
740 void PrintSelf(ostream&, vtkIndent) const;
742 //- Simple memory used debugging information
743 static void printMemory();
748 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
750 } // End namespace Foam
752 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
754 # include "vtkPV3FoamI.H"
758 // ************************************************************************* //