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 #ifndef XPETRA_EPETRACRSGRAPH_HPP
00047 #define XPETRA_EPETRACRSGRAPH_HPP
00048
00049
00050
00051 #include "Xpetra_EpetraConfigDefs.hpp"
00052
00053 #include "Xpetra_CrsGraph.hpp"
00054
00055 #include "Xpetra_EpetraMap.hpp"
00056 #include "Xpetra_EpetraImport.hpp"
00057 #include "Xpetra_EpetraUtils.hpp"
00058
00059 #include <Epetra_CrsGraph.h>
00060
00061 namespace Xpetra {
00062
00063
00064 RCP< const CrsGraph<int, int> > toXpetra(const Epetra_CrsGraph& graph);
00065 const Epetra_CrsGraph & toEpetra(const RCP< const CrsGraph<int, int> > &graph);
00066
00067
00068 #ifndef DOXYGEN_SHOULD_SKIP_THIS
00069
00070 template <class S, class LO, class GO, class N, class SpMatOps>
00071 class CrsMatrix;
00072 #endif
00073
00074 class EpetraCrsGraph
00075 : public CrsGraph<int, int>
00076 {
00077
00078 typedef int LocalOrdinal;
00079 typedef int GlobalOrdinal;
00080 typedef Kokkos::DefaultNode::DefaultNodeType Node;
00082 typedef Map<LocalOrdinal,GlobalOrdinal,Node> map_type;
00083
00084 public:
00085
00087
00088
00090 EpetraCrsGraph(const RCP< const map_type > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const RCP< ParameterList > ¶ms=null);
00091
00093 EpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const RCP< ParameterList > ¶ms=null);
00094
00096 EpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const RCP< ParameterList > ¶ms=null);
00097
00099 EpetraCrsGraph(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const RCP< ParameterList > ¶ms=null);
00100
00102 virtual ~EpetraCrsGraph() { }
00103
00105
00107
00108
00110 void insertGlobalIndices(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &indices);
00111
00113 void insertLocalIndices(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &indices);
00114
00116 void removeLocalIndices(LocalOrdinal localRow) { XPETRA_MONITOR("EpetraCrsGraph::removeLocalIndices"); graph_->RemoveMyIndices(localRow); }
00117
00119
00121
00122
00124 void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > ¶ms=null);
00125
00127 void fillComplete(const RCP< ParameterList > ¶ms=null);
00128
00130
00132
00133
00135 const RCP< const Comm< int > > getComm() const { XPETRA_MONITOR("EpetraCrsGraph::getComm"); return toXpetra(graph_->Comm()); }
00136
00138 const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const { XPETRA_MONITOR("EpetraCrsGraph::getRowMap"); return toXpetra(graph_->RowMap()); }
00139
00141 const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const { XPETRA_MONITOR("EpetraCrsGraph::getColMap"); return toXpetra(graph_->ColMap()); }
00142
00144 const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const { XPETRA_MONITOR("EpetraCrsGraph::getDomainMap"); return toXpetra(graph_->DomainMap()); }
00145
00147 const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const { XPETRA_MONITOR("EpetraCrsGraph::getRangeMap"); return toXpetra(graph_->RangeMap()); }
00148
00150 RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > getImporter() const { XPETRA_MONITOR("EpetraCrsGraph::getImporter"); return toXpetra(graph_->Importer()); }
00151
00153 global_size_t getGlobalNumRows() const { XPETRA_MONITOR("EpetraCrsGraph::getGlobalNumRows"); return graph_->NumGlobalRows(); }
00154
00156 global_size_t getGlobalNumCols() const { XPETRA_MONITOR("EpetraCrsGraph::getGlobalNumCols"); return graph_->NumGlobalCols(); }
00157
00159 size_t getNodeNumRows() const { XPETRA_MONITOR("EpetraCrsGraph::getNodeNumRows"); return graph_->NumMyRows(); }
00160
00162 size_t getNodeNumCols() const { XPETRA_MONITOR("EpetraCrsGraph::getNodeNumCols"); return graph_->NumMyCols(); }
00163
00165 GlobalOrdinal getIndexBase() const { XPETRA_MONITOR("EpetraCrsGraph::getIndexBase"); return graph_->IndexBase(); }
00166
00168 global_size_t getGlobalNumEntries() const { XPETRA_MONITOR("EpetraCrsGraph::getGlobalNumEntries"); return graph_->NumGlobalEntries(); }
00169
00171 size_t getNodeNumEntries() const { XPETRA_MONITOR("EpetraCrsGraph::getNodeNumEntries"); return graph_->NumMyEntries(); }
00172
00174 size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const { XPETRA_MONITOR("EpetraCrsGraph::getNumEntriesInGlobalRow"); return graph_->NumGlobalIndices(globalRow); }
00175
00177 size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const { XPETRA_MONITOR("EpetraCrsGraph::getNumEntriesInLocalRow"); return graph_->NumMyIndices(localRow); }
00178
00180 size_t getNumAllocatedEntriesInGlobalRow(GlobalOrdinal globalRow) const { XPETRA_MONITOR("EpetraCrsGraph::getNumAllocatedEntriesInGlobalRow"); return graph_->NumAllocatedGlobalIndices(globalRow); }
00181
00183 size_t getNumAllocatedEntriesInLocalRow(LocalOrdinal localRow) const { XPETRA_MONITOR("EpetraCrsGraph::getNumAllocatedEntriesInLocalRow"); return graph_->NumAllocatedMyIndices(localRow); }
00184
00186 global_size_t getGlobalNumDiags() const { XPETRA_MONITOR("EpetraCrsGraph::getGlobalNumDiags"); return graph_->NumGlobalDiagonals(); }
00187
00189 size_t getNodeNumDiags() const { XPETRA_MONITOR("EpetraCrsGraph::getNodeNumDiags"); return graph_->NumMyDiagonals(); }
00190
00192 size_t getGlobalMaxNumRowEntries() const { XPETRA_MONITOR("EpetraCrsGraph::getGlobalMaxNumRowEntries"); return graph_->GlobalMaxNumIndices(); }
00193
00195 size_t getNodeMaxNumRowEntries() const { XPETRA_MONITOR("EpetraCrsGraph::getNodeMaxNumRowEntries"); return graph_->MaxNumIndices(); }
00196
00198 bool hasColMap() const { XPETRA_MONITOR("EpetraCrsGraph::hasColMap"); return graph_->HaveColMap(); }
00199
00201 bool isLowerTriangular() const { XPETRA_MONITOR("EpetraCrsGraph::isLowerTriangular"); return graph_->LowerTriangular(); }
00202
00204 bool isUpperTriangular() const { XPETRA_MONITOR("EpetraCrsGraph::isUpperTriangular"); return graph_->UpperTriangular(); }
00205
00207 bool isLocallyIndexed() const { XPETRA_MONITOR("EpetraCrsGraph::isLocallyIndexed"); return graph_->IndicesAreLocal(); }
00208
00210 bool isGloballyIndexed() const { XPETRA_MONITOR("EpetraCrsGraph::isGloballyIndexed"); return graph_->IndicesAreGlobal(); }
00211
00213 bool isFillComplete() const { XPETRA_MONITOR("EpetraCrsGraph::isFillComplete"); return graph_->Filled(); }
00214
00216 bool isStorageOptimized() const { XPETRA_MONITOR("EpetraCrsGraph::isStorageOptimized"); return graph_->StorageOptimized(); }
00217
00219 void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &Indices) const;
00220
00222 void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices) const;
00223
00225
00227
00228
00230 std::string description() const;
00231
00233 void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const;
00234
00236
00238
00239
00241 const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const { XPETRA_MONITOR("EpetraCrsGraph::getMap"); return toXpetra(graph_->Map()); }
00242
00244 void doImport(const DistObject<GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node> &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM);
00245
00247 void doExport(const DistObject<GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node> &dest, const Import< LocalOrdinal, GlobalOrdinal, Node >& importer, CombineMode CM);
00248
00250 void doImport(const DistObject<GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node> &source, const Export< LocalOrdinal, GlobalOrdinal, Node >& exporter, CombineMode CM);
00251
00253 void doExport(const DistObject<GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node> &dest, const Export< LocalOrdinal, GlobalOrdinal, Node >& exporter, CombineMode CM);
00254
00256
00258
00259
00261 EpetraCrsGraph(const Teuchos::RCP<Epetra_CrsGraph> &graph) : graph_(graph) { }
00262
00264 RCP< const Epetra_CrsGraph> getEpetra_CrsGraph() const { return graph_; }
00265
00267
00268 private:
00269
00270 RCP<Epetra_CrsGraph> graph_;
00271
00272 };
00273
00274 }
00275
00276 #endif // XPETRA_EPETRACRSGRAPH_HPP