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_CRSMATRIXFACTORY_HPP
00047 #define XPETRA_CRSMATRIXFACTORY_HPP
00048
00049 #include "Xpetra_ConfigDefs.hpp"
00050
00051 #include "Xpetra_CrsMatrix.hpp"
00052
00053 #ifdef HAVE_XPETRA_TPETRA
00054 #include "Xpetra_TpetraCrsMatrix.hpp"
00055 #endif
00056
00057 #ifdef HAVE_XPETRA_EPETRA
00058 #include "Xpetra_EpetraCrsMatrix.hpp"
00059 #endif
00060
00061 #include "Xpetra_Exceptions.hpp"
00062
00063 namespace Xpetra {
00064
00065 template <class Scalar, class LocalOrdinal = int, class GlobalOrdinal = LocalOrdinal, class Node = Kokkos::DefaultNode::DefaultNodeType, class LocalMatOps = typename Kokkos::DefaultKernels<Scalar,LocalOrdinal,Node>::SparseOps>
00066 class CrsMatrixFactory {
00067
00068 private:
00070 CrsMatrixFactory() {}
00071
00072 public:
00073
00075 static RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps> > Build(const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &rowMap, size_t maxNumEntriesPerRow, Xpetra::ProfileType pftype = Xpetra::DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null) {
00076 XPETRA_MONITOR("CrsMatrixFactory::Build");
00077
00078 #ifdef HAVE_XPETRA_TPETRA
00079 if (rowMap->lib() == UseTpetra)
00080 return rcp( new TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps>(rowMap, maxNumEntriesPerRow, pftype, plist) );
00081 #endif
00082
00083 XPETRA_FACTORY_ERROR_IF_EPETRA(rowMap->lib());
00084 XPETRA_FACTORY_END;
00085 }
00086
00088 static RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps> > Build(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype = Xpetra::DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null) {
00089
00090 #ifdef HAVE_XPETRA_TPETRA
00091 if (rowMap->lib() == UseTpetra)
00092 return rcp( new TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps>(rowMap, NumEntriesPerRowToAlloc, pftype, plist) );
00093 #endif
00094
00095 XPETRA_FACTORY_ERROR_IF_EPETRA(rowMap->lib());
00096 XPETRA_FACTORY_END;
00097 }
00098
00100 RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps> > Build(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null) {
00101 XPETRA_MONITOR("CrsMatrixFactory::Build");
00102
00103 #ifdef HAVE_XPETRA_TPETRA
00104 if (rowMap->lib() == UseTpetra)
00105 return rcp( new TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps>(rowMap, colMap, maxNumEntriesPerRow, pftype, plist) );
00106 #endif
00107
00108 XPETRA_FACTORY_ERROR_IF_EPETRA(rowMap->lib());
00109 XPETRA_FACTORY_END;
00110 }
00111
00113 RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps> > Build(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null) {
00114 XPETRA_MONITOR("CrsMatrixFactory::Build");
00115
00116 #ifdef HAVE_XPETRA_TPETRA
00117 if (rowMap->lib() == UseTpetra)
00118 return rcp( new TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps>(rowMap, colMap, NumEntriesPerRowToAlloc, pftype, plist) );
00119 #endif
00120
00121 XPETRA_FACTORY_ERROR_IF_EPETRA(rowMap->lib());
00122 XPETRA_FACTORY_END;
00123 }
00124
00126 RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps> > Build(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node, LocalMatOps > > &graph, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null) {
00127 XPETRA_MONITOR("CrsMatrixFactory::Build");
00128
00129 #ifdef HAVE_XPETRA_TPETRA
00130 if (graph->getRowMap()->lib() == UseTpetra)
00131 return rcp( new TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps>(graph, plist) );
00132 #endif
00133
00134 XPETRA_FACTORY_ERROR_IF_EPETRA(graph->getRowMap()->lib());
00135 XPETRA_FACTORY_END;
00136 }
00137
00138 };
00139
00140 template <>
00141 class CrsMatrixFactory<double, int, int> {
00142
00143 typedef double Scalar;
00144 typedef int LocalOrdinal;
00145 typedef int GlobalOrdinal;
00146 typedef Kokkos::DefaultNode::DefaultNodeType Node;
00147 typedef Kokkos::DefaultKernels<Scalar,LocalOrdinal,Node>::SparseOps LocalMatOps;
00148
00149 private:
00151 CrsMatrixFactory() {}
00152
00153 public:
00154
00155 static RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps> > Build(const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &rowMap, size_t maxNumEntriesPerRow, Xpetra::ProfileType pftype = Xpetra::DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null) {
00156 XPETRA_MONITOR("CrsMatrixFactory::Build");
00157
00158 #ifdef HAVE_XPETRA_TPETRA
00159 if (rowMap->lib() == UseTpetra)
00160 return rcp( new TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps>(rowMap, maxNumEntriesPerRow, pftype, plist) );
00161 #endif
00162
00163 #ifdef HAVE_XPETRA_EPETRA
00164 if (rowMap->lib() == UseEpetra)
00165 return rcp( new EpetraCrsMatrix(rowMap, maxNumEntriesPerRow, pftype, plist) );
00166 #endif
00167
00168 XPETRA_FACTORY_END;
00169 }
00170
00171 static RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps> > Build(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype = Xpetra::DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null) {
00172 XPETRA_MONITOR("CrsMatrixFactory::Build");
00173
00174 #ifdef HAVE_XPETRA_TPETRA
00175 if (rowMap->lib() == UseTpetra)
00176 return rcp( new TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps>(rowMap, NumEntriesPerRowToAlloc, pftype, plist) );
00177 #endif
00178
00179 #ifdef HAVE_XPETRA_EPETRA
00180 if (rowMap->lib() == UseEpetra)
00181 return rcp( new EpetraCrsMatrix(rowMap, NumEntriesPerRowToAlloc, pftype, plist) );
00182 #endif
00183
00184 XPETRA_FACTORY_END;
00185 }
00186
00188 RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps> > Build(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null) {
00189 XPETRA_MONITOR("CrsMatrixFactory::Build");
00190
00191 #ifdef HAVE_XPETRA_TPETRA
00192 if (rowMap->lib() == UseTpetra)
00193 return rcp( new TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps>(rowMap, colMap, maxNumEntriesPerRow, pftype, plist) );
00194 #endif
00195
00196 #ifdef HAVE_XPETRA_EPETRA
00197 if (rowMap->lib() == UseEpetra)
00198 return rcp( new EpetraCrsMatrix(rowMap, colMap, maxNumEntriesPerRow, pftype, plist) );
00199 #endif
00200
00201 XPETRA_FACTORY_END;
00202 }
00203
00205 RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps> > Build(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null) {
00206 XPETRA_MONITOR("CrsMatrixFactory::Build");
00207
00208 #ifdef HAVE_XPETRA_TPETRA
00209 if (rowMap->lib() == UseTpetra)
00210 return rcp( new TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps>(rowMap, colMap, NumEntriesPerRowToAlloc, pftype, plist) );
00211 #endif
00212
00213 #ifdef HAVE_XPETRA_EPETRA
00214 if (rowMap->lib() == UseEpetra)
00215 return rcp( new EpetraCrsMatrix(rowMap, colMap, NumEntriesPerRowToAlloc, pftype, plist) );
00216 #endif
00217
00218 XPETRA_FACTORY_END;
00219 }
00220
00222 RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps> > Build(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node, LocalMatOps > > &graph, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null) {
00223 XPETRA_MONITOR("CrsMatrixFactory::Build");
00224
00225 #ifdef HAVE_XPETRA_TPETRA
00226 if (graph->getRowMap()->lib() == UseTpetra)
00227 return rcp( new TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps>(graph, plist) );
00228 #endif
00229
00230 #ifdef HAVE_XPETRA_EPETRA
00231 if (graph->getRowMap()->lib() == UseEpetra)
00232 return rcp( new EpetraCrsMatrix(graph, plist) );
00233 #endif
00234
00235 XPETRA_FACTORY_END;
00236 }
00237 };
00238
00239 }
00240
00241 #define XPETRA_CRSMATRIXFACTORY_SHORT
00242 #endif