PetscSolver.cc 4.6 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
/******************************************************************************
 *
 * 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.
 * 
 ******************************************************************************/


22
#ifdef HAVE_SEQ_PETSC
23
24
25

#include "solver/PetscSolver.h"

26
27
using namespace std;

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
namespace AMDiS {  
  
  PetscRunner::PetscRunner(LinearSolver* oem_)
  : oem(*oem_),
    kspPrefix("")
  {
    PetscParameters params;
    bool matSolverPackage = false;
    
    // set the parameter prefix
    Parameters::get(oem.getName() + "->petsc prefix", kspPrefix); // important when more than one solver to configure
    
    // set the solver
    std::string solverName = "petsc";
    Parameters::get(oem.getName(), solverName);
    if (solverName == "petsc") 
      Parameters::get(oem.getName() + "->ksp_type", solverName);
    
    std::string kspSolver = params.solverMap[solverName];
    
    if (params.matSolverPackage[kspSolver]) {
      // direct solvers
      PetscOptionsInsertString(("-" + kspPrefix + "ksp_type preonly").c_str());
      PetscOptionsInsertString(("-" + kspPrefix + "pc_type lu").c_str());
      PetscOptionsInsertString(("-" + kspPrefix + "pc_factor_mat_solver_package " + (kspSolver != "direct" ? kspSolver : "umfpack")).c_str());
      oem.setMaxIterations(1);
      zeroStartVector = true;
      matSolverPackage = true;
    } else if (!params.emptyParam[kspSolver]) {    
      // other solvers
      PetscOptionsInsertString(("-" + kspPrefix + "ksp_type " + kspSolver).c_str());
    }
    
    // set the preconditioner
    string precon = "";
    Parameters::get(oem.getName() + "->pc_type", precon);
    if (!precon.size())
      Parameters::get(oem.getName() + "->left precon", precon);
66
67
    if (!precon.size())
      Parameters::get(oem.getName() + "->right precon", precon);
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
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
    if (!matSolverPackage && !params.emptyParam[precon]) {
      precon = (params.preconMap.find(precon) != params.preconMap.end() ? params.preconMap[precon] : precon);
      PetscOptionsInsertString(("-" + kspPrefix + "pc_type " + precon).c_str());
    }
    
    if (oem.getInfo() >= 20)
      PetscOptionsInsertString(("-" + kspPrefix + "ksp_monitor_true_residual").c_str());
    else if (oem.getInfo() >= 10)
      PetscOptionsInsertString(("-" + kspPrefix + "ksp_monitor").c_str());
    
    // command line string
    string kspString = "";
    Parameters::get(oem.getName() + "->ksp", kspString);
    if (kspString != "")
      PetscOptionsInsertString(kspString.c_str());
  }

  
  void PetscRunner::init(const BlockMatrix& A, const Mat& fullMatrix)
  {
    createSubSolver(ksp, fullMatrix, kspPrefix);
    
    // set tolerance options from LinearSolver-parameters
    KSPSetTolerances(ksp, oem.getRelative(), oem.getTolerance(), PETSC_DEFAULT, oem.getMaxIterations());

    // Do not delete the solution vector, use it for the initial guess.
    if (!zeroStartVector)
      KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
  
    KSPGetPC(ksp, &pc);  
    PCSetFromOptions(pc);
  }

  
  int PetscRunner::solve(const Mat& A, Vec& x, const Vec& b) 
  {
    if (zeroStartVector)
      VecSet(x, 0.0);
    
    PetscErrorCode solverError = KSPSolve(ksp, b, x);
    
    PetscInt nIter = 0;
    KSPGetIterationNumber(ksp, &nIter);
    PetscReal residual_norm = 0.0;
    KSPGetResidualNorm(ksp, &residual_norm);
    
    oem.setErrorCode(solverError);
    oem.setIterations(nIter);
    oem.setResidual(residual_norm);
    
    return solverError;
  }

  
  void PetscRunner::createSubSolver(KSP &ksp_, Mat m, string kspPrefix_)
  {
    KSPCreate(PETSC_COMM_SELF, &ksp_);
    KSPSetOperators(ksp_, m, m, SAME_NONZERO_PATTERN); 
    KSPSetOptionsPrefix(ksp_, kspPrefix_.c_str());
    KSPSetFromOptions(ksp_);
  }


  void PetscRunner::setSolver(KSP ksp_, std::string kspPrefix_,
			      KSPType kspType,  PCType pcType, 
			      PetscReal rtol, PetscReal atol, PetscInt maxIt)
  {
    KSPSetType(ksp_, kspType);
    KSPSetTolerances(ksp_, rtol, atol, PETSC_DEFAULT, maxIt);
    if (kspPrefix_ != "")
      KSPSetOptionsPrefix(ksp_, kspPrefix_.c_str());
    KSPSetFromOptions(ksp_);

    PC pc_;
    KSPGetPC(ksp_, &pc_);
143
    PCSetType(pc_, pcType);
144
145
146
147
    PCSetFromOptions(pc_);
  }
}

148
#endif // HAVE_SEQ_PETSC