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 //
4 // Tpetra: Templated Linear Algebra Services Package
5 // Copyright (2008) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 // @HEADER
41 #ifndef TPETRA_FECRSMATRIX_DECL_HPP
42 #define TPETRA_FECRSMATRIX_DECL_HPP
43 
44 
47 
49 #include "Tpetra_FECrsGraph.hpp"
50 
51 namespace Tpetra {
52 
53 
54 
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>
62 class FECrsMatrix :
63  public CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>
64 {
65  friend class CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
66 public:
68 
69 
73 
84 
87 
90 
93 
96 
99 
106 
109 
112 
115 
118 
121 
125 
128  typedef typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_host_type local_matrix_host_type;
129 
132 
135 
137 
139 
140 
171  explicit FECrsMatrix (const Teuchos::RCP<const fe_crs_graph_type>& graph,
172  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
173 
176 
179 
181  FECrsMatrix&
184  FECrsMatrix&
186 
196  virtual ~FECrsMatrix () = default;
197 
199 
200 
229  void globalAssemble () {endFill();}
230 
232  void endAssembly();
233 
237  void beginAssembly();
238 
240  void endModify();
241 
248  void beginModify();
249 
250  private:
251 
253  void endFill();
254  void beginFill();
255 
260  static const bool useAtomicUpdatesByDefault =
261 #ifdef KOKKOS_ENABLE_SERIAL
262  ! std::is_same<execution_space, Kokkos::Serial>::value;
263 #else
264  true;
265 #endif // KOKKOS_ENABLE_SERIAL
266 
267  protected:
268 
270  LocalOrdinal
272  const crs_graph_type& graph,
273  const RowInfo& rowInfo,
274  const GlobalOrdinal inds[],
275  const impl_scalar_type newVals[],
276  const LocalOrdinal numElts);
277 
278  LocalOrdinal
279  replaceLocalValuesImpl (impl_scalar_type rowVals[],
280  const crs_graph_type& graph,
281  const RowInfo& rowInfo,
282  const LocalOrdinal inds[],
283  const impl_scalar_type newVals[],
284  const LocalOrdinal numElts);
285 
286  LocalOrdinal
287  sumIntoGlobalValuesImpl (impl_scalar_type rowVals[],
288  const crs_graph_type& graph,
289  const RowInfo& rowInfo,
290  const GlobalOrdinal inds[],
291  const impl_scalar_type newVals[],
292  const LocalOrdinal numElts,
293  const bool atomic = useAtomicUpdatesByDefault);
294 
295  LocalOrdinal
296  sumIntoLocalValuesImpl (impl_scalar_type rowVals[],
297  const crs_graph_type& graph,
298  const RowInfo& rowInfo,
299  const LocalOrdinal inds[],
300  const impl_scalar_type newVals[],
301  const LocalOrdinal numElts,
302  const bool atomic = useAtomicUpdatesByDefault);
303 
304  void
305  insertGlobalValuesImpl (crs_graph_type& graph,
306  RowInfo& rowInfo,
307  const GlobalOrdinal gblColInds[],
308  const impl_scalar_type vals[],
309  const size_t numInputEnt);
310 
311  protected:
316 
320 
322  void switchActiveCrsMatrix();
324 
325  private:
326  // Enum for activity
327  enum FEWhichActive
328  {
329  FE_ACTIVE_OWNED,
330  FE_ACTIVE_OWNED_PLUS_SHARED
331  };
332 
333  // The FECrsGraph from construction time
334  Teuchos::RCP<const FECrsGraph<LocalOrdinal, GlobalOrdinal, Node> > feGraph_;
335 
336  // This is whichever multivector isn't currently active
337  Teuchos::RCP<CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > inactiveCrsMatrix_;
338  // This is in RCP to make shallow copies of the FECrsMatrix work correctly
339  Teuchos::RCP<FEWhichActive> activeCrsMatrix_;
340 
341  enum class FillState
342  {
343  open, // matrix is "open". Values can freely summed in to and replaced
344  modify, // matrix is open for modification. *local* values can be replaced
345  closed
346  };
347  Teuchos::RCP<FillState> fillState_;
348 
349 }; // end class FECrsMatrix
350 
351 
352 
353 } // end namespace Tpetra
354 
355 
356 
357 #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...