GEOSX
SiloFile.hpp
Go to the documentation of this file.
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 Total, S.A
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 
19 #ifndef GEOSX_FILEIO_SILO_SILOFILE_HPP_
20 #define GEOSX_FILEIO_SILO_SILOFILE_HPP_
21 
22 #include "common/DataTypes.hpp"
23 #include "constitutive/solid/SolidBase.hpp"
24 #include "mpiCommunications/MpiWrapper.hpp"
26 #include "mesh/CellElementSubRegion.hpp"
29 
30 #include "silo.h"
31 
33 struct _PMPIO_baton_t;
34 
36 typedef _PMPIO_baton_t PMPIO_baton_t;
37 
38 namespace geosx
39 {
40 
41 class DomainPartition;
42 class MeshLevel;
43 
44 namespace constitutive
45 {
46 class ConstitutiveManager;
47 }
48 
49 // *********************************************************************************************************************
50 // *********************************************************************************************************************
55 class SiloFile
56 {
57 
58 public:
59 
61  SiloFile();
62 
64  virtual ~SiloFile();
65 
69  void MakeSiloDirectories();
70 
75  void Initialize( int const numGroups=1 );
76 
80  void Finish();
81 
87  int groupRank( int const i ) const;
88 
99  void WaitForBatonWrite( int const domainNumber,
100  int const cycleNum,
101  integer const eventCounter,
102  bool const isRestart );
103 
109  void WaitForBaton( int const domainNumber, string const & restartFileName );
110 
114  void HandOffBaton();
115 
116 
122  void MakeSubDirectory( string const & subdir, string const & rootdir )
123  {
124  int const rank = MpiWrapper::Comm_rank( MPI_COMM_GEOSX );
125 
126  // char dirname[100];
127  if( rank == 0 )
128  {
129 // DBGetDir(m_dbBaseFilePtr, dirname );
130  DBMkDir( m_dbBaseFilePtr, rootdir.c_str());
131  }
132 
133 // DBGetDir (m_dbFilePtr, dirname );
134  DBMkDir( m_dbFilePtr, subdir.c_str());
135  }
136 
158  void WriteMeshObject( string const & meshName,
159  const localIndex nnodes,
160  real64 * coords[3],
161  const globalIndex * globalNodeNum,
162  char const * const ghostNodeName,
163  char const * const ghostZoneName,
164  int const numShapes,
165  int const * shapecnt,
166  const localIndex * const * const meshConnectivity,
167  const globalIndex * const * const globalElementNum,
168  int const * const shapetype,
169  int const * const shapesize,
170  int const cycleNumber,
171  real64 const problemTime );
172 
192  void WritePolygonMeshObject( const std::string & meshName,
193  const localIndex nnodes,
194  real64 * coords[3],
195  const globalIndex * dummy1,
196  const int numRegions,
197  const int * shapecnt,
198  const localIndex * const * const meshConnectivity,
199  const globalIndex * const * const globalElementNum,
200  const int * const * const dummy2,
201  const int * const shapetype,
202  const int * const shapesize,
203  const int cycleNumber,
204  const real64 problemTime,
205  const int lnodelist );
206 
214  void WriteDomainPartition( DomainPartition const & domain,
215  int const cycleNum,
216  real64 const problemTime,
217  bool const isRestart );
218 
233  void WriteElementMesh( ElementRegionBase const & elementRegion,
234  NodeManager const * const nodeManager,
235  string const & meshName,
236  const localIndex nnodes,
237  real64 * coords[3],
238  globalIndex const * const globalNodeNum,
239  char const * const ghostNodeFlag,
240  int const cycleNumber,
241  real64 const problemTime,
242  bool & writeArbitraryPolygon );
243 
251  void WriteMeshLevel( MeshLevel const * const meshLevel,
252  int const cycleNum,
253  real64 const problemTime,
254  bool const isRestart );
255 
263  void WritePointMesh( string const & meshName,
264  const localIndex numPoints,
265  real64 * coords[3],
266  int const cycleNumber,
267  real64 const problemTime );
268 
279  void WriteBeamMesh( string const & meshName,
280  const localIndex nnodes,
281  real64 * coords[3],
282  const localIndex_array & node1,
283  const localIndex_array & node2,
284  int const cycleNumber,
285  real64 const problemTime );
286 
296  void WriteBeamMesh( string const & meshName,
297  const localIndex nnodes,
298  real64 * coords[3],
299  integer_array & nodelist,
300  int const cycleNumber,
301  real64 const problemTime );
302 
310  void WriteMaterialMapsFullStorage( ElementRegionBase const & elementRegion,
311  string const & meshName,
312  string_array const & regionMaterialList,
313  int const cycleNumber,
314  real64 const problemTime );
326  void WriteGroupSilo( dataRepository::Group const * group,
327  string const & siloDirName,
328  string const & meshname,
329  int const centering,
330  int const cycleNum,
331  real64 const problemTime,
332  bool const isRestart,
333  const localIndex_array & mask );
334 
343  void WriteElementRegionSilo( ElementRegionBase const & elemRegion,
344  string const & siloDirName,
345  string const & meshName,
346  int const cycleNum,
347  real64 const problemTime,
348  bool const isRestart );
361  template< typename OUTPUTTYPE >
362  void WriteWrappersToSilo( string const & meshname,
363  const dataRepository::Group::wrapperMap & wrappers,
364  int const centering,
365  int const cycleNum,
366  real64 const problemTime,
367  bool const isRestart,
368  string const & multiRoot,
369  const localIndex_array & mask );
370 
381  template< typename OUTTYPE, typename TYPE >
382  void WriteDataField( string const & meshName,
383  string const & fieldName,
384  arrayView1d< TYPE const > const & field,
385  int const centering,
386  int const cycleNumber,
387  real64 const problemTime,
388  string const & multiRoot );
389 
400  template< typename OUTTYPE, typename TYPE, int USD >
401  void WriteDataField( string const & meshName,
402  string const & fieldName,
403  arrayView2d< TYPE const, USD > const & field,
404  int const centering,
405  int const cycleNumber,
406  real64 const problemTime,
407  string const & multiRoot );
408 
418  template< typename OUTTYPE, typename TYPE, int USD >
419  void WriteDataField( string const & meshName,
420  string const & fieldName,
421  arrayView3d< TYPE const, USD > const & field,
422  int const centering,
423  int const cycleNumber,
424  real64 const problemTime,
425  string const & multiRoot );
426 
438  template< typename OUTTYPE, typename TYPE, int NDIM, int USD >
439  void WriteDataField( string const & meshName,
440  string const & fieldName,
442  int const siloTensorRank,
443  int const centering,
444  int const cycleNumber,
445  real64 const problemTime,
446  string const & multiRoot );
447 
460  template< typename OUTTYPE, typename TYPE >
461  void WriteMaterialDataField( string const & meshName,
462  string const & fieldName,
463  array1d< array1d< arrayView2d< TYPE const > > > const & field,
464  ElementRegionBase const & elemRegion,
465  int const centering,
466  int const cycleNumber,
467  real64 const problemTime,
468  string const & multiRoot,
469  string_array const & materialNames );
470 
482  template< typename OUTTYPE, typename TYPE >
483  void WriteMaterialDataField2d( string const & meshName,
484  string const & fieldName,
485  ElementRegionBase const & elemRegion,
486  int const centering,
487  int const cycleNumber,
488  real64 const problemTime,
489  string const & multiRoot,
490  string_array const & materialNames );
491 
503  template< typename OUTTYPE, typename TYPE >
504  void WriteMaterialDataField3d( string const & meshName,
505  string const & fieldName,
506  ElementRegionBase const & elemRegion,
507  int const centering,
508  int const cycleNumber,
509  real64 const problemTime,
510  string const & multiRoot,
511  string_array const & materialNames );
512 
524  template< typename OUTTYPE, typename TYPE >
525  void WriteMaterialDataField4d( string const & meshName,
526  string const & fieldName,
527  ElementRegionBase const & elemRegion,
528  int const centering,
529  int const cycleNumber,
530  real64 const problemTime,
531  string const & multiRoot,
532  string_array const & materialNames );
533 
541  void WriteMaterialVarDefinition( string const & subDir,
542  string const & matDir,
543  localIndex const matIndex,
544  string const & fieldName );
545 
550  void WriteStressVarDefinition( string const & MatDir );
551 
557  void WriteVectorVarDefinition( string const & fieldName,
558  string const & subDirectory );
559 
565  int GetMeshType( string const & meshName ) const;
566 
577  template< typename CBF >
578  void WriteMultiXXXX( const DBObjectType type, CBF DBPutMultiCB,
579  int const centering, string const name, int const cycleNumber,
580  string const & multiRoot, const DBoptlist * optlist = nullptr );
581 
582 
591  void ClearEmptiesFromMultiObjects( int const cycleNum );
592 
597  void setNumGroups( int const numGroups )
598  {
599  m_numGroups = numGroups;
600  }
601 
606  void setPlotLevel( int const plotLevel )
607  {
608  m_plotLevel = dataRepository::toPlotLevel( plotLevel );
609  }
610 
615  void setWriteEdgeMesh( int const val )
616  {
617  m_writeEdgeMesh = val;
618  }
619 
624  void setWriteFaceMesh( int const val )
625  {
626  m_writeFaceMesh = val;
627  }
628 
633  void setWriteCellElementMesh( int const val )
634  {
635  m_writeCellElementMesh = val;
636  }
637 
642  void setWriteFaceElementMesh( int const val )
643  {
644  m_writeFaceElementMesh = val;
645  }
646 
651  void setPlotFileRoot( string const & fileRoot )
652  {
653  m_plotFileRoot = fileRoot;
654  }
655 
656 private:
657 
659  DBfile * m_dbFilePtr;
660 
662  DBfile * m_dbBaseFilePtr;
663 
665  int m_numGroups;
666 
668  PMPIO_baton_t *m_baton;
669 
671  int const m_driver;
672 
674  string m_plotFileRoot;
675 
676  string m_restartFileRoot;
677 
678  string const m_siloDirectory = "siloFiles";
679 
680  string const m_siloDataSubDirectory = "data";
681 
682  string m_fileName;
683 
684  string m_baseFileName;
685 
686  std::vector< std::string > m_emptyMeshes;
687 // string_array m_emptyMaterials;
688  std::vector< std::string > m_emptyVariables;
689 
690  integer m_writeEdgeMesh;
691  integer m_writeFaceMesh;
692  integer m_writeCellElementMesh;
693  integer m_writeFaceElementMesh;
694 
695  dataRepository::PlotLevel m_plotLevel;
696 
697  bool m_ghostFlags;
698 };
699 
703 namespace SiloFileUtilities
704 {
713 template< typename OUTTYPE >
714 int DB_TYPE();
715 
723 template< typename TYPE >
725 
730 template< typename TYPE >
731 int GetTensorRank();
732 
740 template< typename OUTTYPE, typename TYPE >
741 OUTTYPE CastField( const TYPE & field, int const i = 0 ); // avoids compiler
742  // warning
749 template<> inline real64 CastField( R1Tensor const & field, int const i )
750 {
751  return field[i];
752 }
753 
761 template<> inline int CastField< int, int >( const int & field, int const dummy )
762 {
763  GEOSX_UNUSED_VAR( dummy );
764  return field;
765 }
766 
774 template<> inline long int CastField< long int, long int >( const long int & field, int const dummy )
775 {
776  GEOSX_UNUSED_VAR( dummy );
777  return field;
778 }
779 
787 template<> inline int CastField< int, long int >( const long int & field, int const dummy )
788 {
789  GEOSX_UNUSED_VAR( dummy );
790  return LvArray::integerConversion< int >( field );
791 }
792 
800 template<> inline long long int CastField< long long int, long long int >( const long long int & field, int const dummy )
801 {
802  GEOSX_UNUSED_VAR( dummy );
803  return field;
804 }
805 
813 template<> inline int CastField< int, long long int >( const long long int & field, int const dummy )
814 {
815  GEOSX_UNUSED_VAR( dummy );
816  return LvArray::integerConversion< int >( field );
817 }
818 
826 template<> inline real64 CastField< real64, real64 >( const real64 & field, int const dummy )
827 {
828  GEOSX_UNUSED_VAR( dummy );
829  return field;
830 }
831 
839 template<> inline float CastField< float, real64 >( const real64 & field, int const dummy )
840 {
841  GEOSX_UNUSED_VAR( dummy );
842  return static_cast< float >(field);
843 }
844 
855 template< typename TYPE >
856 void SetVariableNames( string const & fieldName, string_array & varnamestring, char const * varnames[] );
857 
858 
859 }
860 
861 
862 
863 }
864 #endif /* GEOSX_FILEIO_SILO_SILOFILE_HPP_ */
void setNumGroups(int const numGroups)
Sets the number of individual Silo files to generate.
Definition: SiloFile.hpp:597
void setPlotFileRoot(string const &fileRoot)
Sets root of the filename that will be read/written.
Definition: SiloFile.hpp:651
long long int globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:128
void setPlotLevel(int const plotLevel)
Sets the plot level option.
Definition: SiloFile.hpp:606
void MakeSubDirectory(string const &subdir, string const &rootdir)
Make a subdirectory within the silo file.
Definition: SiloFile.hpp:122
Class facilitating the representation of a multi-level discretization of a MeshBody.
Definition: MeshLevel.hpp:38
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data...
Definition: NodeManager.hpp:47
_PMPIO_baton_t PMPIO_baton_t
Type alias for _PMPIO_baton_t struct.
Definition: SiloFile.hpp:33
long int CastField< long int, long int >(const long int &field, int const dummy)
Definition: SiloFile.hpp:774
This class serves to provide a "view" of a multidimensional array.
Definition: ArrayView.hpp:67
long long int CastField< long long int, long long int >(const long long int &field, int const dummy)
Definition: SiloFile.hpp:800
double real64
64-bit floating point type.
Definition: DataTypes.hpp:136
int CastField< int, int >(const int &field, int const dummy)
Definition: SiloFile.hpp:761
int CastField< int, long int >(const long int &field, int const dummy)
Definition: SiloFile.hpp:787
void setWriteCellElementMesh(int const val)
Sets the cell element mesh output option.
Definition: SiloFile.hpp:633
void setWriteEdgeMesh(int const val)
Sets the edge mesh output option.
Definition: SiloFile.hpp:615
Lightweight wrapper around a c-array.
Definition: Tensor.hpp:29
int MPI_COMM_GEOSX
Global MPI communicator used by GEOSX.
real64 CastField< real64, real64 >(const real64 &field, int const dummy)
Definition: SiloFile.hpp:826
float CastField< float, real64 >(const real64 &field, int const dummy)
Definition: SiloFile.hpp:839
PlotLevel toPlotLevel(int const val)
Function to get a PlotLevel enum from an int.
std::int32_t integer
Signed integer type.
Definition: DataTypes.hpp:122
The ElementRegionBase is the base class to manage the data stored at the element level.
#define GEOSX_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:78
void SetVariableNames(string const &fieldName, string_array &varnamestring, char const *varnames[])
std::ptrdiff_t localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:125
real64 CastField(R1Tensor const &field, int const i)
Specialization for R1Tensor.
Definition: SiloFile.hpp:749
void setWriteFaceElementMesh(int const val)
Sets the face element mesh output option.
Definition: SiloFile.hpp:642
std::string string
String type.
Definition: DataTypes.hpp:131
int CastField< int, long long int >(const long long int &field, int const dummy)
Definition: SiloFile.hpp:813
This class provides a fixed dimensional resizeable array interface in addition to an interface simila...
Definition: Array.hpp:55
Partition of the decomposed physical domain. It also manages the connexion information to its neighbo...
void setWriteFaceMesh(int const val)
Sets the face mesh output option.
Definition: SiloFile.hpp:624