KrylovPreconditioner.h 5.17 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
/******************************************************************************
 *
 * AMDiS - Adaptive multidimensional simulations
 *
 * Copyright (C) 2013 Dresden University of Technology. All Rights Reserved.
 * Web: https://fusionforge.zih.tu-dresden.de/projects/amdis
 *
 * Authors: 
 * Simon Vey, Thomas Witkowski, Andreas Naumann, Simon Praetorius, et al.
 *
 * This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
 * WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
 *
 *
 * This file is part of AMDiS
 *
 * See also license.opensource.txt in the distribution.
 * 
 ******************************************************************************/


/** \file KrylovPreconditioner.h */

#ifndef AMDIS_KRYLOV_PRECONDITIONER_H
#define AMDIS_KRYLOV_PRECONDITIONER_H

#include "solver/ITL_Preconditioner.h"
#ifdef HAVE_PARALLEL_MTL4
#include "parallel/PITL_Solver.h"
#endif

namespace AMDiS {

  /**
   * \ingroup Solver
   * 
   * \brief Implements an inner solver for (flexible) krylov subspace solvers
   * 
   * Any AMDiS solver can be used as preconditioner for any other AMDiS solver.
   * Initfile-syntax:<br />
   * [...]->solver: any-krylov-solver<br />
   * [...]->solver->left precon: [krylov|solver]<br />
   * [...]->solver->left precon->solver: any-amdis-solver<br />
   * [...]->solver->left precon->solver->parameter1: ...<br />
   * [...]->solver->left precon->solver->parameter2: ...<br />
   * [...]->solver->left precon->solver->parameter3: ...
   * 
   * This Krylov-Preconditioner can also be used as preconditioner for the inner solver, so
   * that a sequence of iner-inner-solver can be loaded.
   **/

  template< typename MatrixType, typename VectorType >
53
  struct KrylovPreconditioner : ITL_PreconditionerBase< MatrixType, VectorType >
54
  {
55
56
    typedef ITL_PreconditionerBase< MatrixType, VectorType >  precon_base;
    typedef KrylovPreconditioner< MatrixType, VectorType >    self;
57
    
58
    class Creator : public CreatorInterfaceName<precon_base>
59
60
61
62
    {
    public:
      virtual ~Creator() {}
      
63
64
      precon_base* create() { 
	return new self(this->name);
65
66
67
      }
    };
    
68
69
70
71
    KrylovPreconditioner(std::string name) 
      : fullMatrix(NULL), 
	solver(NULL), 
	runner(NULL)
72
73
    {

74
#if defined HAVE_PARALLEL_PETSC
75
      std::string backend("p_petsc");
76
#elif defined HAVE_PARALLEL_MTL
77
      std::string backend("p_mtl");
78
#elif defined HAVE_PETSC || defined HAVE_SEQ_PETSC
79
      std::string backend("petsc");
80
#else
81
      std::string backend("mtl");
82
#endif
83
84
      
      // === read backend-name ===
85
      std::string initFileStr = name + "->solver";    
86
87
88
      Parameters::get(initFileStr + "->backend", backend);

      // === read solver-name ===
89
      std::string solverType("0");
90
91
92
93
94
95
      Parameters::get(initFileStr, solverType);
      
      if (backend != "0" && backend != "no" && backend != "")
	solverType = backend + "_" + solverType;
    
      LinearSolverCreator *solverCreator = 
96
	dynamic_cast<LinearSolverCreator*>(CreatorMap<LinearSolverInterface>::getCreator(solverType, initFileStr));
97
98
99
100
101
      TEST_EXIT(solverCreator)
	("No valid solver type found in parameter \"%s\"\n", initFileStr.c_str());
      solverCreator->setName(initFileStr);
      solver = solverCreator->create();
      
102
      runner = dynamic_cast<RunnerBase< MatrixType, VectorType >*>(solver->getRunner());
103
104
105
106
107
108
109
    }
    
    virtual ~KrylovPreconditioner()
    {
      delete solver;
    }
    
110
    /// Implementation of \ref ITL_PreconditionerBase::init()
111
    void init(const SolverMatrix<Matrix<DOFMatrix*> >& A, const MatrixType& fullMatrix_) override
112
113
114
115
116
    {
      fullMatrix = &fullMatrix_;
      runner->init(A, fullMatrix_);
    }
    
117
    /// Implementation of \ref PreconditionerInterface::init()
118
    void exit() override
119
120
121
122
    {
      runner->exit();
    }
    
123
    /// Implementation of \ref ITL_PreconditionerBase::solve()
124
    void solve(const VectorType& b, VectorType& x) const override
125
126
127
128
129
    {
      initVector(x);
      runner->solve(*fullMatrix, x, b);
    }
    
130
    /// Implementation of \ref ITL_PreconditionerBase::adjoint_solve()
131
    void adjoint_solve(const VectorType& b, VectorType& x) const override
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
    {
      initVector(x);
      runner->adjoint_solve(*fullMatrix, x, b);
    }
    
  protected: // methods
    
    template<typename VectorT>
    typename boost::enable_if< mtl::traits::is_distributed<VectorT>, void >::type
    initVector(VectorT& x) const
    {
#ifdef HAVE_PARALLEL_MTL4
      x.init_distribution(col_distribution(*fullMatrix), num_cols(*fullMatrix));
#endif
      set_to_zero(x);
    }
    
    template<typename VectorT>
    typename boost::disable_if< mtl::traits::is_distributed<VectorT>, void >::type
    initVector(VectorT& x) const
    {
      x.change_dim(num_cols(*fullMatrix));
      set_to_zero(x);
    }
    
  protected: // member variables

    const MatrixType* fullMatrix;
    
161
162
    LinearSolverInterface* solver;
    RunnerBase< MatrixType, VectorType >* runner;
163
164
165
166
167
168
169
170
171
172
173
174
  };
  
  
#ifdef HAVE_PARALLEL_MTL4
  typedef KrylovPreconditioner<MTLTypes::PMTLMatrix, MTLTypes::PMTLVector> KrylovPreconditionerParallel;
#endif
  typedef KrylovPreconditioner<MTLTypes::MTLMatrix, MTLTypes::MTLVector> KrylovPreconditionerSeq;

  
} // end namespace AMDiS

#endif // AMDIS_KRYLOV_PRECONDITIONER_H