41 #ifndef TPETRA_FECRSMATRIX_DECL_HPP
42 #define TPETRA_FECRSMATRIX_DECL_HPP
49 #include "Tpetra_FECrsGraph.hpp"
58 template<
class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type,
60 class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type,
61 class Node = ::Tpetra::Details::DefaultTypes::node_type>
63 public CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>
65 friend class CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
72 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::scalar_type scalar_type;
83 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::impl_scalar_type impl_scalar_type;
86 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_ordinal_type
local_ordinal_type;
89 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::global_ordinal_type global_ordinal_type;
92 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::node_type node_type;
95 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::device_type device_type;
98 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::execution_space execution_space;
105 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::mag_type mag_type;
108 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::map_type map_type;
111 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::import_type import_type;
114 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::export_type export_type;
117 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::crs_graph_type crs_graph_type;
120 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_graph_device_type local_graph_device_type;
124 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_device_type local_matrix_device_type;
128 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_host_type local_matrix_host_type;
131 typedef CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> crs_matrix_type;
134 typedef FECrsGraph<LocalOrdinal, GlobalOrdinal, Node> fe_crs_graph_type;
171 explicit FECrsMatrix (
const Teuchos::RCP<const fe_crs_graph_type>& graph,
172 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
175 FECrsMatrix (
const FECrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>&) =
delete;
178 FECrsMatrix (FECrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>&&) =
delete;
182 operator= (
const FECrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>&) =
delete;
185 operator= (FECrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>&&) =
delete;
196 virtual ~FECrsMatrix () =
default;
229 void globalAssemble () {endFill();}
235 void beginAssembly();
253 static const bool useAtomicUpdatesByDefault =
254 #ifdef KOKKOS_ENABLE_SERIAL
255 ! std::is_same<execution_space, Kokkos::Serial>::value;
258 #endif // KOKKOS_ENABLE_SERIAL
265 const crs_graph_type& graph,
266 const RowInfo& rowInfo,
267 const GlobalOrdinal inds[],
268 const impl_scalar_type newVals[],
269 const LocalOrdinal numElts);
273 const crs_graph_type& graph,
274 const RowInfo& rowInfo,
275 const LocalOrdinal inds[],
276 const impl_scalar_type newVals[],
277 const LocalOrdinal numElts);
281 const crs_graph_type& graph,
282 const RowInfo& rowInfo,
283 const GlobalOrdinal inds[],
284 const impl_scalar_type newVals[],
285 const LocalOrdinal numElts,
286 const bool atomic = useAtomicUpdatesByDefault);
290 const crs_graph_type& graph,
291 const RowInfo& rowInfo,
292 const LocalOrdinal inds[],
293 const impl_scalar_type newVals[],
294 const LocalOrdinal numElts,
295 const bool atomic = useAtomicUpdatesByDefault);
300 const GlobalOrdinal gblColInds[],
301 const impl_scalar_type vals[],
302 const size_t numInputEnt);
315 void switchActiveCrsMatrix();
323 FE_ACTIVE_OWNED_PLUS_SHARED
327 Teuchos::RCP<const FECrsGraph<LocalOrdinal, GlobalOrdinal, Node> > feGraph_;
330 Teuchos::RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > inactiveCrsMatrix_;
332 Teuchos::RCP<FEWhichActive> activeCrsMatrix_;
340 Teuchos::RCP<FillState> fillState_;
350 #endif // TPETRA_FECRSMATRIX_DECL_HPP
virtual LocalOrdinal replaceGlobalValuesImpl(impl_scalar_type rowVals[], const crs_graph_type &graph, const RowInfo &rowInfo, const GlobalOrdinal inds[], const impl_scalar_type newVals[], const LocalOrdinal numElts)
Implementation detail of replaceGlobalValues.
virtual void insertGlobalValuesImpl(crs_graph_type &graph, RowInfo &rowInfo, const GlobalOrdinal gblColInds[], const impl_scalar_type vals[], const size_t numInputEnt)
Common implementation detail of insertGlobalValues and insertGlobalValuesFiltered.
Declaration of the Tpetra::CrsMatrix class.
int local_ordinal_type
Default value of Scalar template parameter.
virtual LocalOrdinal sumIntoLocalValuesImpl(impl_scalar_type rowVals[], const crs_graph_type &graph, const RowInfo &rowInfo, const LocalOrdinal inds[], const impl_scalar_type newVals[], const LocalOrdinal numElts, const bool atomic=useAtomicUpdatesByDefault)
Implementation detail of sumIntoLocalValues.
CombineMode
Rule for combining data in an Import or Export.
virtual LocalOrdinal sumIntoGlobalValuesImpl(impl_scalar_type rowVals[], const crs_graph_type &graph, const RowInfo &rowInfo, const GlobalOrdinal inds[], const impl_scalar_type newVals[], const LocalOrdinal numElts, const bool atomic=useAtomicUpdatesByDefault)
Implementation detail of sumIntoGlobalValues.
virtual LocalOrdinal replaceLocalValuesImpl(impl_scalar_type rowVals[], const crs_graph_type &graph, const RowInfo &rowInfo, const LocalOrdinal inds[], const impl_scalar_type newVals[], const LocalOrdinal numElts)
Implementation detail of replaceLocalValues.
CrsMatrix(const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &)=default
Copy constructor.