00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046 #include "Xpetra_EpetraCrsGraph.hpp"
00047
00048 #include "Xpetra_Exceptions.hpp"
00049 #include "Xpetra_Utils.hpp"
00050 #include "Xpetra_EpetraExport.hpp"
00051 #include "Xpetra_EpetraImport.hpp"
00052
00053 namespace Xpetra {
00054
00055
00056 const Epetra_CrsGraph & toEpetra(const RCP< const CrsGraph<int, int> > &graph) {
00057 XPETRA_RCP_DYNAMIC_CAST(const EpetraCrsGraph, graph, epetraGraph, "toEpetra");
00058 return *(epetraGraph->getEpetra_CrsGraph());
00059 }
00060
00061 EpetraCrsGraph::EpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype, const Teuchos::RCP< Teuchos::ParameterList > &plist)
00062 : graph_(Teuchos::rcp(new Epetra_CrsGraph(Copy, toEpetra(rowMap), maxNumEntriesPerRow, toEpetra(pftype)))) { }
00063
00064
00065
00066
00067
00068 EpetraCrsGraph::EpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype, const Teuchos::RCP< Teuchos::ParameterList > &plist)
00069 : graph_(Teuchos::rcp(new Epetra_CrsGraph(Copy, toEpetra(rowMap), toEpetra(colMap), maxNumEntriesPerRow, toEpetra(pftype)))) { }
00070
00071
00072
00073
00074
00075 void EpetraCrsGraph::insertGlobalIndices(int globalRow, const ArrayView<const int> &indices) {
00076 XPETRA_MONITOR("EpetraCrsGraph::insertGlobalIndices");
00077
00078 int* indices_rawPtr = const_cast<int*>(indices.getRawPtr());
00079 XPETRA_ERR_CHECK(graph_->InsertGlobalIndices(globalRow, indices.size(), indices_rawPtr));
00080 }
00081
00082 void EpetraCrsGraph::insertLocalIndices(int localRow, const ArrayView<const int> &indices) {
00083 XPETRA_MONITOR("EpetraCrsGraph::insertLocalIndices");
00084
00085 int* indices_rawPtr = const_cast<int*>(indices.getRawPtr());
00086 XPETRA_ERR_CHECK(graph_->InsertMyIndices(localRow, indices.size(), indices_rawPtr));
00087 }
00088
00089 void EpetraCrsGraph::getGlobalRowView(int GlobalRow, ArrayView<const int> &Indices) const {
00090 XPETRA_MONITOR("EpetraCrsGraph::getGlobalRowView");
00091
00092 int numEntries;
00093 int * eIndices;
00094
00095 XPETRA_ERR_CHECK(graph_->ExtractGlobalRowView(GlobalRow, numEntries, eIndices));
00096 if (numEntries == 0) { eIndices = NULL; }
00097
00098 Indices = ArrayView<const int>(eIndices, numEntries);
00099 }
00100
00101 void EpetraCrsGraph::getLocalRowView(int LocalRow, ArrayView<const int> &indices) const {
00102 XPETRA_MONITOR("EpetraCrsGraph::getLocalRowView");
00103
00104 int numEntries;
00105 int * eIndices;
00106
00107 XPETRA_ERR_CHECK(graph_->ExtractMyRowView(LocalRow, numEntries, eIndices));
00108 if (numEntries == 0) { eIndices = NULL; }
00109
00110 indices = ArrayView<const int>(eIndices, numEntries);
00111 }
00112
00113 void EpetraCrsGraph::fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > ¶ms){
00114 XPETRA_MONITOR("EpetraCrsGraph::fillComplete");
00115
00116 graph_->FillComplete(toEpetra(domainMap), toEpetra(rangeMap));
00117 bool doOptimizeStorage = true;
00118 if (params != null && params->get("Optimize Storage",true) == false) doOptimizeStorage = false;
00119 if (doOptimizeStorage) graph_->OptimizeStorage();
00120 }
00121
00122 void EpetraCrsGraph::fillComplete(const RCP< ParameterList > ¶ms) {
00123 XPETRA_MONITOR("EpetraCrsGraph::fillComplete");
00124
00125 graph_->FillComplete();
00126 bool doOptimizeStorage = true;
00127 if (params != null && params->get("Optimize Storage",true) == false) doOptimizeStorage = false;
00128 if (doOptimizeStorage) graph_->OptimizeStorage();
00129 }
00130
00131 std::string EpetraCrsGraph::description() const { XPETRA_MONITOR("EpetraCrsGraph::description"); return "NotImplemented"; }
00132
00133 void EpetraCrsGraph::describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const {
00134 XPETRA_MONITOR("EpetraCrsGraph::describe");
00135
00136 out << "EpetraCrsGraph::describe : Warning, verbosity level is ignored by this method." << std::endl;
00137 const Epetra_BlockMap rowmap = graph_->RowMap();
00138 if (rowmap.Comm().MyPID() == 0) out << "** EpetraCrsGraph **\n\nrowmap" << std::endl;
00139 rowmap.Print(out);
00140 graph_->Print(out);
00141 }
00142
00143
00144 RCP< const CrsGraph<int, int> > toXpetra(const Epetra_CrsGraph &g) {
00145 RCP<const Epetra_CrsGraph> const_graph = rcp(new Epetra_CrsGraph(g));
00146
00147 RCP<Epetra_CrsGraph> graph = Teuchos::rcp_const_cast<Epetra_CrsGraph>(const_graph);
00148 return rcp( new Xpetra::EpetraCrsGraph(graph) );
00149 }
00150
00151
00152
00153 void EpetraCrsGraph::doImport(const DistObject<int, int, int> &source,
00154 const Import<int, int> &importer, CombineMode CM) {
00155 XPETRA_MONITOR("EpetraCrsGraph::doImport");
00156
00157 XPETRA_DYNAMIC_CAST(const EpetraCrsGraph, source, tSource, "Xpetra::EpetraCrsGraph::doImport only accept Xpetra::EpetraCrsGraph as input arguments.");
00158 XPETRA_DYNAMIC_CAST(const EpetraImport, importer, tImporter, "Xpetra::EpetraCrsGraph::doImport only accept Xpetra::EpetraImport as input arguments.");
00159
00160 RCP<const Epetra_CrsGraph> v = tSource.getEpetra_CrsGraph();
00161 int err = graph_->Import(*v, *tImporter.getEpetra_Import(), toEpetra(CM));
00162 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, "Catch error code returned by Epetra.");
00163 }
00164
00165 void EpetraCrsGraph::doExport(const DistObject<int, int, int> &dest,
00166 const Import<int, int>& importer, CombineMode CM) {
00167 XPETRA_MONITOR("EpetraCrsGraph::doExport");
00168
00169 XPETRA_DYNAMIC_CAST(const EpetraCrsGraph, dest, tDest, "Xpetra::EpetraCrsGraph::doImport only accept Xpetra::EpetraCrsGraph as input arguments.");
00170 XPETRA_DYNAMIC_CAST(const EpetraImport, importer, tImporter, "Xpetra::EpetraCrsGraph::doImport only accept Xpetra::EpetraImport as input arguments.");
00171
00172 RCP<const Epetra_CrsGraph> v = tDest.getEpetra_CrsGraph();
00173 int err = graph_->Export(*v, *tImporter.getEpetra_Import(), toEpetra(CM));
00174 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, "Catch error code returned by Epetra.");
00175 }
00176
00177 void EpetraCrsGraph::doImport(const DistObject<int, int, int> &source,
00178 const Export<int, int>& exporter, CombineMode CM) {
00179 XPETRA_MONITOR("EpetraCrsGraph::doImport");
00180
00181 XPETRA_DYNAMIC_CAST(const EpetraCrsGraph, source, tSource, "Xpetra::EpetraCrsGraph::doImport only accept Xpetra::EpetraCrsGraph as input arguments.");
00182 XPETRA_DYNAMIC_CAST(const EpetraExport, exporter, tExporter, "Xpetra::EpetraCrsGraph::doImport only accept Xpetra::EpetraImport as input arguments.");
00183
00184 RCP<const Epetra_CrsGraph> v = tSource.getEpetra_CrsGraph();
00185 int err = graph_->Import(*v, *tExporter.getEpetra_Export(), toEpetra(CM));
00186 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, "Catch error code returned by Epetra.");
00187
00188 }
00189
00190 void EpetraCrsGraph::doExport(const DistObject<int, int, int> &dest,
00191 const Export<int, int>& exporter, CombineMode CM) {
00192 XPETRA_MONITOR("EpetraCrsGraph::doExport");
00193
00194 XPETRA_DYNAMIC_CAST(const EpetraCrsGraph, dest, tDest, "Xpetra::EpetraCrsGraph::doImport only accept Xpetra::EpetraCrsGraph as input arguments.");
00195 XPETRA_DYNAMIC_CAST(const EpetraExport, exporter, tExporter, "Xpetra::EpetraCrsGraph::doImport only accept Xpetra::EpetraImport as input arguments.");
00196
00197 RCP<const Epetra_CrsGraph> v = tDest.getEpetra_CrsGraph();
00198 int err = graph_->Export(*v, *tExporter.getEpetra_Export(), toEpetra(CM));
00199 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, "Catch error code returned by Epetra.");
00200 }
00201
00202 }