GEOS
SiloFile.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2016-2024 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2024 Total, S.A
7  * Copyright (c) 2018-2024 The Board of Trustees of the Leland Stanford Junior University
8  * Copyright (c) 2023-2024 Chevron
9  * Copyright (c) 2019- GEOS/GEOSX Contributors
10  * All rights reserved
11  *
12  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
13  * ------------------------------------------------------------------------------------------------------------
14  */
15 
20 #ifndef GEOS_FILEIO_SILO_SILOFILE_HPP_
21 #define GEOS_FILEIO_SILO_SILOFILE_HPP_
22 
23 #include "common/DataTypes.hpp"
24 #include "constitutive/solid/SolidBase.hpp"
25 #include "common/MpiWrapper.hpp"
27 #include "mesh/CellElementSubRegion.hpp"
30 
31 #include "silo.h"
32 
34 struct _PMPIO_baton_t;
35 
37 typedef _PMPIO_baton_t PMPIO_baton_t;
38 
39 namespace geos
40 {
41 
42 class DomainPartition;
43 class MeshLevel;
44 
45 namespace constitutive
46 {
47 class ConstitutiveManager;
48 }
49 
50 // *********************************************************************************************************************
51 // *********************************************************************************************************************
56 class SiloFile
57 {
58 
59 public:
60 
63 
65  virtual ~SiloFile();
66 
71 
76  void initialize( int const numGroups=1 );
77 
81  void finish();
82 
88  int groupRank( int const i ) const;
89 
100  void waitForBatonWrite( int const domainNumber,
101  int const cycleNum,
102  integer const eventCounter,
103  bool const isRestart );
104 
110  void waitForBaton( int const domainNumber, string const & restartFileName );
111 
115  void handOffBaton();
116 
117 
123  void makeSubDirectory( string const & subdir, string const & rootdir )
124  {
125  int const rank = MpiWrapper::commRank( MPI_COMM_GEOS );
126 
127  // char dirname[100];
128  if( rank == 0 )
129  {
130 // DBGetDir(m_dbBaseFilePtr, dirname );
131  DBMkDir( m_dbBaseFilePtr, rootdir.c_str());
132  }
133 
134 // DBGetDir (m_dbFilePtr, dirname );
135  DBMkDir( m_dbFilePtr, subdir.c_str());
136  }
137 
159  void writeMeshObject( string const & meshName,
160  const localIndex nnodes,
161  real64 * coords[3],
162  const globalIndex * globalNodeNum,
163  char const * const ghostNodeName,
164  char const * const ghostZoneName,
165  int const numShapes,
166  int const * shapecnt,
167  const localIndex * const * const meshConnectivity,
168  const globalIndex * const * const globalElementNum,
169  int const * const shapetype,
170  int const * const shapesize,
171  int const cycleNumber,
172  real64 const problemTime );
173 
193  void writePolygonMeshObject( const string & meshName,
194  const localIndex nnodes,
195  real64 * coords[3],
196  const globalIndex * dummy1,
197  const int numRegions,
198  const int * shapecnt,
199  const localIndex * const * const meshConnectivity,
200  const globalIndex * const * const globalElementNum,
201  const int * const * const dummy2,
202  const int * const shapetype,
203  const int * const shapesize,
204  const int cycleNumber,
205  const real64 problemTime,
206  const int lnodelist );
207 
215  void writeDomainPartition( DomainPartition const & domain,
216  int const cycleNum,
217  real64 const problemTime,
218  bool const isRestart );
219 
234  void writeElementMesh( ElementRegionBase const & elementRegion,
235  NodeManager const & nodeManager,
236  string const & meshName,
237  const localIndex nnodes,
238  real64 * coords[3],
239  globalIndex const * const globalNodeNum,
240  char const * const ghostNodeFlag,
241  int const cycleNumber,
242  real64 const problemTime,
243  bool & writeArbitraryPolygon );
244 
252  void writeMeshLevel( MeshLevel const & meshLevel,
253  int const cycleNum,
254  real64 const problemTime,
255  bool const isRestart );
256 
264  void writePointMesh( string const & meshName,
265  const localIndex numPoints,
266  real64 * coords[3],
267  int const cycleNumber,
268  real64 const problemTime );
269 
280  void writeBeamMesh( string const & meshName,
281  const localIndex nnodes,
282  real64 * coords[3],
283  const localIndex_array & node1,
284  const localIndex_array & node2,
285  int const cycleNumber,
286  real64 const problemTime );
287 
297  void writeBeamMesh( string const & meshName,
298  const localIndex nnodes,
299  real64 * coords[3],
300  integer_array & nodelist,
301  int const cycleNumber,
302  real64 const problemTime );
303 
311  void writeMaterialMapsFullStorage( ElementRegionBase const & elementRegion,
312  string const & meshName,
313  string_array const & regionMaterialList,
314  int const cycleNumber,
315  real64 const problemTime );
328  string const & siloDirName,
329  string const & meshname,
330  int const centering,
331  int const cycleNum,
332  real64 const problemTime,
333  bool const isRestart,
334  const localIndex_array & mask );
335 
344  void writeElementRegionSilo( ElementRegionBase const & elemRegion,
345  string const & siloDirName,
346  string const & meshName,
347  int const cycleNum,
348  real64 const problemTime,
349  bool const isRestart );
362  template< typename OUTPUTTYPE >
363  void writeWrappersToSilo( string const & meshname,
364  const dataRepository::Group::wrapperMap & wrappers,
365  int const centering,
366  int const cycleNum,
367  real64 const problemTime,
368  bool const isRestart,
369  string const & multiRoot,
370  const localIndex_array & mask );
371 
382  template< typename OUTTYPE, typename TYPE >
383  void writeDataField( string const & meshName,
384  string const & fieldName,
385  arrayView1d< TYPE const > const & field,
386  int const centering,
387  int const cycleNumber,
388  real64 const problemTime,
389  string const & multiRoot );
390 
401  template< typename OUTTYPE, typename TYPE, int USD >
402  void writeDataField( string const & meshName,
403  string const & fieldName,
404  arrayView2d< TYPE const, USD > const & field,
405  int const centering,
406  int const cycleNumber,
407  real64 const problemTime,
408  string const & multiRoot );
409 
419  template< typename OUTTYPE, typename TYPE, int USD >
420  void writeDataField( string const & meshName,
421  string const & fieldName,
422  arrayView3d< TYPE const, USD > const & field,
423  int const centering,
424  int const cycleNumber,
425  real64 const problemTime,
426  string const & multiRoot );
427 
439  template< typename OUTTYPE, typename TYPE, int NDIM, int USD >
440  void writeDataField( string const & meshName,
441  string const & fieldName,
443  int const siloTensorRank,
444  int const centering,
445  int const cycleNumber,
446  real64 const problemTime,
447  string const & multiRoot );
448 
461  template< typename OUTTYPE, typename TYPE >
462  void writeMaterialDataField( string const & meshName,
463  string const & fieldName,
464  array1d< array1d< arrayView2d< TYPE const > > > const & field,
465  ElementRegionBase const & elemRegion,
466  int const centering,
467  int const cycleNumber,
468  real64 const problemTime,
469  string const & multiRoot,
470  string_array const & materialNames );
471 
483  template< typename OUTTYPE, typename TYPE >
484  void writeMaterialDataField2d( string const & meshName,
485  string const & fieldName,
486  ElementRegionBase const & elemRegion,
487  int const centering,
488  int const cycleNumber,
489  real64 const problemTime,
490  string const & multiRoot,
491  string_array const & materialNames );
492 
504  template< typename OUTTYPE, typename TYPE >
505  void writeMaterialDataField3d( string const & meshName,
506  string const & fieldName,
507  ElementRegionBase const & elemRegion,
508  int const centering,
509  int const cycleNumber,
510  real64 const problemTime,
511  string const & multiRoot,
512  string_array const & materialNames );
513 
525  template< typename OUTTYPE, typename TYPE >
526  void writeMaterialDataField4d( string const & meshName,
527  string const & fieldName,
528  ElementRegionBase const & elemRegion,
529  int const centering,
530  int const cycleNumber,
531  real64 const problemTime,
532  string const & multiRoot,
533  string_array const & materialNames );
534 
542  void writeMaterialVarDefinition( string const & subDir,
543  string const & matDir,
544  localIndex const matIndex,
545  string const & fieldName );
546 
551  void writeStressVarDefinition( string const & MatDir );
552 
558  void writeVectorVarDefinition( string const & fieldName,
559  string const & subDirectory );
560 
566  int getMeshType( string const & meshName ) const;
567 
578  template< typename CBF >
579  void writeMultiXXXX( const DBObjectType type, CBF DBPutMultiCB,
580  int const centering, string const name, int const cycleNumber,
581  string const & multiRoot, const DBoptlist * optlist = nullptr );
582 
583 
592  void clearEmptiesFromMultiObjects( int const cycleNum );
593 
598  void setNumGroups( int const numGroups )
599  {
600  m_numGroups = numGroups;
601  }
602 
607  void setPlotLevel( int const plotLevel )
608  {
609  m_plotLevel = dataRepository::toPlotLevel( plotLevel );
610  }
611 
616  void setWriteEdgeMesh( int const val )
617  {
618  m_writeEdgeMesh = val;
619  }
620 
625  void setWriteFaceMesh( int const val )
626  {
627  m_writeFaceMesh = val;
628  }
629 
634  void setWriteCellElementMesh( int const val )
635  {
636  m_writeCellElementMesh = val;
637  }
638 
643  void setWriteFaceElementMesh( int const val )
644  {
645  m_writeFaceElementMesh = val;
646  }
647 
652  void setPlotFileRoot( string const & fileRoot )
653  {
654  m_plotFileRoot = fileRoot;
655  }
656 
661  void setOutputDirectory( string const & path )
662  {
663  m_siloDirectory = joinPath( path, "siloFiles" );
664  }
665 
670  void setOnlyPlotSpecifiedFieldNamesFlag( integer const onlyPlotSpecifiedFieldNames )
671  {
672  m_onlyPlotSpecifiedFieldNames = onlyPlotSpecifiedFieldNames;
673  }
674 
679  void setFieldNames( arrayView1d< string const > const & fieldNames )
680  {
681  m_fieldNames.insert( fieldNames.begin(), fieldNames.end() );
682  }
683 
684 private:
685 
691  bool isFieldPlotEnabled( dataRepository::WrapperBase const & wrapper ) const;
692 
694  DBfile * m_dbFilePtr;
695 
697  DBfile * m_dbBaseFilePtr;
698 
700  int m_numGroups;
701 
703  PMPIO_baton_t *m_baton;
704 
706  int const m_driver;
707 
709  string m_plotFileRoot;
710 
711  string m_restartFileRoot;
712 
713  string m_siloDirectory;
714 
715  string const m_siloDataSubDirectory = "data";
716 
717  string m_fileName;
718 
719  string m_baseFileName;
720 
721  std::vector< string > m_emptyMeshes;
722 // string_array m_emptyMaterials;
723  std::vector< string > m_emptyVariables;
724 
725  integer m_writeEdgeMesh;
726  integer m_writeFaceMesh;
727  integer m_writeCellElementMesh;
728  integer m_writeFaceElementMesh;
729 
730  dataRepository::PlotLevel m_plotLevel;
731 
733  integer m_onlyPlotSpecifiedFieldNames;
734 
736  bool m_requireFieldRegistrationCheck;
737 
739  std::set< string > m_fieldNames;
740 
741 
742  bool m_ghostFlags;
743 };
744 
748 namespace siloFileUtilities
749 {
758 template< typename OUTTYPE >
759 int DB_TYPE();
760 
768 template< typename TYPE >
770 
775 template< typename TYPE >
777 
785 template< typename OUTTYPE, typename TYPE >
786 OUTTYPE CastField( const TYPE & field, int const i = 0 ); // avoids compiler
787  // warning
794 template<> inline real64 CastField( R1Tensor const & field, int const i )
795 {
796  return field[i];
797 }
798 
806 template<> inline int CastField< int, int >( const int & field, int const dummy )
807 {
808  GEOS_UNUSED_VAR( dummy );
809  return field;
810 }
811 
819 template<> inline long int CastField< long int, long int >( const long int & field, int const dummy )
820 {
821  GEOS_UNUSED_VAR( dummy );
822  return field;
823 }
824 
832 template<> inline int CastField< int, long int >( const long int & field, int const dummy )
833 {
834  GEOS_UNUSED_VAR( dummy );
835  return LvArray::integerConversion< int >( field );
836 }
837 
845 template<> inline long long int CastField< long long int, long long int >( const long long int & field, int const dummy )
846 {
847  GEOS_UNUSED_VAR( dummy );
848  return field;
849 }
850 
858 template<> inline int CastField< int, long long int >( const long long int & field, int const dummy )
859 {
860  GEOS_UNUSED_VAR( dummy );
861  return LvArray::integerConversion< int >( field );
862 }
863 
871 template<> inline real64 CastField< real64, real64 >( const real64 & field, int const dummy )
872 {
873  GEOS_UNUSED_VAR( dummy );
874  return field;
875 }
876 
884 template<> inline float CastField< float, real64 >( const real64 & field, int const dummy )
885 {
886  GEOS_UNUSED_VAR( dummy );
887  return static_cast< float >(field);
888 }
889 
900 template< typename TYPE >
901 void SetVariableNames( string const & fieldName, string_array & varnamestring, char const * varnames[] );
902 
903 
904 }
905 
906 
907 
908 }
909 #endif /* GEOS_FILEIO_SILO_SILOFILE_HPP_ */
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:84
_PMPIO_baton_t PMPIO_baton_t
Type alias for _PMPIO_baton_t struct.
Definition: SiloFile.hpp:34
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.
Class facilitating the representation of a multi-level discretization of a MeshBody.
Definition: MeshLevel.hpp:42
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
Definition: NodeManager.hpp:46
void clearEmptiesFromMultiObjects(int const cycleNum)
void writeMaterialDataField4d(string const &meshName, string const &fieldName, ElementRegionBase const &elemRegion, int const centering, int const cycleNumber, real64 const problemTime, string const &multiRoot, string_array const &materialNames)
void setNumGroups(int const numGroups)
Sets the number of individual Silo files to generate.
Definition: SiloFile.hpp:598
void writeDomainPartition(DomainPartition const &domain, int const cycleNum, real64 const problemTime, bool const isRestart)
write a domain parititon out to silo file
void writeDataField(string const &meshName, string const &fieldName, arrayView2d< TYPE const, USD > const &field, int const centering, int const cycleNumber, real64 const problemTime, string const &multiRoot)
void writeDataField(string const &meshName, string const &fieldName, ArrayView< TYPE const, NDIM, USD > const &field, int const siloTensorRank, int const centering, int const cycleNumber, real64 const problemTime, string const &multiRoot)
void writeWrappersToSilo(string const &meshname, const dataRepository::Group::wrapperMap &wrappers, int const centering, int const cycleNum, real64 const problemTime, bool const isRestart, string const &multiRoot, const localIndex_array &mask)
void makeSiloDirectories()
function to setup directories where silo files will be written
void writeMaterialVarDefinition(string const &subDir, string const &matDir, localIndex const matIndex, string const &fieldName)
void writeStressVarDefinition(string const &MatDir)
void writeMaterialDataField(string const &meshName, string const &fieldName, array1d< array1d< arrayView2d< TYPE const > > > const &field, ElementRegionBase const &elemRegion, int const centering, int const cycleNumber, real64 const problemTime, string const &multiRoot, string_array const &materialNames)
int getMeshType(string const &meshName) const
void setOutputDirectory(string const &path)
Sets the top-level output directory, under which Silo's output dir is nested.
Definition: SiloFile.hpp:661
void writeBeamMesh(string const &meshName, const localIndex nnodes, real64 *coords[3], const localIndex_array &node1, const localIndex_array &node2, int const cycleNumber, real64 const problemTime)
write a beam mesh to silo file
void writeDataField(string const &meshName, string const &fieldName, arrayView3d< TYPE const, USD > const &field, int const centering, int const cycleNumber, real64 const problemTime, string const &multiRoot)
void setPlotFileRoot(string const &fileRoot)
Sets root of the filename that will be read/written.
Definition: SiloFile.hpp:652
void waitForBatonWrite(int const domainNumber, int const cycleNum, integer const eventCounter, bool const isRestart)
Wait for the Baton when writing using PMPIO.
SiloFile()
Default Constructor.
void finish()
finishes/closes up the silo interface
void setWriteCellElementMesh(int const val)
Sets the cell element mesh output option.
Definition: SiloFile.hpp:634
void waitForBaton(int const domainNumber, string const &restartFileName)
Wait for the Baton when reading using PMPIO.
void setPlotLevel(int const plotLevel)
Sets the plot level option.
Definition: SiloFile.hpp:607
void writeMeshObject(string const &meshName, const localIndex nnodes, real64 *coords[3], const globalIndex *globalNodeNum, char const *const ghostNodeName, char const *const ghostZoneName, int const numShapes, int const *shapecnt, const localIndex *const *const meshConnectivity, const globalIndex *const *const globalElementNum, int const *const shapetype, int const *const shapesize, int const cycleNumber, real64 const problemTime)
Write out a single silo mesh object.
void writeMeshLevel(MeshLevel const &meshLevel, int const cycleNum, real64 const problemTime, bool const isRestart)
write a mesh level out to the silo file
void setWriteEdgeMesh(int const val)
Sets the edge mesh output option.
Definition: SiloFile.hpp:616
void makeSubDirectory(string const &subdir, string const &rootdir)
Make a subdirectory within the silo file.
Definition: SiloFile.hpp:123
int groupRank(int const i) const
obtain the group number of the calling processor, indexed from zero.
void setWriteFaceElementMesh(int const val)
Sets the face element mesh output option.
Definition: SiloFile.hpp:643
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...
Definition: SiloFile.hpp:670
void writePolygonMeshObject(const string &meshName, const localIndex nnodes, real64 *coords[3], const globalIndex *dummy1, const int numRegions, const int *shapecnt, const localIndex *const *const meshConnectivity, const globalIndex *const *const globalElementNum, const int *const *const dummy2, const int *const shapetype, const int *const shapesize, const int cycleNumber, const real64 problemTime, const int lnodelist)
void writeGroupSilo(dataRepository::Group const &group, string const &siloDirName, string const &meshname, int const centering, int const cycleNum, real64 const problemTime, bool const isRestart, const localIndex_array &mask)
void writeBeamMesh(string const &meshName, const localIndex nnodes, real64 *coords[3], integer_array &nodelist, int const cycleNumber, real64 const problemTime)
void setFieldNames(arrayView1d< string const > const &fieldNames)
Set the names of the fields to output.
Definition: SiloFile.hpp:679
void setWriteFaceMesh(int const val)
Sets the face mesh output option.
Definition: SiloFile.hpp:625
void writePointMesh(string const &meshName, const localIndex numPoints, real64 *coords[3], int const cycleNumber, real64 const problemTime)
void initialize(int const numGroups=1)
Initializes silo for input/output.
void writeElementRegionSilo(ElementRegionBase const &elemRegion, string const &siloDirName, string const &meshName, int const cycleNum, real64 const problemTime, bool const isRestart)
void writeElementMesh(ElementRegionBase const &elementRegion, NodeManager const &nodeManager, string const &meshName, const localIndex nnodes, real64 *coords[3], globalIndex const *const globalNodeNum, char const *const ghostNodeFlag, int const cycleNumber, real64 const problemTime, bool &writeArbitraryPolygon)
void writeMultiXXXX(const DBObjectType type, CBF DBPutMultiCB, int const centering, string const name, int const cycleNumber, string const &multiRoot, const DBoptlist *optlist=nullptr)
void handOffBaton()
Hand off the Baton when done writing to file.
void writeVectorVarDefinition(string const &fieldName, string const &subDirectory)
void writeMaterialDataField2d(string const &meshName, string const &fieldName, ElementRegionBase const &elemRegion, int const centering, int const cycleNumber, real64 const problemTime, string const &multiRoot, string_array const &materialNames)
void writeDataField(string const &meshName, string const &fieldName, arrayView1d< TYPE const > const &field, int const centering, int const cycleNumber, real64 const problemTime, string const &multiRoot)
void writeMaterialMapsFullStorage(ElementRegionBase const &elementRegion, string const &meshName, string_array const &regionMaterialList, int const cycleNumber, real64 const problemTime)
void writeMaterialDataField3d(string const &meshName, string const &fieldName, ElementRegionBase const &elemRegion, int const centering, int const cycleNumber, real64 const problemTime, string const &multiRoot, string_array const &materialNames)
virtual ~SiloFile()
Destructor.
Base class for all wrappers containing common operations.
Definition: WrapperBase.hpp:56
PlotLevel toPlotLevel(int const val)
Function to get a PlotLevel enum from an int.
long long int CastField< long long int, long long int >(const long long int &field, int const dummy)
Definition: SiloFile.hpp:845
float CastField< float, real64 >(const real64 &field, int const dummy)
Definition: SiloFile.hpp:884
int CastField< int, long int >(const long int &field, int const dummy)
Definition: SiloFile.hpp:832
int CastField< int, long long int >(const long long int &field, int const dummy)
Definition: SiloFile.hpp:858
real64 CastField< real64, real64 >(const real64 &field, int const dummy)
Definition: SiloFile.hpp:871
OUTTYPE CastField(const TYPE &field, int const i=0)
void SetVariableNames(string const &fieldName, string_array &varnamestring, char const *varnames[])
int CastField< int, int >(const int &field, int const dummy)
Definition: SiloFile.hpp:806
long int CastField< long int, long int >(const long int &field, int const dummy)
Definition: SiloFile.hpp:819
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:180
array1d< localIndex > localIndex_array
A 1-dimensional array of geos::localIndex types.
Definition: DataTypes.hpp:398
int MPI_COMM_GEOS
Global MPI communicator used by GEOSX.
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:88
array1d< string > string_array
A 1-dimensional array of geos::string types.
Definition: DataTypes.hpp:392
array1d< integer > integer_array
A 1-dimensional array of geos::integer types.
Definition: DataTypes.hpp:383
double real64
64-bit floating point type.
Definition: DataTypes.hpp:99
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:85
std::int32_t integer
Signed integer type.
Definition: DataTypes.hpp:82
std::string joinPath(ARGS const &... args)
Join parts of a path.
Definition: Path.hpp:189
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Definition: DataTypes.hpp:196
Array< T, 1 > array1d
Alias for 1D array.
Definition: DataTypes.hpp:176
LvArray::ArrayView< T, NDIM, USD, localIndex, LvArray::ChaiBuffer > ArrayView
Multidimensional array view type. See LvArray:ArrayView for details.
Definition: DataTypes.hpp:148
ArrayView< T, 3, USD > arrayView3d
Alias for 3D array view.
Definition: DataTypes.hpp:212