[Trilinos-Users] RCP and loading CrsGraph from HDF5 file
Bartlett, Roscoe A
rabartl at sandia.gov
Fri Sep 26 05:44:43 MDT 2008
Davood,
What you have will not work and will most likely cause a segfault.
When working with RCP, it is critical that you understand exactly what is going on and how it interacts with the C++ type system.
For a basic introduction on Teuchos::RCP, see:
http://www.cs.sandia.gov/~rabartl/RefCountPtrBeginnersGuideSAND.pdf
For some more detailed information, see:
http://www.cs.sandia.gov/~rabartl/rabartl_TUG_RCP_Talk.pdf
Note that all of Trilinos is not updated to use Teuchos::RCP so in cases like with the class EpetraExt::HDF5, you have to be very careful to get them into an RCP correctly.
>From looking at the code in EpetraExt_HDF5.cpp, it seems that this is how you would successfully read into RCP wrapped objects:
using Teuchos::RCP;
using Teuchos::rcp;
HDF5.Open("CRSGraph.h5");
Epetra_Map *tmpRowMap = 0;
HDF5.Read("Map", tmpRowMap);
RCP<Epetra_Map> rowMap = rcp(tmpRowMap);
tmpRowMap = 0;
Epetra_CrsGraph *tmpCrsGraph = 0;
HDF5.Read("Graph", tmpCrsGraph) ;
RCP<Epetra_CrsGraph> crsGraph = rcp(tmpCrsGraph);
tmpCrsGraph = 0;
HDF5.Close();
I have not compiled and tested this code so I don't know if this is 100% correct but it should be very close.
If you look through the documents above and through the file EpetraExt_HDF5.cpp, hopefully you should be able to see what is going on here.
Cheers,
- Ross
P.S. You should almost never employ a C-style cast like (Epetra *&). Such an attempt is almost always going to silence the C++ compiler but will result in wrong code that will segfault your program. In the future, if you are tempted to resort to a C-style cast like this then you need to stop and investigate what is really going on.
---------------------------------------------------------
Dr. Roscoe A. Bartlett
Senior Member of the Technical Staff
Sandia National Laboratories
Phone: (505) 275-6147
________________________________
From: trilinos-users-bounces at software.sandia.gov [mailto:trilinos-users-bounces at software.sandia.gov] On Behalf Of Davood Ansari
Sent: Friday, September 26, 2008 1:45 AM
To: trilinos-users at software.sandia.gov
Subject: [Trilinos-Users] RCP and loading CrsGraph from HDF5 file
Hi all
This should seem a bit dumb.
I am not good at OOp so I put it like this.
In the following piece, I had to use the type cast - (Epetra_Map*&) - or it would not compile.
My intention is to load a EpetraCrsGraph and a Epetra_Map into RCP pointed objects.
Would the following code handle the garbage collection right? Or should I write it differently?
Kindly Advise
Davood
Teuchos::RCP<Epetra_Map::Epetra_Map> RowMap ;
Teuchos::RCP<Epetra_CrsGraph::Epetra_CrsGraph> CrsGraph ;
HDF5.Open("CRSGraph.h5");
HDF5.Read("Map", (Epetra_Map*&)RowMap) ;
HDF5.Read("Graph", (Epetra_CrsGraph*&)CrsGraph) ;
HDF5.Close();
-------------- next part --------------
An HTML attachment was scrubbed...
URL: http://software.sandia.gov/mailman/private/trilinos-users/attachments/20080926/eae2ca6e/attachment-0001.html
More information about the Trilinos-Users
mailing list