HypreSolver.h 2.65 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
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
/******************************************************************************
 *
 * 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 HypreSolver.h */


#ifndef AMDIS_HYPRE_SOLVER_H
#define AMDIS_HYPRE_SOLVER_H

#ifdef MTL_HAS_HYPRE

#include "solver/MTL4Solver.h"
#include "solver/ITL_Preconditioner.h"
#include "MTL4Types.h"
#include "solver/itl/hypre.hpp"

#include <boost/numeric/itl/itl.hpp>
#include <boost/numeric/mtl/mtl.hpp>

namespace AMDiS {

  struct Hypre_Runner : MTL4Runner< MTLTypes::MTLMatrix, MTLTypes::MTLVector >
  {       
    typedef MTLTypes::MTLMatrix MatrixType;
    typedef MTLTypes::MTLVector VectorType;
    typedef MTL4Runner< MatrixType, VectorType > super;
    
    Hypre_Runner(LinearSolver* oemPtr)
    : oem(*oemPtr),
      solver(NULL)
    { 
      int cycleMode = 1, interpolation = 0;
      Parameters::get(oem.getName() + "->cycle mode", cycleMode);
      Parameters::get(oem.getName() + "->interpolation type", interpolation);
      
      config.maxIter(oem.getMaxIterations());
      config.interpolation(interpolation);
      config.cycle(cycleMode);
      config.tolerance(oem.getTolerance());
      config.printLevel(oem.getInfo());
    }
    
    ~Hypre_Runner()
    {
      if (solver) {
	delete solver;
	solver = NULL;
      }
    }      
    
    virtual void init(const typename super::BlockMatrix& A, const MatrixType& fullMatrix)
    {
      solver = new mtl::HypreSolverWrapper(fullMatrix);
    }
        
    virtual int solve(const MatrixType& A , VectorType& x, const VectorType& b)
    {      
      int error = (*solver)(x, b, config);      
      return error;
    }
    
    virtual void exit()
    {
      if (solver) {
	delete solver;
	solver = NULL;
      }
    }
    
  protected:
    LinearSolver& oem;
    mtl::HypreSolverWrapper *solver;
    mtl::AMGConfigurator config;
  };
  

  /**
   * \ingroup Solver
   * \class AMDiS::HypreSolver
   * \brief
   * Wrapper for the external HYPRE-AMG solver
   */
  typedef MTL4Solver< MTLTypes::MTLMatrix, MTLTypes::MTLVector, Hypre_Runner > HypreSolver;

} // namespace AMDiS

#endif // MTL_HAS_HYPRE

#endif // AMDIS_HYPRE_SOLVER_H