UmfPackSolver.h 2.99 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
/******************************************************************************
 *
 * 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 UmfPackSolver.h */

#ifndef AMDIS_UMFPACKSOLVER_H
#define AMDIS_UMFPACKSOLVER_H

#if defined HAVE_UMFPACK && defined MTL_HAS_UMFPACK

#include <iostream>
#include <boost/numeric/mtl/operation/two_norm.hpp>
#include <boost/numeric/mtl/interface/umfpack_solve.hpp>
#include "solver/MTL4Solver.h"

namespace AMDiS {
  
  template< typename MatrixType, typename VectorType >
  struct Umfpack_Runner : MTL4Runner<MatrixType, VectorType>
  {
    typedef MTL4Runner<MatrixType, VectorType> super;
    
    Umfpack_Runner(LinearSolver* oem_):
      oem(*oem_),
      solver(NULL),
      store_symbolic(0),
      symmetric_strategy(0),
      alloc_init(0.7)
    {
      Parameters::get(oem.getName() + "->store symbolic", store_symbolic); // ?
      Parameters::get(oem.getName() + "->symmetric strategy", symmetric_strategy);
      Parameters::get(oem.getName() + "->alloc init", alloc_init);
    }

    void init(const typename super::BlockMatrix& A, const MatrixType& fullMatrix)
    {
      if (solver != NULL) {
	delete solver;
	solver = NULL;
      }

      solver = new mtl::matrix::umfpack::solver<MatrixType>(fullMatrix, symmetric_strategy, alloc_init);
    }

    int solve(const MatrixType& A, VectorType& x, const VectorType& b) 
    { 
      FUNCNAME("Umfpack_Runner::solve()");
      TEST_EXIT(solver != NULL)("The umfpack solver was not initialized\n");

      int code = (*solver)(x, b);
      
Praetorius, Simon's avatar
Praetorius, Simon committed
70
71
72
73
74
75
      VectorType r(b); 
      r -= A * x; 
      double residual = two_norm(r);
      oem.setResidual(residual);
      oem.setErrorCode(code);
      
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
108
109
110
111
112
113
114
115
116
117
118
      return code;
    }

    void exit() {}

    ~Umfpack_Runner() 
    {
      if (solver != NULL) {	
	delete solver;
	solver = NULL;
      }
    }
    
  public:
    LinearSolver& oem;

  private:
    mtl::matrix::umfpack::solver<MatrixType> *solver;

    int store_symbolic;

    int symmetric_strategy;

    double alloc_init;
  };


  /**
   * \ingroup Solver
   * \class AMDiS::UmfPackSolver
   * \brief \implements MTL4Solver
   * Wrapper for the external UMFPACK solver:
   *   http://www.cise.ufl.edu/research/sparse/umfpack/
   *
   * This is a direct solver for large sparse matrices.
   */
  typedef MTL4Solver<MTLTypes::MTLMatrix, MTLTypes::MTLVector, Umfpack_Runner< MTLTypes::MTLMatrix, MTLTypes::MTLVector > > UmfPackSolver;
  
}

#endif // HAVE_UMFPACK

#endif // AMDIS_UMFPACKSOLVER_H