OGS
PETScNonlinearSolver.cpp
Go to the documentation of this file.
1 
11 #ifdef USE_PETSC
12 
13 #include "PETScNonlinearSolver.h"
14 
15 #include <petscmat.h>
16 #include <petscvec.h>
17 
18 #include "BaseLib/RunTime.h"
19 
20 namespace
21 {
23 {
26  std::vector<GlobalVector*>& x;
27  std::vector<GlobalVector*> const& x_prev;
30  int const process_id;
31 };
32 
33 PetscErrorCode updateResidual(SNES /*snes*/, Vec x, Vec petsc_r,
34  void* petsc_context)
35 {
36  auto context = static_cast<PetscContext*>(petsc_context);
37 
38  DBUG("PETScNonlinearSolver: residual callback called.");
39 
40  VecCopy(x, context->x[context->process_id]->getRawVector());
41 
42  // Assemble in ogs context.
43  BaseLib::RunTime time_assembly;
44  time_assembly.start();
45  context->system->assemble(context->x, context->x_prev, context->process_id);
46 
47  INFO("[time] Assembly took {} s.", time_assembly.elapsed());
48  context->system->getResidual(*context->x[context->process_id],
49  *context->x_prev[context->process_id],
50  *context->r);
51  context->r->finalizeAssembly();
52  context->J->finalizeAssembly();
53 
54  context->system->getJacobian(*context->J);
55  context->system->applyKnownSolutionsNewton(
56  *context->J, *context->r, *context->x[context->process_id]);
57 
58  VecCopy(context->r->getRawVector(), petsc_r);
59 
60  return 0;
61 }
62 
63 PetscErrorCode updateJacobian(SNES /*snes*/, Vec /*x*/, Mat J,
64  Mat /*same as J*/, void* petsc_context)
65 {
66  DBUG("PETScNonlinearSolver: jacobian callback called.");
67  // Assume the system is already assembled.
68  // Copy the results into petsc vectors.
69 
70  auto context = static_cast<PetscContext*>(petsc_context);
71  MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);
72  MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);
73  MatCopy(context->J->getRawMatrix(), J, DIFFERENT_NONZERO_PATTERN);
74  return 0;
75 }
76 } // namespace
77 
78 namespace NumLib
79 {
81  GlobalLinearSolver& /*linear_solver*/, int const maxiter,
82  std::string prefix)
83 {
84  SNESCreate(PETSC_COMM_WORLD, &_snes_solver);
85  if (!prefix.empty())
86  {
87  prefix = prefix + "_";
88  SNESSetOptionsPrefix(_snes_solver, prefix.c_str());
89  }
90  // force SNESSolve() to take at least one iteration regardless of the
91  // initial residual norm
92  SNESSetForceIteration(_snes_solver, PETSC_TRUE);
93 
94  // Set the maximum iterations.
95  PetscReal atol, rtol, stol;
96  PetscInt maxf;
97  SNESGetTolerances(_snes_solver, &atol, &rtol, &stol, nullptr, &maxf);
98  SNESSetTolerances(_snes_solver, atol, rtol, stol, maxiter, maxf);
99 
100  SNESSetFromOptions(_snes_solver);
101 #ifndef NDEBUG
102  PetscOptionsView(nullptr, PETSC_VIEWER_STDOUT_WORLD);
103 #endif // NDEBUG
104 }
105 
107  ConvergenceCriterion& conv_crit)
108 {
109  _equation_system = &eq;
110  _convergence_criterion = &conv_crit;
111 }
112 
114  bool const value)
115 {
117 }
118 
120  std::vector<GlobalVector*> const& /*x*/,
121  std::vector<GlobalVector*> const& /*x_prev*/, int const /*process_id*/)
122 {
124  {
125  return;
126  }
127 
128  OGS_FATAL(
129  "Non-equilibrium initial residuum is not implemented for the "
130  "PETScNonlinearSolver.");
131 }
132 
134  std::vector<GlobalVector*>& x,
135  std::vector<GlobalVector*> const& x_prev,
136  std::function<void(
137  int,
138  std::vector<GlobalVector*> const&)> const& /*postIterationCallback*/,
139  int const process_id)
140 {
141  DBUG("PETScNonlinearSolver: solve()");
142  using TimeDiscretizedSystem =
145 
146  auto* system = static_cast<TimeDiscretizedSystem*>(_equation_system);
147 
148  DBUG("PETScNonlinearSolver: create vectors");
149  // r and J on which the ogs assembly operates.
151  system->getMatrixSpecifications(process_id), _jacobian_id);
153  system->getMatrixSpecifications(process_id), _residual_id);
154 
155  // Vectors and matrices used by SNES for solutions. These will be locked
156  // during the SNES' solve call.
158  system->getMatrixSpecifications(process_id), _petsc_jacobian_id);
159  BaseLib::RunTime timer_dirichlet;
160  double time_dirichlet = 0.0;
161 
162  timer_dirichlet.start();
163  system->computeKnownSolutions(*x[process_id], process_id);
164  system->applyKnownSolutions(*x[process_id]);
165  time_dirichlet += timer_dirichlet.elapsed();
166  INFO("[time] Applying Dirichlet BCs took {} s.", time_dirichlet);
167 
169  system->getMatrixSpecifications(process_id), _petsc_residual_id);
171  system->getMatrixSpecifications(process_id), _petsc_x_id);
172  MathLib::LinAlg::copy(*x[process_id], x_snes); // Initial guess.
173 
174  PetscContext petsc_context{_equation_system, x, x_prev, &r, &J, process_id};
175 
176  DBUG("PETScNonlinearSolver: set function");
177  SNESSetFunction(_snes_solver, r_snes.getRawVector(), updateResidual,
178  &petsc_context);
179 
180  DBUG("PETScNonlinearSolver: set jacobian");
181  // The jacobian and the preconditioner matrices point to the same location.
182  SNESSetJacobian(_snes_solver, J_snes.getRawMatrix(), J_snes.getRawMatrix(),
183  updateJacobian, &petsc_context);
184 
185  std::unique_ptr<GlobalVector> xl = nullptr;
186  std::unique_ptr<GlobalVector> xu = nullptr;
187 
188  SNESType snes_type;
189  SNESGetType(_snes_solver, &snes_type);
190  if ((std::strcmp(snes_type, SNESVINEWTONRSLS) == 0) ||
191  (std::strcmp(snes_type, SNESVINEWTONSSLS) == 0))
192  {
193  // Set optional constraints via callback.
194  DBUG("PETScNonlinearSolver: set constraints");
196  system->getMatrixSpecifications(process_id));
198  system->getMatrixSpecifications(process_id));
199 
200  system->updateConstraints(*xl, *xu, process_id);
203 
204  SNESVISetVariableBounds(_snes_solver, xl->getRawVector(),
205  xu->getRawVector());
206  }
207 
208  DBUG("PETScNonlinearSolver: call SNESSolve");
209  SNESSolve(_snes_solver, nullptr, x_snes.getRawVector());
210 
211  SNESConvergedReason reason = SNES_CONVERGED_ITERATING;
212  SNESGetConvergedReason(_snes_solver, &reason);
213  INFO("PETSsSNES convergence reason {}.", reason);
214 
215  PetscInt iterations;
216  SNESGetIterationNumber(_snes_solver, &iterations);
217  INFO("PETScSNES used {} iterations.", iterations);
218 
219  // Copy back the solution.
220  MathLib::LinAlg::copy(x_snes, *x[process_id]);
221 
227 
228  return {reason >= 0, iterations};
229 }
230 
231 } // namespace NumLib
232 #endif // USE_PETSC
#define OGS_FATAL(...)
Definition: Error.h:26
void INFO(char const *fmt, Args const &... args)
Definition: Logging.h:32
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
Definition of the RunTime class.
Count the running time.
Definition: RunTime.h:29
double elapsed() const
Get the elapsed time in seconds.
Definition: RunTime.h:42
void start()
Start the timer.
Definition: RunTime.h:32
Global vector based on Eigen vector.
Definition: EigenVector.h:26
virtual GlobalMatrix & getMatrix(std::size_t &id)=0
Get an uninitialized matrix with the given id.
virtual void releaseMatrix(GlobalMatrix const &A)=0
ConvergenceCriterion * _convergence_criterion
NonlinearSolverStatus solve(std::vector< GlobalVector * > &x, std::vector< GlobalVector * > const &x_prev, std::function< void(int, std::vector< GlobalVector * > const &)> const &postIterationCallback, int const process_id) override
void setEquationSystem(System &eq, ConvergenceCriterion &conv_crit)
Set the nonlinear equation system that will be solved.
PETScNonlinearSolver(GlobalLinearSolver &linear_solver, int maxiter, std::string prefix)
std::size_t _jacobian_id
ID of the Jacobian matrix.
void compensateNonEquilibriumInitialResiduum(bool const value)
std::size_t _residual_id
ID of the residual vector.
void calculateNonEquilibriumInitialResiduum(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id) override
virtual GlobalVector & getVector(std::size_t &id)=0
Get an uninitialized vector with the given id.
virtual void releaseVector(GlobalVector const &x)=0
void copy(PETScVector const &x, PETScVector &y)
Definition: LinAlg.cpp:37
void finalizeVectorAssembly(VEC_T &)
General function to finalize the vector assembly.
static const double r
PetscErrorCode updateResidual(SNES, Vec x, Vec petsc_r, void *petsc_context)
PetscErrorCode updateJacobian(SNES, Vec, Mat J, Mat, void *petsc_context)
static NUMLIB_EXPORT MatrixProvider & provider
static NUMLIB_EXPORT VectorProvider & provider
Status of the non-linear solver.