Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ifpack2_BlockTriDiContainer_def.hpp
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 
43 #ifndef IFPACK2_BLOCKTRIDICONTAINER_DEF_HPP
44 #define IFPACK2_BLOCKTRIDICONTAINER_DEF_HPP
45 
47 
48 #include <Tpetra_Distributor.hpp>
49 #include <Tpetra_BlockMultiVector.hpp>
50 
51 #include <Kokkos_ArithTraits.hpp>
52 #include <KokkosBatched_Util.hpp>
53 #include <KokkosBatched_Vector.hpp>
54 #include <KokkosBatched_AddRadial_Decl.hpp>
55 #include <KokkosBatched_AddRadial_Impl.hpp>
56 #include <KokkosBatched_Gemm_Decl.hpp>
57 #include <KokkosBatched_Gemm_Serial_Impl.hpp>
58 #include <KokkosBatched_Gemv_Decl.hpp>
59 #include <KokkosBatched_Gemv_Serial_Impl.hpp>
60 #include <KokkosBatched_Trsm_Decl.hpp>
61 #include <KokkosBatched_Trsm_Serial_Impl.hpp>
62 #include <KokkosBatched_Trsv_Decl.hpp>
63 #include <KokkosBatched_Trsv_Serial_Impl.hpp>
64 #include <KokkosBatched_LU_Decl.hpp>
65 #include <KokkosBatched_LU_Serial_Impl.hpp>
66 
68 #include "Ifpack2_BlockTriDiContainer_impl.hpp"
69 
70 #include <memory>
71 
72 
73 namespace Ifpack2 {
74 
78 
79  template <typename MatrixType>
80  void
81  BlockTriDiContainer<MatrixType, BlockTriDiContainerDetails::ImplSimdTag>
82  ::initInternal (const Teuchos::RCP<const row_matrix_type>& matrix,
84  const Teuchos::RCP<const import_type>& importer,
85  const bool overlapCommAndComp,
86  const bool useSeqMethod)
87  {
88  // create pointer of impl
89  impl_ = Teuchos::rcp(new BlockTriDiContainerDetails::ImplObject<MatrixType>());
90 
91  using impl_type = BlockTriDiContainerDetails::ImplType<MatrixType>;
92  using block_crs_matrix_type = typename impl_type::tpetra_block_crs_matrix_type;
93 
94  impl_->A = Teuchos::rcp_dynamic_cast<const block_crs_matrix_type>(matrix);
96  (impl_->A.is_null(), "BlockTriDiContainer currently supports Tpetra::BlockCrsMatrix only.");
97 
98  impl_->tpetra_importer = Teuchos::null;
99  impl_->async_importer = Teuchos::null;
100 
101  if (importer.is_null()) // there is no given importer, then create one
102  if (useSeqMethod)
103  impl_->tpetra_importer = BlockTriDiContainerDetails::createBlockCrsTpetraImporter<MatrixType>(impl_->A);
104  else
105  impl_->async_importer = BlockTriDiContainerDetails::createBlockCrsAsyncImporter<MatrixType>(impl_->A);
106  else
107  impl_->tpetra_importer = importer; // if there is a given importer, use it
108 
109  // as a result, there are
110  // 1) tpetra_importer is null , async_importer is null (no need for importer)
111  // 2) tpetra_importer is NOT null , async_importer is null (sequential method is used)
112  // 3) tpetra_importer is null , async_importer is NOT null (async method is used)
113 
114  // temporary disabling
115  impl_->overlap_communication_and_computation = overlapCommAndComp;
116 
117  impl_->Z = typename impl_type::tpetra_multivector_type();
118  impl_->W = typename impl_type::impl_scalar_type_1d_view();
119 
120  impl_->part_interface = BlockTriDiContainerDetails::createPartInterface<MatrixType>(impl_->A, partitions);
121  impl_->block_tridiags = BlockTriDiContainerDetails::createBlockTridiags<MatrixType>(impl_->part_interface);
122  impl_->norm_manager = BlockTriDiContainerDetails::NormManager<MatrixType>(impl_->A->getComm());
123  }
124 
125  template <typename MatrixType>
126  void
127  BlockTriDiContainer<MatrixType, BlockTriDiContainerDetails::ImplSimdTag>
128  ::clearInternal ()
129  {
130  using impl_type = BlockTriDiContainerDetails::ImplType<MatrixType>;
131  using part_interface_type = BlockTriDiContainerDetails::PartInterface<MatrixType>;
132  using block_tridiags_type = BlockTriDiContainerDetails::BlockTridiags<MatrixType>;
133  using amd_type = BlockTriDiContainerDetails::AmD<MatrixType>;
134  using norm_manager_type = BlockTriDiContainerDetails::NormManager<MatrixType>;
135 
136  impl_->A = Teuchos::null;
137  impl_->tpetra_importer = Teuchos::null;
138  impl_->async_importer = Teuchos::null;
139 
140  impl_->Z = typename impl_type::tpetra_multivector_type();
141  impl_->W = typename impl_type::impl_scalar_type_1d_view();
142 
143  impl_->part_interface = part_interface_type();
144  impl_->block_tridiags = block_tridiags_type();
145  impl_->a_minus_d = amd_type();
146  impl_->work = typename impl_type::vector_type_1d_view();
147  impl_->norm_manager = norm_manager_type();
148 
149  impl_ = Teuchos::null;
150  }
151 
152  template <typename MatrixType>
156  const Teuchos::RCP<const import_type>& importer,
157  bool pointIndexed)
158  : Container<MatrixType>(matrix, partitions, pointIndexed)
159  {
160  const bool useSeqMethod = false;
161  const bool overlapCommAndComp = false;
162  initInternal(matrix, partitions, importer, overlapCommAndComp, useSeqMethod);
163  }
164 
165  template <typename MatrixType>
169  const bool overlapCommAndComp, const bool useSeqMethod)
170  : Container<MatrixType>(matrix, partitions, false)
171  {
172  initInternal(matrix, partitions, Teuchos::null, overlapCommAndComp, useSeqMethod);
173  }
174 
175  template <typename MatrixType>
178  {
179  }
180 
181  template <typename MatrixType>
182  void
185  {
186  // the solver doesn't currently take any parameters
187  }
188 
189  template <typename MatrixType>
190  void
193  {
194  this->IsInitialized_ = true;
195  // We assume that if you called this method, you intend to recompute
196  // everything.
197  this->IsComputed_ = false;
198  TEUCHOS_ASSERT(!impl_->A.is_null()); // when initInternal is called, A_ must be set
199  {
200  BlockTriDiContainerDetails::performSymbolicPhase<MatrixType>
201  (impl_->A,
202  impl_->part_interface, impl_->block_tridiags,
203  impl_->a_minus_d,
204  impl_->overlap_communication_and_computation);
205  }
206  }
207 
208  template <typename MatrixType>
209  void
212  {
213  this->IsComputed_ = false;
214  if (!this->isInitialized())
215  this->initialize();
216  {
217  BlockTriDiContainerDetails::performNumericPhase<MatrixType>
218  (impl_->A,
219  impl_->part_interface, impl_->block_tridiags,
220  Kokkos::ArithTraits<magnitude_type>::zero());
221  }
222  this->IsComputed_ = true;
223  }
224 
225  template <typename MatrixType>
226  void
229  {
230  clearInternal();
231  this->IsInitialized_ = false;
232  this->IsComputed_ = false;
234  }
235 
236  template <typename MatrixType>
237  void
238  BlockTriDiContainer<MatrixType, BlockTriDiContainerDetails::ImplSimdTag>
239  ::applyInverseJacobi (const mv_type& X, mv_type& Y, scalar_type dampingFactor,
240  bool zeroStartingSolution, int numSweeps) const
241  {
242  const magnitude_type tol = Kokkos::ArithTraits<magnitude_type>::zero();
243  const int check_tol_every = 1;
244 
245  BlockTriDiContainerDetails::applyInverseJacobi<MatrixType>
246  (impl_->A,
247  impl_->tpetra_importer,
248  impl_->async_importer,
249  impl_->overlap_communication_and_computation,
250  X, Y, impl_->Z, impl_->W,
251  impl_->part_interface, impl_->block_tridiags, impl_->a_minus_d,
252  impl_->work,
253  impl_->norm_manager,
254  dampingFactor,
255  zeroStartingSolution,
256  numSweeps,
257  tol,
258  check_tol_every);
259  }
260 
261  template <typename MatrixType>
262  typename BlockTriDiContainer<MatrixType, BlockTriDiContainerDetails::ImplSimdTag>::ComputeParameters
265  {
266  return ComputeParameters();
267  }
268 
269  template <typename MatrixType>
270  void
272  ::compute (const ComputeParameters& in)
273  {
274  this->IsComputed_ = false;
275  if (!this->isInitialized())
276  this->initialize();
277  {
278  BlockTriDiContainerDetails::performNumericPhase<MatrixType>
279  (impl_->A,
280  impl_->part_interface, impl_->block_tridiags,
281  in.addRadiallyToDiagonal);
282  }
283  this->IsComputed_ = true;
284  }
285 
286  template <typename MatrixType>
290  {
291  ApplyParameters in;
292  in.dampingFactor = Teuchos::ScalarTraits<scalar_type>::one();
293  return in;
294  }
295 
296  template <typename MatrixType>
297  int
299  ::applyInverseJacobi (const mv_type& X, mv_type& Y,
300  const ApplyParameters& in) const
301  {
302  int r_val = 0;
303  {
304  r_val = BlockTriDiContainerDetails::applyInverseJacobi<MatrixType>
305  (impl_->A,
306  impl_->tpetra_importer,
307  impl_->async_importer,
308  impl_->overlap_communication_and_computation,
309  X, Y, impl_->Z, impl_->W,
310  impl_->part_interface, impl_->block_tridiags, impl_->a_minus_d,
311  impl_->work,
312  impl_->norm_manager,
313  in.dampingFactor,
314  in.zeroStartingSolution,
315  in.maxNumSweeps,
316  in.tolerance,
317  in.checkToleranceEvery);
318  }
319  return r_val;
320  }
321 
322  template <typename MatrixType>
325  ::getNorms0 () const {
326  return impl_->norm_manager.getNorms0();
327  }
328 
329  template <typename MatrixType>
333  return impl_->norm_manager.getNormsFinal();
334  }
335 
336  template <typename MatrixType>
337  void
339  ::apply (HostView /* X */, HostView /* Y */, int /* blockIndex */, Teuchos::ETransp /* mode */,
340  scalar_type /* alpha */, scalar_type /* beta */) const
341  {
342  TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "BlockTriDiContainer::apply is not implemented. You may have reached this message "
343  << "because you want to use this container's performance-portable Jacobi iteration. In "
344  << "that case, set \"relaxation: type\" to \"MT Split Jacobi\" rather than \"Jacobi\".");
345  }
346 
347  template <typename MatrixType>
348  void
350  ::weightedApply (HostView /* X */, HostView /* Y */, HostView /* D */, int /* blockIndex */,
351  Teuchos::ETransp /* mode */, scalar_type /* alpha */, scalar_type /* beta */) const
352  {
353  TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "BlockTriDiContainer::weightedApply is not implemented.");
354  }
355 
356  template <typename MatrixType>
357  std::ostream&
359  ::print (std::ostream& os) const
360  {
361  Teuchos::FancyOStream fos(Teuchos::rcp(&os,false));
362  fos.setOutputToRootOnly(0);
363  describe(fos);
364  return os;
365  }
366 
367  template <typename MatrixType>
368  std::string
371  {
372  std::ostringstream oss;
374  if (this->isInitialized()) {
375  if (this->isComputed()) {
376  oss << "{status = initialized, computed";
377  }
378  else {
379  oss << "{status = initialized, not computed";
380  }
381  }
382  else {
383  oss << "{status = not initialized, not computed";
384  }
385 
386  oss << "}";
387  return oss.str();
388  }
389 
390  template <typename MatrixType>
391  void
394  const Teuchos::EVerbosityLevel verbLevel) const
395  {
396  using std::endl;
397  if(verbLevel==Teuchos::VERB_NONE) return;
398  os << "================================================================================" << endl
399  << "Ifpack2::BlockTriDiContainer" << endl
400  << "Number of blocks = " << this->numBlocks_ << endl
401  << "isInitialized() = " << this->IsInitialized_ << endl
402  << "isComputed() = " << this->IsComputed_ << endl
403  << "================================================================================" << endl
404  << endl;
405  }
406 
407  template <typename MatrixType>
408  std::string
410  ::getName() { return "Ifpack2::BlockTriDiContainer::ImplSimdTag"; }
411 
412 #define IFPACK2_BLOCKTRIDICONTAINER_INSTANT(S,LO,GO,N) \
413  template class Ifpack2::BlockTriDiContainer< Tpetra::RowMatrix<S, LO, GO, N> >;
414 }
415 #endif
Ifpack2::BlockTriDiContainer class declaration.
Store and solve local block tridiagonal linear problems.
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:134
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
#define TEUCHOS_TEST_FOR_EXCEPT_MSG(throw_exception_test, msg)
BlockTriDiContainer(const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< local_ordinal_type > > &partitions, const Teuchos::RCP< const import_type > &importer, bool pointIndexed)
Constructor.
Definition: Ifpack2_BlockTriDiContainer_def.hpp:154
virtual std::string description() const
basic_FancyOStream & setOutputToRootOnly(const int rootRank)
Interface for creating and solving a set of local linear problems.
Definition: Ifpack2_Container_decl.hpp:112
ComputeParameters createDefaultComputeParameters() const
Create a ComputeParameters struct with default values.
Definition: Ifpack2_BlockTriDiContainer_def.hpp:264
#define TEUCHOS_ASSERT(assertion_test)
bool is_null() const