10 #ifndef TPETRA_FECRSMATRIX_DECL_HPP
11 #define TPETRA_FECRSMATRIX_DECL_HPP
17 #include "Tpetra_ConfigDefs.hpp"
19 #include "Tpetra_FECrsGraph.hpp"
28 template<
class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type,
30 class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type,
31 class Node = ::Tpetra::Details::DefaultTypes::node_type>
33 public CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>
35 friend class CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
98 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_host_type
local_matrix_host_type;
141 explicit FECrsMatrix (
const Teuchos::RCP<const fe_crs_graph_type>& graph,
142 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
230 static const bool useAtomicUpdatesByDefault =
231 #ifdef KOKKOS_ENABLE_SERIAL
232 ! std::is_same<execution_space, Kokkos::Serial>::value;
235 #endif // KOKKOS_ENABLE_SERIAL
244 const GlobalOrdinal inds[],
246 const LocalOrdinal numElts);
252 const LocalOrdinal inds[],
254 const LocalOrdinal numElts);
260 const GlobalOrdinal inds[],
262 const LocalOrdinal numElts,
263 const bool atomic = useAtomicUpdatesByDefault);
269 const LocalOrdinal inds[],
271 const LocalOrdinal numElts,
272 const bool atomic = useAtomicUpdatesByDefault);
277 const GlobalOrdinal gblColInds[],
279 const size_t numInputEnt);
298 Teuchos::RCP<const FECrsGraph<LocalOrdinal, GlobalOrdinal, Node> > feGraph_;
301 Teuchos::RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > inactiveCrsMatrix_;
303 Teuchos::RCP<FE::WhichActive> activeCrsMatrix_;
305 Teuchos::RCP<FE::FillState> fillState_;
315 #endif // TPETRA_FECRSMATRIX_DECL_HPP
Communication plan for data redistribution from a uniquely-owned to a (possibly) multiply-owned distr...
Sparse matrix that presents a row-oriented interface that lets users read or modify entries...
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
Sparse matrix that presents a row-oriented interface that lets users read or modify entries...
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::device_type device_type
The Kokkos device type.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type local_ordinal_type
This class' second template parameter; the type of local indices.
typename device_type::execution_space execution_space
The Kokkos execution space.
typename Kokkos::ArithTraits< impl_scalar_type >::mag_type mag_type
Type of a norm result.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > crs_matrix_type
Parent CrsMatrix type using the same scalars.
Declaration of the Tpetra::CrsMatrix class.
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)
Overloads of modification methods.
GlobalOrdinal global_ordinal_type
The type of each global index in the matrix.
FECrsGraph< LocalOrdinal, GlobalOrdinal, Node > fe_crs_graph_type
The CrsGraph specialization suitable for this CrsMatrix specialization.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scalar_type scalar_type
This class' first template parameter; the type of each entry in the matrix.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::map_type map_type
The Map specialization suitable for this CrsMatrix specialization.
Scalar scalar_type
The type of each entry in the matrix.
Allocation information for a locally owned row in a CrsGraph or CrsMatrix.
int local_ordinal_type
Default value of Scalar template parameter.
void doOwnedToOwnedPlusShared(const CombineMode CM=Tpetra::ADD)
Migrate data from the owned to the owned+shared matrix Precondition: Must be FE_ACTIVE_OWNED mode...
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::impl_scalar_type impl_scalar_type
The type used internally in place of Scalar.
KokkosSparse::CrsMatrix< impl_scalar_type, local_ordinal_type, device_type, void, typename local_graph_device_type::size_type > local_matrix_device_type
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI pro...
void globalAssemble()
Communicate nonlocal contributions to other processes.
virtual ~FECrsMatrix()=default
Destructor (virtual for memory safety of derived classes).
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::crs_graph_type crs_graph_type
The CrsGraph specialization suitable for this CrsMatrix specialization.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::import_type import_type
The Import specialization suitable for this CrsMatrix specialization.
Communication plan for data redistribution from a (possibly) multiply-owned to a uniquely-owned distr...
void switchActiveCrsMatrix()
Switches which CrsGraph is active (without migrating data)
Node node_type
This class' Kokkos Node type.
CombineMode
Rule for combining data in an Import or Export.
void beginAssembly()
Activates the owned+shared mode for assembly.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::node_type node_type
This class' fourth template parameter; the Kokkos device type.
FECrsMatrix & operator=(const FECrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &)=delete
Copy assignment (forbidden).
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_device_type local_matrix_device_type
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix for each MPI pr...
void endModify()
Closes modification phase.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mag_type mag_type
Type of a norm result.
typename Node::device_type device_type
The Kokkos device type.
void doOwnedPlusSharedToOwned(const CombineMode CM=Tpetra::ADD)
Migrate data from the owned+shared to the owned matrix Since this is non-unique -> unique...
FECrsMatrix(const Teuchos::RCP< const fe_crs_graph_type > &graph, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor specifying one or two previously constructed graphs.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::export_type export_type
The Export specialization suitable for this CrsMatrix specialization.
void endAssembly()
Migrates data to the owned mode.
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
typename row_matrix_type::impl_scalar_type impl_scalar_type
The type used internally in place of Scalar.
A parallel distribution of indices over processes.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::execution_space execution_space
The Kokkos execution space.
void beginModify()
Activates the owned mode for modifying local values.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_graph_device_type local_graph_device_type
The part of the sparse matrix's graph on each MPI process.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::global_ordinal_type global_ordinal_type
This class' third template parameter; the type of global indices.
LocalOrdinal local_ordinal_type
The type of each local index in the matrix.
typename crs_graph_type::local_graph_device_type local_graph_device_type
The part of the sparse matrix's graph on each MPI process.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_host_type local_matrix_host_type
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix for each MPI pr...