Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ifpack2_Experimental_RBILUK_decl.hpp
Go to the documentation of this file.
1 /*@HEADER
2 // ***********************************************************************
3 //
4 // Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
5 // Copyright (2009) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
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 */
42 
45 
46 #ifndef IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP
47 #define IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP
48 
49 #include <Tpetra_BlockCrsMatrix.hpp>
50 
51 #include <Ifpack2_RILUK.hpp>
52 
53 namespace Ifpack2 {
54 
55 namespace Experimental {
56 
127 template<class MatrixType>
128 class RBILUK : virtual public Ifpack2::RILUK< Tpetra::RowMatrix< typename MatrixType::scalar_type,
129  typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type> >
130 {
131  public:
132 
134 
135  typedef typename MatrixType::scalar_type scalar_type;
137 
138  //typedef typename MatrixType::impl_scalar_type impl_scalar_type;
139  typedef typename Kokkos::ArithTraits<typename MatrixType::scalar_type>::val_type impl_scalar_type;
140 
142  typedef typename MatrixType::local_ordinal_type local_ordinal_type;
143  typedef typename MatrixType::local_ordinal_type LO;
144 
146  typedef typename MatrixType::global_ordinal_type global_ordinal_type;
147  typedef typename MatrixType::global_ordinal_type GO;
148 
150  typedef typename MatrixType::node_type node_type;
151 
154 
156  typedef Tpetra::RowMatrix<scalar_type,
160 
162  typedef Tpetra::CrsMatrix<scalar_type,
166 
167  using crs_graph_type = Tpetra::CrsGraph<local_ordinal_type,
169  node_type>;
170 
171  typedef Tpetra::BlockCrsMatrix<scalar_type,
174  node_type> block_crs_matrix_type;
175 
176  template <class NewMatrixType> friend class RBILUK;
177 
179 
181 
182  typedef typename crs_matrix_type::local_matrix_device_type local_matrix_device_type;
183  typedef typename local_matrix_device_type::StaticCrsGraphType::row_map_type lno_row_view_t;
184  typedef typename local_matrix_device_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
185  typedef typename local_matrix_device_type::values_type scalar_nonzero_view_t;
186  typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space TemporaryMemorySpace;
187  typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::memory_space PersistentMemorySpace;
188  typedef typename local_matrix_device_type::StaticCrsGraphType::device_type::execution_space HandleExecSpace;
189  typedef typename KokkosKernels::Experimental::KokkosKernelsHandle
190  <typename lno_row_view_t::const_value_type, typename lno_nonzero_view_t::const_value_type, typename scalar_nonzero_view_t::value_type,
191  HandleExecSpace, TemporaryMemorySpace,PersistentMemorySpace > kk_handle_type;
192  //typedef typename KokkosKernels::Experimental::KokkosKernelsHandle
193  // <typename lno_row_view_t::non_const_value_type, typename lno_nonzero_view_t::non_const_value_type, typename scalar_nonzero_view_t::value_type,
194  // HandleExecSpace, TemporaryMemorySpace,PersistentMemorySpace > kk_handle_type;//test
195  Teuchos::RCP<kk_handle_type> KernelHandle_;
196 
198 
200 
205 
210 
211  private:
214  RBILUK (const RBILUK<MatrixType> & src);
215 
216  public:
217 
219  virtual ~RBILUK ();
221 
223  void initialize ();
224 
233  void compute ();
234 
236 
237 
238  // Declare that we intend to overload RILUK::setMatrix, not hide it.
239  // This avoids build warnings that the method below "hides
240  // overloaded virtual function" (e.g., Clang 3.5).
241  //
242  // NOTE: If the base class of this class changes, e.g., if its
243  // template parameter changes, then be sure to change the code below
244  // to refer to the proper base class.
245  using RILUK<Tpetra::RowMatrix<typename MatrixType::scalar_type,
246  typename MatrixType::local_ordinal_type,
247  typename MatrixType::global_ordinal_type,
248  typename MatrixType::node_type> >::setMatrix;
249 
272  void
274 
276 
278 
280  std::string description () const;
281 
283 
285 
315  void
316  apply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
317  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
322 
323 public:
324 
327 
329  const block_crs_matrix_type& getLBlock () const;
330 
332  const block_crs_matrix_type& getDBlock () const;
333 
335  const block_crs_matrix_type& getUBlock () const;
336 
337 private:
338  typedef Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> MV;
341  typedef typename block_crs_matrix_type::little_block_type little_block_type;
342  typedef typename block_crs_matrix_type::little_block_host_type little_block_host_type;
343  typedef typename block_crs_matrix_type::little_vec_type little_vec_type;
344  typedef typename block_crs_matrix_type::little_host_vec_type little_host_vec_type;
345  typedef typename block_crs_matrix_type::const_host_little_vec_type const_host_little_vec_type;
346 
347  using local_inds_host_view_type = typename block_crs_matrix_type::local_inds_host_view_type;
348  using values_host_view_type = typename block_crs_matrix_type::values_host_view_type;
349  using local_inds_device_view_type = typename block_crs_matrix_type::local_inds_device_view_type;
350  using values_device_view_type = typename block_crs_matrix_type::values_device_view_type;
351 
352  void allocate_L_and_U_blocks();
353  void initAllValues ();
354 
356  local_ordinal_type blockSize_;
357 
364 
366  Teuchos::RCP<block_crs_matrix_type> D_block_inverse_;
367 };
368 
369 
370 } // namepsace Experimental
371 
372 } // namespace Ifpack2
373 
374 #endif /* IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP */
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:344
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:146
virtual ~RBILUK()
Destructor (declared virtual for memory safety).
Definition: Ifpack2_Experimental_RBILUK_def.hpp:133
const block_crs_matrix_type & getUBlock() const
Return the U factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:190
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:136
void apply(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &X, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const
Apply the (inverse of the) incomplete factorization to X, resulting in Y.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:920
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:153
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition: Ifpack2_RILUK_decl.hpp:245
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:150
Tpetra::CrsMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > crs_matrix_type
Tpetra::CrsMatrix specialization used by this class for representing L and U.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:165
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:142
const block_crs_matrix_type & getDBlock() const
Return the diagonal entries of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:176
void compute()
Compute the (numeric) incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:578
Teuchos::RCP< const block_crs_matrix_type > getBlockMatrix() const
Get the input matrix.
void setMatrix(const Teuchos::RCP< const block_crs_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:138
const block_crs_matrix_type & getLBlock() const
Return the L factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:162
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Tpetra::RowMatrix specialization used by this class.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:159
std::string description() const
A one-line description of this object.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:1072
ILU(k) factorization of a given Tpetra::BlockCrsMatrix.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:128