Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
HDF5Import_def.hpp
1#ifndef HDF5IMPORT_DEF_hpp
2#define HDF5IMPORT_DEF_hpp
3
4#include "HDF5Import_decl.hpp"
5
15
16namespace FEDD {
17
18template<class SC,class LO,class GO,class NO>
19HDF5Import<SC,LO,GO,NO>::HDF5Import(MapConstPtr_Type readMap, std::string inputFilename):
21comm_(),
22commEpetra_()
23{
24 Teuchos::RCP<const Teuchos::MpiComm<int> > mpiComm = Teuchos::rcp_dynamic_cast<const Teuchos::MpiComm<int> >( readMap->getComm() );
25 commEpetra_.reset( new Epetra_MpiComm( *mpiComm->getRawMpiComm() ) ); // Communicator for epetra
26
27 // We convert the current map to an eptra map
28 Teuchos::ArrayView< const GO > indices = readMap->getNodeElementList(); // Global Ids of map on this processor
29 int* intGlobIDs = new int[indices.size()];
30 for (int i=0; i<indices.size(); i++) {
31 intGlobIDs[i] = (int) indices[i];
32 }
33
34 int nmbPointsGlob = readMap->getGlobalNumElements(); // Number of global points
35
36 EpetraMapPtr_Type mapEpetra = Teuchos::rcp(new Epetra_Map((int)nmbPointsGlob,indices.size(),intGlobIDs,0,*commEpetra_)); // Building epetra map
37
38 readMap_ = mapEpetra; // Defining the read map. All different variables that might be contained in the .h5 file must have the same map
39
40 hdf5importer_.reset( new HDF5_Type(*commEpetra_) ); // Build HDF5 importer as Epetra_Ext HDF5 Type
41
42 inputFilename_ = inputFilename; // Name of input file
43 hdf5importer_->Open(inputFilename_+".h5"); // We 'open' the file and connect to HDF5 iXpetramporter
44
45 u_import_Tpetra_.reset(new MultiVector_Type(readMap)); // The general import vector is defined via readMap
46
47
48}
49
50template<class SC,class LO,class GO,class NO>
51typename HDF5Import<SC, LO, GO, NO>::MultiVectorPtr_Type HDF5Import<SC,LO,GO,NO>::readVariablesHDF5(std::string varName){
52 // Testing if requested vairable is contained in file
53 TEUCHOS_TEST_FOR_EXCEPTION( !hdf5importer_->IsContained(varName), std::logic_error, "Requested varName: " << varName << " not contained in hdf file "<< inputFilename_ << ".h5.");
54 hdf5importer_->Read(varName,*readMap_,u_import_Epetra_); // Reading the variable 'varName' from the file and distribute according to readMap_ to u_import_
55
56 // Now we write the contents of the eptera multivector u_import to u_import_mv
57 Teuchos::ArrayRCP<SC> tmpData = u_import_Tpetra_->getDataNonConst(0);
58 for (int i=0; i<u_import_Tpetra_->getLocalLength(); i++) {
59 tmpData[i] = u_import_Epetra_->Values()[i]; // this is how you access epetra values of multi vectors
60 }
61
62 hdf5importer_->Flush();
63
64 return u_import_Tpetra_;
65
66}
67template<class SC,class LO,class GO,class NO>
68void HDF5Import<SC,LO,GO,NO>::closeImporter(){
69 hdf5importer_->Close();
70}
71
72}
73#endif
EpetraMapPtr_Type readMap_
Name of Map of import multivector.
Definition HDF5Import_decl.hpp:80
MultiVectorPtr_Type readVariablesHDF5(std::string varName)
Reading a variable 'varName' from the inputFile with inputFilename of file type HDF5.
Definition HDF5Import_def.hpp:51
std::string inputFilename_
Name of input file.
Definition HDF5Import_decl.hpp:78
Epetra_MultiVector * u_import_Epetra_
Imported MultiVector in Epetra format.
Definition HDF5Import_decl.hpp:82
MultiVectorPtr_Type u_import_Tpetra_
Imported file in Xpetra format.
Definition HDF5Import_decl.hpp:84
HDF5Import(MapConstPtr_Type readMap, std::string inputFilename)
Constructor of HDF import. Here the general setting are defined. An epetra map build based on the rea...
Definition HDF5Import_def.hpp:19
HDF5Ptr_Type hdf5importer_
HDF5 importer based on EpetraExt HDF5 importer.
Definition HDF5Import_decl.hpp:70
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement.cpp:5