KrylovPreconditioner.h 5.09 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
/******************************************************************************
 *
 * 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 {
33
  
34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53

  /**
   * \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 >
54
  struct KrylovPreconditioner : public ITL_BasePreconditioner< MatrixType, VectorType >
55 56 57 58 59 60 61 62 63 64 65 66 67
  {
    typedef ITL_BasePreconditioner< MatrixType, VectorType > super;
    
    class Creator : public CreatorInterfaceName< super >
    {
    public:
      virtual ~Creator() {}
      
      super* create() { 
	return new KrylovPreconditioner(this->name);
      }
    };
    
68
    KrylovPreconditioner(std::string name) : fullMatrix(nullptr), solver(nullptr), runner(nullptr)
69 70
    {

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

      // === read solver-name ===
86
      std::string solverType("0");
87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106
      Parameters::get(initFileStr, solverType);
      
      if (backend != "0" && backend != "no" && backend != "")
	solverType = backend + "_" + solverType;
    
      LinearSolverCreator *solverCreator = 
	dynamic_cast<LinearSolverCreator*>(CreatorMap<LinearSolver>::getCreator(solverType, initFileStr));
      TEST_EXIT(solverCreator)
	("No valid solver type found in parameter \"%s\"\n", initFileStr.c_str());
      solverCreator->setName(initFileStr);
      solver = solverCreator->create();
      
      runner = dynamic_cast<MTL4Runner< MatrixType, VectorType >*>(solver->getRunner());
    }
    
    virtual ~KrylovPreconditioner()
    {
      delete solver;
    }
    
107
    /// Implementation of \ref ITL_BasePreconditioner::init()
108
    void init(const SolverMatrix<Matrix<DOFMatrix*> >& A, const MatrixType& fullMatrix_) override
109 110 111 112 113
    {
      fullMatrix = &fullMatrix_;
      runner->init(A, fullMatrix_);
    }
    
114 115
    /// Implementation of \ref OEMPreconditioner::init()
    void exit() override
116 117 118 119
    {
      runner->exit();
    }
    
120 121
    /// Implementation of \ref ITL_BasePreconditioner::solve()
    void solve(const VectorType& b, VectorType& x) const override
122 123 124 125 126
    {
      initVector(x);
      runner->solve(*fullMatrix, x, b);
    }
    
127 128
    /// Implementation of \ref ITL_BasePreconditioner::adjoint_solve()
    void adjoint_solve(const VectorType& b, VectorType& x) const override
129 130 131 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 161 162 163 164 165 166 167 168 169 170 171
    {
      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;
    
    LinearSolver* solver;
    MTL4Runner< MatrixType, VectorType >* runner;
  };
  
  
#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