Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_FECrsMatrix_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Tpetra: Templated Linear Algebra Services Package
4//
5// Copyright 2008 NTESS and the Tpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef TPETRA_FECRSMATRIX_DECL_HPP
11#define TPETRA_FECRSMATRIX_DECL_HPP
12
15
16#include "Tpetra_ConfigDefs.hpp"
18#include "Tpetra_FECrsGraph.hpp"
19
20namespace Tpetra {
21
25template <class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type,
27 class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type,
28 class Node = ::Tpetra::Details::DefaultTypes::node_type>
29class FECrsMatrix : public CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
30 friend class CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
31
32 public:
34
35
39
50
53
56
59
62
65
72
75
78
81
84
87
91
94 typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_host_type local_matrix_host_type;
95
98
101
103
105
136 explicit FECrsMatrix(const Teuchos::RCP<const fe_crs_graph_type>& graph,
137 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
138
141
144
151
161 virtual ~FECrsMatrix() = default;
162
164
165
194 void globalAssemble() { endFill(); }
195
197 void endAssembly();
198
202 void beginAssembly();
203
205 void endModify();
206
213 void beginModify();
214
215 private:
217 void endFill();
218 void beginFill();
219
224 static const bool useAtomicUpdatesByDefault =
225#ifdef KOKKOS_ENABLE_SERIAL
226 !std::is_same<execution_space, Kokkos::Serial>::value;
227#else
228 true;
229#endif // KOKKOS_ENABLE_SERIAL
230
231 protected:
235 const crs_graph_type& graph,
236 const RowInfo& rowInfo,
237 const GlobalOrdinal inds[],
239 const LocalOrdinal numElts);
240
243 const crs_graph_type& graph,
244 const RowInfo& rowInfo,
245 const LocalOrdinal inds[],
247 const LocalOrdinal numElts);
248
251 const crs_graph_type& graph,
252 const RowInfo& rowInfo,
253 const GlobalOrdinal inds[],
255 const LocalOrdinal numElts,
256 const bool atomic = useAtomicUpdatesByDefault);
257
260 const crs_graph_type& graph,
261 const RowInfo& rowInfo,
262 const LocalOrdinal inds[],
264 const LocalOrdinal numElts,
265 const bool atomic = useAtomicUpdatesByDefault);
266
267 void
271 const impl_scalar_type vals[],
272 const size_t numInputEnt);
273
274 protected:
279
283
287
288 private:
289 // The FECrsGraph from construction time
290 Teuchos::RCP<const FECrsGraph<LocalOrdinal, GlobalOrdinal, Node> > feGraph_;
291
292 // This is whichever multivector isn't currently active
293 Teuchos::RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > inactiveCrsMatrix_;
294 // This is in RCP to make shallow copies of the FECrsMatrix work correctly
295 Teuchos::RCP<FE::WhichActive> activeCrsMatrix_;
296
297 Teuchos::RCP<FE::FillState> fillState_;
298
299}; // end class FECrsMatrix
300
301} // end namespace Tpetra
302
303#endif // TPETRA_FECRSMATRIX_DECL_HPP
Declaration of the Tpetra::CrsMatrix class.
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
Struct that holds views of the contents of a CrsMatrix.
Sparse matrix that presents a row-oriented interface that lets users read or modify entries.
typename device_type::execution_space execution_space
The Kokkos execution space.
Node node_type
This class' Kokkos Node type.
typename Node::device_type device_type
The Kokkos device type.
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...
typename Kokkos::ArithTraits< impl_scalar_type >::mag_type mag_type
Type of a norm result.
typename crs_graph_type::local_graph_device_type local_graph_device_type
The part of the sparse matrix's graph on each MPI process.
typename row_matrix_type::impl_scalar_type impl_scalar_type
The type used internally in place of Scalar.
Communication plan for data redistribution from a (possibly) multiply-owned to a uniquely-owned distr...
Sparse matrix that presents a row-oriented interface that lets users read or modify entries.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > crs_matrix_type
Parent CrsMatrix type using the same scalars.
FECrsMatrix(const FECrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &)=delete
Copy constructor (forbidden).
void beginAssembly()
Activates the owned+shared mode for assembly.
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.
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.
FECrsGraph< LocalOrdinal, GlobalOrdinal, Node > fe_crs_graph_type
The CrsGraph specialization suitable for this CrsMatrix specialization.
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.
FECrsMatrix(FECrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &&)=delete
Move constructor (forbidden).
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.
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.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::node_type node_type
This class' fourth template parameter; the Kokkos device type.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::impl_scalar_type impl_scalar_type
The type used internally in place of Scalar.
void endAssembly()
Migrates data to the owned mode.
void doOwnedPlusSharedToOwned(const CombineMode CM=Tpetra::ADD)
Migrate data from the owned+shared to the owned matrix Since this is non-unique -> unique,...
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::execution_space execution_space
The Kokkos execution space.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mag_type mag_type
Type of a norm result.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::map_type map_type
The Map specialization suitable for this CrsMatrix specialization.
FECrsMatrix & operator=(FECrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &&)=delete
Move assignment (forbidden).
void endModify()
Closes modification phase.
virtual ~FECrsMatrix()=default
Destructor (virtual for memory safety of derived classes).
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::export_type export_type
The Export specialization suitable for this CrsMatrix specialization.
void globalAssemble()
Communicate nonlocal contributions to other processes.
void beginModify()
Activates the owned mode for modifying local values.
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...
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.
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.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::global_ordinal_type global_ordinal_type
This class' third template parameter; the type of global indices.
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 >::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...
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::device_type device_type
The Kokkos device type.
FECrsMatrix & operator=(const FECrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &)=delete
Copy assignment (forbidden).
void switchActiveCrsMatrix()
Switches which CrsGraph is active (without migrating data)
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type local_ordinal_type
This class' second template parameter; the type of local indices.
void doOwnedToOwnedPlusShared(const CombineMode CM=Tpetra::ADD)
Migrate data from the owned to the owned+shared matrix Precondition: Must be FE_ACTIVE_OWNED mode.
Communication plan for data redistribution from a uniquely-owned to a (possibly) multiply-owned distr...
A parallel distribution of indices over processes.
int local_ordinal_type
Default value of Scalar template parameter.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
CombineMode
Rule for combining data in an Import or Export.
@ ADD
Sum new values.
Allocation information for a locally owned row in a CrsGraph or CrsMatrix.