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 
17 #include "Tpetra_ConfigDefs.hpp"
19 #include "Tpetra_FECrsGraph.hpp"
20 
21 namespace Tpetra {
22 
23 
24 
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>
32 class FECrsMatrix :
33  public CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>
34 {
35  friend class CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
36 public:
38 
39 
43 
54 
57 
60 
63 
66 
69 
76 
79 
82 
85 
88 
91 
95 
98  typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_host_type local_matrix_host_type;
99 
102 
105 
107 
109 
110 
141  explicit FECrsMatrix (const Teuchos::RCP<const fe_crs_graph_type>& graph,
142  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
143 
146 
149 
151  FECrsMatrix&
154  FECrsMatrix&
156 
166  virtual ~FECrsMatrix () = default;
167 
169 
170 
199  void globalAssemble () {endFill();}
200 
202  void endAssembly();
203 
207  void beginAssembly();
208 
210  void endModify();
211 
218  void beginModify();
219 
220  private:
221 
223  void endFill();
224  void beginFill();
225 
230  static const bool useAtomicUpdatesByDefault =
231 #ifdef KOKKOS_ENABLE_SERIAL
232  ! std::is_same<execution_space, Kokkos::Serial>::value;
233 #else
234  true;
235 #endif // KOKKOS_ENABLE_SERIAL
236 
237  protected:
238 
240  LocalOrdinal
242  const crs_graph_type& graph,
243  const RowInfo& rowInfo,
244  const GlobalOrdinal inds[],
245  const impl_scalar_type newVals[],
246  const LocalOrdinal numElts);
247 
248  LocalOrdinal
249  replaceLocalValuesImpl (impl_scalar_type rowVals[],
250  const crs_graph_type& graph,
251  const RowInfo& rowInfo,
252  const LocalOrdinal inds[],
253  const impl_scalar_type newVals[],
254  const LocalOrdinal numElts);
255 
256  LocalOrdinal
257  sumIntoGlobalValuesImpl (impl_scalar_type rowVals[],
258  const crs_graph_type& graph,
259  const RowInfo& rowInfo,
260  const GlobalOrdinal inds[],
261  const impl_scalar_type newVals[],
262  const LocalOrdinal numElts,
263  const bool atomic = useAtomicUpdatesByDefault);
264 
265  LocalOrdinal
266  sumIntoLocalValuesImpl (impl_scalar_type rowVals[],
267  const crs_graph_type& graph,
268  const RowInfo& rowInfo,
269  const LocalOrdinal inds[],
270  const impl_scalar_type newVals[],
271  const LocalOrdinal numElts,
272  const bool atomic = useAtomicUpdatesByDefault);
273 
274  void
275  insertGlobalValuesImpl (crs_graph_type& graph,
276  RowInfo& rowInfo,
277  const GlobalOrdinal gblColInds[],
278  const impl_scalar_type vals[],
279  const size_t numInputEnt);
280 
281  protected:
286 
290 
292  void switchActiveCrsMatrix();
294 
295  private:
296 
297  // The FECrsGraph from construction time
298  Teuchos::RCP<const FECrsGraph<LocalOrdinal, GlobalOrdinal, Node> > feGraph_;
299 
300  // This is whichever multivector isn't currently active
301  Teuchos::RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > inactiveCrsMatrix_;
302  // This is in RCP to make shallow copies of the FECrsMatrix work correctly
303  Teuchos::RCP<FE::WhichActive> activeCrsMatrix_;
304 
305  Teuchos::RCP<FE::FillState> fillState_;
306 
307 }; // end class FECrsMatrix
308 
309 
310 
311 } // end namespace Tpetra
312 
313 
314 
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&#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...