Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ifpack2_BlockTriDiContainer_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 
43 #ifndef IFPACK2_BLOCKTRIDICONTAINER_DECL_HPP
44 #define IFPACK2_BLOCKTRIDICONTAINER_DECL_HPP
45 
48 
49 #include "Ifpack2_config.h"
50 #include "Ifpack2_Container.hpp"
51 #include "Tpetra_MultiVector.hpp"
52 #include "Tpetra_Map.hpp"
53 #include "Tpetra_RowMatrix.hpp"
54 #include "Tpetra_BlockCrsMatrix_decl.hpp"
55 #include <type_traits>
56 #include <string>
57 
58 namespace Ifpack2 {
59 
107 
111  namespace BlockTriDiContainerDetails {
115  struct ImplNotAvailTag {};
116  struct ImplSimdTag {};
117  struct ImplSacadoTag {};
118 
119  template<typename T> struct ImplTag { typedef ImplNotAvailTag type; };
120  template<> struct ImplTag<float> { typedef ImplSimdTag type; };
121  template<> struct ImplTag<double> { typedef ImplSimdTag type; };
122  template<> struct ImplTag<std::complex<float> > { typedef ImplSimdTag type; };
123  template<> struct ImplTag<std::complex<double> > { typedef ImplSimdTag type; };
124 
126  template<typename MatrixType> struct ImplObject;
127  }
128 
132  template <typename MatrixType,
135 
141  template <typename MatrixType>
142  class BlockTriDiContainer<MatrixType,BlockTriDiContainerDetails::ImplSimdTag>
143  : public Container<MatrixType> {
145 
146  private:
152  typedef MatrixType matrix_type;
153 
155  typedef typename Container<MatrixType>::scalar_type scalar_type;
157  typedef typename Kokkos::ArithTraits<scalar_type>::magnitudeType magnitude_type;
159  typedef typename Container<MatrixType>::local_ordinal_type local_ordinal_type;
161  typedef typename Container<MatrixType>::global_ordinal_type global_ordinal_type;
163  typedef typename Container<MatrixType>::node_type node_type;
164 
165  typedef typename Container<MatrixType>::mv_type mv_type;
166  typedef typename Container<MatrixType>::map_type map_type;
167  typedef typename Container<MatrixType>::vector_type vector_type;
169  typedef typename Container<MatrixType>::import_type import_type;
170 
171  typedef typename Container<MatrixType>::HostView host_view_type;
172  typedef host_view_type HostView;
173  //typedef Tpetra::MultiVector<local_scalar_type, local_ordinal_type, global_ordinal_type, node_type> local_mv_type;
174  //typedef typename Kokkos::View<local_scalar_type**, Kokkos::HostSpace> HostViewLocal;
175 
176  typedef Tpetra::Experimental::BlockCrsMatrix
177  <scalar_type,local_ordinal_type,global_ordinal_type,node_type> block_crs_matrix_type;
178 
185  typedef typename Container<MatrixType>::row_matrix_type row_matrix_type;
186 
187  static_assert (std::is_same<MatrixType, row_matrix_type>::value,
188  "Ifpack2::BlockTriDiContainer: MatrixType must be a Tpetra::RowMatrix specialization.");
190  public:
192 
193 
211  const Teuchos::RCP<const import_type>& importer,
212  int OverlapLevel,
213  scalar_type DampingFactor);
214 
231  bool overlapCommAndComp = false, bool useSequentialMethod = false);
232 
234  ~BlockTriDiContainer () override;
235 
236  struct ComputeParameters {
243  magnitude_type addRadiallyToDiagonal = Kokkos::ArithTraits<magnitude_type>::zero();
244  };
245 
247  struct ApplyParameters {
250  bool zeroStartingSolution = false;
252  scalar_type dampingFactor = Kokkos::ArithTraits<scalar_type>::one();
255  int maxNumSweeps = 1;
265  magnitude_type tolerance = Kokkos::ArithTraits<magnitude_type>::zero();
273  int checkToleranceEvery = 1;
274  };
275 
277 
279 
281  bool isInitialized() const override;
282 
284  bool isComputed() const override;
285 
287  void setParameters(const Teuchos::ParameterList& List) override;
288 
289  void clearBlocks() override;
290 
292 
294 
296  void initialize () override;
297 
299  void compute () override;
300 
301  // \brief Compute <tt>Y := D^{-1} (X - R*Y)</tt>.
302  void applyInverseJacobi (const mv_type& X, mv_type& Y,
303  bool zeroStartingSolution = false,
304  int numSweeps = 1) const override;
305 
307  ComputeParameters createDefaultComputeParameters () const;
308 
320  void compute (const ComputeParameters& input);
321 
323  ApplyParameters createDefaultApplyParameters () const;
324 
331  int applyInverseJacobi (const mv_type& X, mv_type& Y,
332  const ApplyParameters& input) const;
333 
337  const Teuchos::ArrayRCP<const magnitude_type> getNorms0 () const;
338 
342  const Teuchos::ArrayRCP<const magnitude_type> getNormsFinal () const;
343 
346  void
347  apply (host_view_type& X,
348  host_view_type& Y,
349  int blockIndex,
350  int stride,
352  scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
353  scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override;
354 
357  void
358  weightedApply (host_view_type& X,
359  host_view_type& Y,
360  host_view_type& W,
361  int blockIndex,
362  int stride,
364  scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
365  scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override;
366 
368 
370 
374  std::ostream& print (std::ostream& os) const override;
375 
377 
379 
381  std::string description () const override;
382 
384  void
385  describe (Teuchos::FancyOStream &out,
386  const Teuchos::EVerbosityLevel verbLevel =
388 
390 
392  static std::string getName();
393 
394  private:
397 
399  bool IsInitialized_;
400 
402  bool IsComputed_;
403 
404  // hide details of impl using ImplObj; finally I understand why AMB did that way.
406 
407  // initialize distributed and local objects
408  void initInternal (const Teuchos::RCP<const row_matrix_type>& matrix,
410  const Teuchos::RCP<const import_type> &importer,
411  const bool overlapCommAndComp,
412  const bool useSeqMethod);
413 
414  void clearInternal();
415  };
416 
424  template <typename MatrixType>
425  class BlockTriDiContainer<MatrixType,BlockTriDiContainerDetails::ImplNotAvailTag>
426  : public Container<MatrixType> {
427  private:
428  typedef typename Container<MatrixType>::scalar_type scalar_type;
429  typedef typename Kokkos::ArithTraits<scalar_type>::magnitudeType magnitude_type;
430  typedef typename Container<MatrixType>::local_ordinal_type local_ordinal_type;
431  typedef typename Container<MatrixType>::global_ordinal_type global_ordinal_type;
432 
433  typedef typename Container<MatrixType>::mv_type mv_type;
435  typedef typename Container<MatrixType>::import_type import_type;
436 
437  typedef typename Container<MatrixType>::HostView host_view_type;
438  typedef typename Container<MatrixType>::row_matrix_type row_matrix_type;
439 
440  static_assert (std::is_same<MatrixType, row_matrix_type>::value,
441  "Ifpack2::BlockTriDiContainer: MatrixType must be a Tpetra::RowMatrix specialization.");
442  public:
443 
446  const Teuchos::RCP<const import_type>& importer,
447  int OverlapLevel,
448  scalar_type DampingFactor)
449  : Container<MatrixType>(matrix, partitions, importer, OverlapLevel, DampingFactor) {
450  TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "Error: BlockTriDiContainer is not available for this scalar_type");
451  }
452 
453  bool isInitialized() const override { return 0; }
454  bool isComputed() const override { return 0; }
455  void setParameters(const Teuchos::ParameterList& List) override {}
456  void clearBlocks() override {}
457 
458  void initialize () override {}
459  void compute () override {}
460  void applyInverseJacobi (const mv_type& X, mv_type& Y,
461  bool zeroStartingSolution = false,
462  int numSweeps = 1) const override {}
463 
464  void
465  apply (host_view_type& X,
466  host_view_type& Y,
467  int blockIndex,
468  int stride,
470  scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
471  scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override {}
472 
473  void
474  weightedApply (host_view_type& X,
475  host_view_type& Y,
476  host_view_type& W,
477  int blockIndex,
478  int stride,
480  scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
481  scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const override {}
482 
483  std::ostream& print (std::ostream& os) const override {
484  return os << "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
485  }
486 
487  std::string description () const override {
488  return "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
489  }
490 
491  void
492  describe (Teuchos::FancyOStream &out,
493  const Teuchos::EVerbosityLevel verbLevel =
495  out << "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
496  }
497 
498  static std::string getName() {
499  return "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
500  }
501  };
502 
503 
504 } // namespace Ifpack2
505 
506 #endif // IFPACK2_BLOCKTRIDICONTAINER_DECL_HPP
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:115
Ifpack2::Container class declaration.
void compute() override
Extract the local diagonal block and prepare the solver.
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:459
Store and solve local block tridiagonal linear problems.
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:134
void weightedApply(host_view_type &X, host_view_type &Y, host_view_type &W, int blockIndex, int stride, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const override
Compute Y := alpha * diag(D) * M^{-1} (diag(D) * X) + beta*Y.
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:474
Ifpack2::Partitioner:
Definition: Ifpack2_Partitioner.hpp:179
void initialize() override
Do all set-up operations that only require matrix structure.
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:458
#define TEUCHOS_TEST_FOR_EXCEPT_MSG(throw_exception_test, msg)
forward declaration
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:126
std::ostream & print(std::ostream &os) const override
Print basic information about the container to os.
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:483
void setParameters(const Teuchos::ParameterList &List) override
Set parameters.
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:455
static const EVerbosityLevel verbLevel_default
void apply(host_view_type &X, host_view_type &Y, int blockIndex, int stride, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const override
Compute Y := alpha * M^{-1} X + beta*Y.
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:465
Interface for creating and solving a local linear problem.
Definition: Ifpack2_Container.hpp:114
void applyInverseJacobi(const mv_type &X, mv_type &Y, bool zeroStartingSolution=false, int numSweeps=1) const override
Compute Y := (1 - a) Y + a D^{-1} (X - R*Y).
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:460
bool isInitialized() const override
Return true if the container has been successfully initialized.
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:453
bool isComputed() const override
Return true if the container has been successfully computed.
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:454