GEOSX
VTKPolyDataWriterInterface.hpp
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2018-2020 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2020 The Board of Trustees of the Leland Stanford Junior University
7  * Copyright (c) 2018-2020 TotalEnergies
8  * Copyright (c) 2019- GEOSX Contributors
9  * All rights reserved
10  *
11  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
12  * ------------------------------------------------------------------------------------------------------------
13  */
14 
15 #ifndef GEOS_FILEIO_VTK_VTKPOLYDATAWRITERINTERFACE_HPP_
16 #define GEOS_FILEIO_VTK_VTKPOLYDATAWRITERINTERFACE_HPP_
17 
18 #include "common/DataTypes.hpp"
21 #include "fileIO/vtk/VTKPVDWriter.hpp"
22 #include "fileIO/vtk/VTKVTMWriter.hpp"
23 #include "codingUtilities/EnumStrings.hpp"
24 
25 class vtkUnstructuredGrid;
26 class vtkPointData;
27 class vtkCellData;
28 
29 namespace geos
30 {
31 
32 class DomainPartition;
33 class ElementRegionBase;
34 class ParticleRegionBase;
35 class EmbeddedSurfaceNodeManager;
36 class ElementRegionManager;
37 class NodeManager;
38 class ParticleManager;
39 
40 namespace vtk
41 {
42 
43 enum struct VTKOutputMode
44 {
45  BINARY,
46  ASCII
47 };
48 
49 enum struct VTKRegionTypes
50 {
51  CELL,
52  WELL,
53  SURFACE,
54  PARTICLE,
55  ALL
56 };
57 
59 ENUM_STRINGS( VTKOutputMode,
60  "binary",
61  "ascii" );
62 
64 ENUM_STRINGS( VTKRegionTypes,
65  "cell",
66  "well",
67  "surface",
68  "particle",
69  "all" );
70 
75 {
76 public:
81  explicit VTKPolyDataWriterInterface( string outputName );
82 
87  void setWriteGhostCells( bool writeGhostCells )
88  {
89  m_writeGhostCells = writeGhostCells;
90  }
91 
98  void setPlotLevel( integer plotLevel )
99  {
100  m_plotLevel = dataRepository::toPlotLevel( plotLevel );
101  }
102 
107  void setOutputMode( VTKOutputMode mode )
108  {
109  m_outputMode = mode;
110  }
111 
116  void setOutputRegionType( VTKRegionTypes regionType )
117  {
118  m_outputRegionType = regionType;
119  }
120 
126  void setOutputLocation( string outputDir, string outputName )
127  {
128  m_outputDir = std::move( outputDir );
129  m_outputName = std::move( outputName );
130  m_pvd.setFileName( joinPath( m_outputDir, m_outputName ) + ".pvd" );
131  }
132 
137  void setOnlyPlotSpecifiedFieldNamesFlag( integer const onlyPlotSpecifiedFieldNames )
138  {
139  m_onlyPlotSpecifiedFieldNames = onlyPlotSpecifiedFieldNames;
140  }
141 
146  void setFieldNames( arrayView1d< string const > const & fieldNames )
147  {
148  m_fieldNames.insert( fieldNames.begin(), fieldNames.end() );
149  }
150 
155  void setLevelNames( arrayView1d< string const > const & levelNames )
156  {
157  m_levelNames.insert( levelNames.begin(), levelNames.end() );
158  }
159 
193  void write( real64 time, integer cycle, DomainPartition const & domain );
194 
199  void clearData();
200 
201 
202 private:
203 
209  bool isFieldPlotEnabled( dataRepository::WrapperBase const & wrapper ) const;
210 
221  void writeCellElementRegions( real64 time,
222  ElementRegionManager const & elemManager,
223  NodeManager const & nodeManager,
224  string const & path ) const;
225 
226  void writeParticleRegions( real64 const time,
227  ParticleManager const & particleManager,
228  string const & path ) const;
229 
238  void writeWellElementRegions( real64 time,
239  ElementRegionManager const & elemManager,
240  NodeManager const & nodeManager,
241  string const & path ) const;
242 
253  void writeSurfaceElementRegions( real64 time,
254  ElementRegionManager const & elemManager,
255  NodeManager const & nodeManager,
256  EmbeddedSurfaceNodeManager const & embSurfNodeManager,
257  string const & path ) const;
258 
267  void writeVtmFile( integer const cycle,
268  DomainPartition const & domain,
269  VTKVTMWriter const & vtmWriter ) const;
276  void writeNodeFields( NodeManager const & nodeManager,
277  arrayView1d< localIndex const > const & nodeIndices,
278  vtkPointData * pointData ) const;
279 
285  void writeElementFields( ElementRegionBase const & subRegion,
286  vtkCellData * cellData ) const;
287 
288  void writeParticleFields( ParticleRegionBase const & region,
289  vtkCellData * cellData ) const;
290 
299  void writeUnstructuredGrid( string const & path,
300  vtkUnstructuredGrid * ug ) const;
301 
302 private:
303 
305  string m_outputDir;
306 
308  string m_outputName;
309 
312  VTKPVDWriter m_pvd;
313 
315  bool m_writeGhostCells;
316 
318  dataRepository::PlotLevel m_plotLevel;
319 
321  integer m_onlyPlotSpecifiedFieldNames;
322 
324  bool m_requireFieldRegistrationCheck;
325 
327  std::set< string > m_fieldNames;
328 
330  std::set< string > m_levelNames;
331 
333  integer m_previousCycle;
334 
336  VTKOutputMode m_outputMode;
337 
339  VTKRegionTypes m_outputRegionType;
340 };
341 
342 } // namespace vtk
343 } // namespace geos
344 
345 #endif
Partition of the decomposed physical domain. It also manages the connexion information to its neighbo...
The ElementRegionBase is the base class to manage the data stored at the element level.
The ElementRegionManager class provides an interface to ObjectManagerBase in order to manage ElementR...
The EmbeddedSurfaceNodeManager class provides an interface to ObjectManagerBase in order to manage no...
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
Definition: NodeManager.hpp:45
The ParticleManager class provides an interface to ObjectManagerBase in order to manage ParticleRegio...
The ParticleRegionBase is the base class to manage the data stored at the particle level.
Base class for all wrappers containing common operations.
Definition: WrapperBase.hpp:55
VTK PVD Writer class.
void setFileName(string fileName)
Set the output file name.
Encapsulate output methods for vtk.
void setOutputLocation(string outputDir, string outputName)
Set the output directory name.
void setWriteGhostCells(bool writeGhostCells)
Defines if the vtk outputs should contain the ghost cells.
void write(real64 time, integer cycle, DomainPartition const &domain)
Main method of this class. Write all the files for one time step.
void setFieldNames(arrayView1d< string const > const &fieldNames)
Set the names of the fields to output.
void setPlotLevel(integer plotLevel)
Sets the plot level.
void setOnlyPlotSpecifiedFieldNamesFlag(integer const onlyPlotSpecifiedFieldNames)
Set the flag to decide whether we only plot the fields specified by fieldNames, or if we also plot fi...
void setOutputRegionType(VTKRegionTypes regionType)
Set the output region type.
VTKPolyDataWriterInterface(string outputName)
Constructor.
void setOutputMode(VTKOutputMode mode)
Set the binary mode.
void clearData()
Clears the datasets accumulated in the pvd writer.
void setLevelNames(arrayView1d< string const > const &levelNames)
Set the names of the mesh levels to output.
VTM Writer class.
PlotLevel toPlotLevel(int const val)
Function to get a PlotLevel enum from an int.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:220
ENUM_STRINGS(LinearSolverParameters::SolverType, "direct", "cg", "gmres", "fgmres", "bicgstab", "preconditioner")
Declare strings associated with enumeration values.
double real64
64-bit floating point type.
Definition: DataTypes.hpp:139
std::int32_t integer
Signed integer type.
Definition: DataTypes.hpp:122
std::string joinPath(ARGS const &... args)
Join parts of a path.
Definition: Path.hpp:188