Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
HDF5Export_def.hpp
1#ifndef HDF5EXPORT_DEF_hpp
2#define HDF5EXPORT_DEF_hpp
3
4#include "HDF5Export_decl.hpp"
5
6
7namespace FEDD {
8
9template<class SC,class LO,class GO,class NO>
10HDF5Export<SC,LO,GO,NO>::HDF5Export(MapConstPtr_Type writeMap, std::string outputFilename):
11hdf5exporter_(),
12comm_(),
13commEpetra_()
14{
15
16 Teuchos::RCP<const Teuchos::MpiComm<int> > mpiComm = Teuchos::rcp_dynamic_cast<const Teuchos::MpiComm<int> >( writeMap->getComm() );
17 commEpetra_.reset( new Epetra_MpiComm( *mpiComm->getRawMpiComm() ) );
18
19 // We convert the current map to an eptra map
20 Teuchos::ArrayView< const GO > indices = writeMap->getNodeElementList();
21 int* intGlobIDs = new int[indices.size()];
22 for (int i=0; i<indices.size(); i++) {
23 intGlobIDs[i] = (int) indices[i];
24 }
25
26 int nmbPointsGlob = writeMap->getGlobalNumElements(); // Number of global points
27
28 EpetraMapPtr_Type mapEpetra = Teuchos::rcp(new Epetra_Map((int)nmbPointsGlob,indices.size(),intGlobIDs,0,*commEpetra_));
29
30 writeMap_ = mapEpetra; // Defining the read map. All different variables that might be written to the .h5 file have the same map
31
32 hdf5exporter_.reset( new HDF5_Type(*commEpetra_) ); // Building HDF5 Exporter
33
34 hdf5exporter_->Create(outputFilename+".h5"); // Creating output file with the 'outoutFilename'
35
36 outputFilename_ = outputFilename;
37}
38
39template<class SC,class LO,class GO,class NO>
40void HDF5Export<SC,LO,GO,NO>::writeVariablesHDF5(std::string varName,MultiVectorConstPtr_Type writeVector){
41
42 EpetraMVPtr_Type u_export(new Epetra_MultiVector(*(writeMap_),1)); // Epetra export vector
43
44 TEUCHOS_TEST_FOR_EXCEPTION( std::abs(writeMap_->NumMyElements() - writeVector->getLocalLength()) > 1e-12, std::logic_error, " The local length of map does not match the local mv length. Map and MultiVector are not compatible");
45
46 // We need to write the contents of the writeVector into the Epetra export vector: Convert Xpetra -> Epetra
47 Teuchos::ArrayRCP<const SC> tmpData = writeVector->getData(0);
48 for (int i=0; i<writeVector->getLocalLength(); i++) {
49 u_export->ReplaceMyValue( i, 0, tmpData[i] );
50 }
51
52 hdf5exporter_->Write(varName,*u_export); // Writing u_export as variable 'varName' in file
53
54 if(writeVector->getMap()->getComm()->getRank() == 0 )
55 std::cout << " HDF5_Export:: Exporting to file " << outputFilename_ << " with variable name " << varName << std::endl;
56
57 hdf5exporter_->Flush();
58
59}
60
61template<class SC,class LO,class GO,class NO>
63 hdf5exporter_->Close();
64}
65
66}
67#endif
void closeExporter()
Closing Exporter.
Definition HDF5Export_def.hpp:62
HDF5Export(MapConstPtr_Type writeMap, std::string outputFilename)
Constructor for HDF5 Exporter.
Definition HDF5Export_def.hpp:10
void writeVariablesHDF5(std::string varName, MultiVectorConstPtr_Type writeVector)
Exporting MultiVector writeVector as HDF5 File with the variable name varName.
Definition HDF5Export_def.hpp:40
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement.cpp:5