Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 
13 
16 
18 #include "Tpetra_FECrsGraph.hpp"
19 
20 namespace Tpetra {
21 
22 
23 
27 template<class Scalar = ::Tpetra::Details::DefaultTypes::scalar_type,
29  class GlobalOrdinal = ::Tpetra::Details::DefaultTypes::global_ordinal_type,
30  class Node = ::Tpetra::Details::DefaultTypes::node_type>
31 class FECrsMatrix :
32  public CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>
33 {
34  friend class CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
35 public:
37 
38 
42 
53 
56 
59 
62 
65 
68 
75 
78 
81 
84 
87 
90 
94 
97  typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_host_type local_matrix_host_type;
98 
101 
104 
106 
108 
109 
140  explicit FECrsMatrix (const Teuchos::RCP<const fe_crs_graph_type>& graph,
141  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
142 
145 
148 
150  FECrsMatrix&
153  FECrsMatrix&
155 
165  virtual ~FECrsMatrix () = default;
166 
168 
169 
198  void globalAssemble () {endFill();}
199 
201  void endAssembly();
202 
206  void beginAssembly();
207 
209  void endModify();
210 
217  void beginModify();
218 
219  private:
220 
222  void endFill();
223  void beginFill();
224 
229  static const bool useAtomicUpdatesByDefault =
230 #ifdef KOKKOS_ENABLE_SERIAL
231  ! std::is_same<execution_space, Kokkos::Serial>::value;
232 #else
233  true;
234 #endif // KOKKOS_ENABLE_SERIAL
235 
236  protected:
237 
239  LocalOrdinal
241  const crs_graph_type& graph,
242  const RowInfo& rowInfo,
243  const GlobalOrdinal inds[],
244  const impl_scalar_type newVals[],
245  const LocalOrdinal numElts);
246 
247  LocalOrdinal
248  replaceLocalValuesImpl (impl_scalar_type rowVals[],
249  const crs_graph_type& graph,
250  const RowInfo& rowInfo,
251  const LocalOrdinal inds[],
252  const impl_scalar_type newVals[],
253  const LocalOrdinal numElts);
254 
255  LocalOrdinal
256  sumIntoGlobalValuesImpl (impl_scalar_type rowVals[],
257  const crs_graph_type& graph,
258  const RowInfo& rowInfo,
259  const GlobalOrdinal inds[],
260  const impl_scalar_type newVals[],
261  const LocalOrdinal numElts,
262  const bool atomic = useAtomicUpdatesByDefault);
263 
264  LocalOrdinal
265  sumIntoLocalValuesImpl (impl_scalar_type rowVals[],
266  const crs_graph_type& graph,
267  const RowInfo& rowInfo,
268  const LocalOrdinal inds[],
269  const impl_scalar_type newVals[],
270  const LocalOrdinal numElts,
271  const bool atomic = useAtomicUpdatesByDefault);
272 
273  void
274  insertGlobalValuesImpl (crs_graph_type& graph,
275  RowInfo& rowInfo,
276  const GlobalOrdinal gblColInds[],
277  const impl_scalar_type vals[],
278  const size_t numInputEnt);
279 
280  protected:
285 
289 
291  void switchActiveCrsMatrix();
293 
294  private:
295  // Enum for activity
296  enum FEWhichActive
297  {
298  FE_ACTIVE_OWNED,
299  FE_ACTIVE_OWNED_PLUS_SHARED
300  };
301 
302  // The FECrsGraph from construction time
303  Teuchos::RCP<const FECrsGraph<LocalOrdinal, GlobalOrdinal, Node> > feGraph_;
304 
305  // This is whichever multivector isn't currently active
306  Teuchos::RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > inactiveCrsMatrix_;
307  // This is in RCP to make shallow copies of the FECrsMatrix work correctly
308  Teuchos::RCP<FEWhichActive> activeCrsMatrix_;
309 
310  enum class FillState
311  {
312  open, // matrix is "open". Values can freely summed in to and replaced
313  modify, // matrix is open for modification. *local* values can be replaced
314  closed
315  };
316  Teuchos::RCP<FillState> fillState_;
317 
318 }; // end class FECrsMatrix
319 
320 
321 
322 } // end namespace Tpetra
323 
324 
325 
326 #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&#39; 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&#39; 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&#39; Kokkos Node type.
CombineMode
Rule for combining data in an Import or Export.
void beginAssembly()
Activates the owned+shared mode for assembly.
Sum new values.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::node_type node_type
This class&#39; 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 -&gt; unique...
FECrsMatrix(const Teuchos::RCP< const fe_crs_graph_type > &graph, const Teuchos::RCP< Teuchos::ParameterList > &params=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&#39;s graph on each MPI process.
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::global_ordinal_type global_ordinal_type
This class&#39; 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&#39;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...