petsc_nonlinear_solver.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #include "libmesh/libmesh_common.h"
21 
22 #ifdef LIBMESH_HAVE_PETSC
23 
24 // C++ includes
25 
26 // Local Includes
30 #include "libmesh/petsc_vector.h"
31 #include "libmesh/petsc_matrix.h"
32 #include "libmesh/dof_map.h"
33 #include "libmesh/preconditioner.h"
35 
36 /* DMlibMesh include. */
37 #include "libmesh/petscdmlibmesh.h"
38 
39 // Pick the right version of the petsc line search getter function name
40 #if PETSC_VERSION_LESS_THAN(3,4,0)
41 # define SNESGETLINESEARCH SNESGetSNESLineSearch
42 #else
43 # define SNESGETLINESEARCH SNESGetLineSearch
44 #endif
45 
46 namespace libMesh
47 {
49 {
50 public:
52  PetscErrorCode ierr_in) :
53  solver(solver_in),
54  sys(sys_in),
55  ierr(ierr_in)
56  {}
57 
60  PetscErrorCode ierr;
61 };
62 
64 libmesh_petsc_snes_residual_helper (SNES snes, Vec x, void * ctx)
65 {
66  LOG_SCOPE("residual()", "PetscNonlinearSolver");
67 
68  PetscErrorCode ierr = 0;
69 
70  libmesh_assert(x);
71  libmesh_assert(ctx);
72 
73  // No way to safety-check this cast, since we got a void *...
75  static_cast<PetscNonlinearSolver<Number> *> (ctx);
76 
77  // Get the current iteration number from the snes object,
78  // store it in the PetscNonlinearSolver object for possible use
79  // by the user's residual function.
80  {
81  PetscInt n_iterations = 0;
82  ierr = SNESGetIterationNumber(snes, &n_iterations);
83  CHKERRABORT(solver->comm().get(),ierr);
84  solver->_current_nonlinear_iteration_number = cast_int<unsigned>(n_iterations);
85  }
86 
87  NonlinearImplicitSystem & sys = solver->system();
88 
89  PetscVector<Number> & X_sys = *cast_ptr<PetscVector<Number> *>(sys.solution.get());
90 
91  PetscVector<Number> X_global(x, sys.comm());
92 
93  // Use the system's update() to get a good local version of the
94  // parallel solution. This operation does not modify the incoming
95  // "x" vector, it only localizes information from "x" into
96  // sys.current_local_solution.
97  X_global.swap(X_sys);
98  sys.update();
99  X_global.swap(X_sys);
100 
101  // Enforce constraints (if any) exactly on the
102  // current_local_solution. This is the solution vector that is
103  // actually used in the computation of the residual below, and is
104  // not locked by debug-enabled PETSc the way that "x" is.
106 
107  return ResidualContext(solver, sys, ierr);
108 }
109 
110 //--------------------------------------------------------------------
111 // Functions with C linkage to pass to PETSc. PETSc will call these
112 // methods as needed.
113 //
114 // Since they must have C linkage they have no knowledge of a namespace.
115 // Give them an obscure name to avoid namespace pollution.
116 extern "C"
117 {
118 
119  //-------------------------------------------------------------------
120  // this function is called by PETSc at the end of each nonlinear step
121  PetscErrorCode
122  libmesh_petsc_snes_monitor (SNES, PetscInt its, PetscReal fnorm, void *)
123  {
124  //PetscErrorCode ierr=0;
125 
126  //if (its > 0)
127  libMesh::out << " NL step "
128  << std::setw(2) << its
129  << std::scientific
130  << ", |residual|_2 = " << fnorm
131  << std::endl;
132 
133  //return ierr;
134  return 0;
135  }
136 
137 #ifdef LIBMESH_ENABLE_DEPRECATED
138  PetscErrorCode
139  __libmesh_petsc_snes_monitor (SNES, PetscInt its, PetscReal fnorm, void *)
140  {
141  libmesh_deprecated();
142  return libmesh_petsc_snes_monitor(nullptr, its, fnorm, nullptr);
143  }
144 #endif
145 
146 
147  //---------------------------------------------------------------
148  // this function is called by PETSc to evaluate the residual at X
149  PetscErrorCode
150  libmesh_petsc_snes_residual (SNES snes, Vec x, Vec r, void * ctx)
151  {
153 
154  libmesh_assert(r);
155  PetscVector<Number> R(r, rc.sys.comm());
156 
157  if (rc.solver->_zero_out_residual)
158  R.zero();
159 
160  //-----------------------------------------------------------------------------
161  // if the user has provided both function pointers and objects only the pointer
162  // will be used, so catch that as an error
163  if (rc.solver->residual && rc.solver->residual_object)
164  libmesh_error_msg("ERROR: cannot specify both a function and object to compute the Residual!");
165 
167  libmesh_error_msg("ERROR: cannot specify both a function and object to compute the combined Residual & Jacobian!");
168 
169  if (rc.solver->residual != nullptr)
170  rc.solver->residual(*rc.sys.current_local_solution.get(), R, rc.sys);
171 
172  else if (rc.solver->residual_object != nullptr)
174 
175  else if (rc.solver->matvec != nullptr)
176  rc.solver->matvec (*rc.sys.current_local_solution.get(), &R, nullptr, rc.sys);
177 
178  else if (rc.solver->residual_and_jacobian_object != nullptr)
180 
181  else
182  libmesh_error_msg("Error! Unable to compute residual and/or Jacobian!");
183 
184  R.close();
185 
186  return rc.ierr;
187  }
188 
189 #ifdef LIBMESH_ENABLE_DEPRECATED
190  PetscErrorCode
191  __libmesh_petsc_snes_residual (SNES snes, Vec x, Vec r, void * ctx)
192  {
193  libmesh_deprecated();
194  return libmesh_petsc_snes_residual(snes, x, r, ctx);
195  }
196 #endif
197 
198  //-----------------------------------------------------------------------------------------
199  // this function is called by PETSc to approximate the Jacobian at X via finite differences
200  PetscErrorCode
201  libmesh_petsc_snes_fd_residual (SNES snes, Vec x, Vec r, void * ctx)
202  {
204 
205  libmesh_assert(r);
206  PetscVector<Number> R(r, rc.sys.comm());
207 
208  if (rc.solver->_zero_out_residual)
209  R.zero();
210 
211  if (rc.solver->fd_residual_object != nullptr)
213 
214  else if (rc.solver->residual_object != nullptr)
216 
217  else
218  libmesh_error_msg("Error! Unable to compute residual for forming finite difference Jacobian!");
219 
220  R.close();
221 
222  return rc.ierr;
223  }
224 
225 #ifdef LIBMESH_ENABLE_DEPRECATED
226  PetscErrorCode
227  __libmesh_petsc_snes_fd_residual (SNES snes, Vec x, Vec r, void * ctx)
228  {
229  libmesh_deprecated();
230  return libmesh_petsc_snes_fd_residual(snes, x, r, ctx);
231  }
232 #endif
233 
234  //----------------------------------------------------------------
235  // this function is called by PETSc to approximate Jacobian-vector
236  // products at X via finite differences
237  PetscErrorCode
238  libmesh_petsc_snes_mffd_residual (SNES snes, Vec x, Vec r, void * ctx)
239  {
241 
242  libmesh_assert(r);
243  PetscVector<Number> R(r, rc.sys.comm());
244 
245  if (rc.solver->_zero_out_residual)
246  R.zero();
247 
248  if (rc.solver->mffd_residual_object != nullptr)
250 
251  else if (rc.solver->residual_object != nullptr)
253 
254  else
255  libmesh_error_msg("Error! Unable to compute residual for forming finite differenced"
256  "Jacobian-vector products!");
257 
258  R.close();
259 
260  return rc.ierr;
261  }
262 
263  //----------------------------------------------------------
264  // this function serves an interface between the petsc layer
265  // and the actual mffd residual computing routine
266  PetscErrorCode
267  libmesh_petsc_snes_mffd_interface (void * ctx, Vec x, Vec r)
268  {
269  // No way to safety-check this cast, since we got a void *...
271  static_cast<PetscNonlinearSolver<Number> *> (ctx);
272 
273  return libmesh_petsc_snes_mffd_residual(solver->snes(), x, r, ctx);
274  }
275 
276 #ifdef LIBMESH_ENABLE_DEPRECATED
277  PetscErrorCode
278  __libmesh_petsc_snes_mffd_interface (void * ctx, Vec x, Vec r)
279  {
280  libmesh_deprecated();
281  return libmesh_petsc_snes_mffd_interface(ctx, x, r);
282  }
283 #endif
284 
285  //---------------------------------------------------------------
286  // this function is called by PETSc to evaluate the Jacobian at X
287  PetscErrorCode
289 #if PETSC_RELEASE_LESS_THAN(3,5,0)
290  SNES snes, Vec x, Mat * jac, Mat * pc, MatStructure * msflag, void * ctx
291 #else
292  SNES snes, Vec x, Mat jac, Mat pc, void * ctx
293 #endif
294  )
295  {
296  LOG_SCOPE("jacobian()", "PetscNonlinearSolver");
297 
298  PetscErrorCode ierr=0;
299 
300  libmesh_assert(ctx);
301 
302  // No way to safety-check this cast, since we got a void *...
304  static_cast<PetscNonlinearSolver<Number> *> (ctx);
305 
306  // Get the current iteration number from the snes object,
307  // store it in the PetscNonlinearSolver object for possible use
308  // by the user's Jacobian function.
309  {
310  PetscInt n_iterations = 0;
311  ierr = SNESGetIterationNumber(snes, &n_iterations);
312  CHKERRABORT(solver->comm().get(),ierr);
313  solver->_current_nonlinear_iteration_number = cast_int<unsigned>(n_iterations);
314  }
315 
316  NonlinearImplicitSystem & sys = solver->system();
317 #if PETSC_RELEASE_LESS_THAN(3,5,0)
318  PetscMatrix<Number> PC(*pc, sys.comm());
319  PetscMatrix<Number> Jac(*jac, sys.comm());
320 #else
321  PetscMatrix<Number> PC(pc, sys.comm());
322  PetscMatrix<Number> Jac(jac, sys.comm());
323 #endif
324  PetscVector<Number> & X_sys = *cast_ptr<PetscVector<Number> *>(sys.solution.get());
325  PetscVector<Number> X_global(x, sys.comm());
326 
327  // Set the dof maps
328  PC.attach_dof_map(sys.get_dof_map());
329  Jac.attach_dof_map(sys.get_dof_map());
330 
331  // Use the systems update() to get a good local version of the parallel solution
332  X_global.swap(X_sys);
333  sys.update();
334  X_global.swap(X_sys);
335 
336  // Enforce constraints (if any) exactly on the
337  // current_local_solution. This is the solution vector that is
338  // actually used in the computation of the residual below, and is
339  // not locked by debug-enabled PETSc the way that "x" is.
341 
342  if (solver->_zero_out_jacobian)
343  PC.zero();
344 
345  //-----------------------------------------------------------------------------
346  // if the user has provided both function pointers and objects only the pointer
347  // will be used, so catch that as an error
348  if (solver->jacobian && solver->jacobian_object)
349  libmesh_error_msg("ERROR: cannot specify both a function and object to compute the Jacobian!");
350 
351  if (solver->matvec && solver->residual_and_jacobian_object)
352  libmesh_error_msg("ERROR: cannot specify both a function and object to compute the combined Residual & Jacobian!");
353 
354  if (solver->jacobian != nullptr)
355  solver->jacobian(*sys.current_local_solution.get(), PC, sys);
356 
357  else if (solver->jacobian_object != nullptr)
358  solver->jacobian_object->jacobian(*sys.current_local_solution.get(), PC, sys);
359 
360  else if (solver->matvec != nullptr)
361  solver->matvec(*sys.current_local_solution.get(), nullptr, &PC, sys);
362 
363  else if (solver->residual_and_jacobian_object != nullptr)
364  solver->residual_and_jacobian_object->residual_and_jacobian (*sys.current_local_solution.get(), nullptr, &PC, sys);
365 
366  else
367  libmesh_error_msg("Error! Unable to compute residual and/or Jacobian!");
368 
369  PC.close();
370  Jac.close();
371 #if PETSC_RELEASE_LESS_THAN(3,5,0)
372  *msflag = SAME_NONZERO_PATTERN;
373 #endif
374 
375  return ierr;
376  }
377 
378  // This function gets called by PETSc in place of the standard Petsc line searches
379  // if a linesearch object is supplied to the PetscNonlinearSolver class. It wraps
380  // the lineserach algorithm implemented on the linesearch object.
381  // * "linesearch" is an object that can be used to access the non-linear and linear solution
382  // vectors as well as the residual and SNES object
383  // * "ctx" is the PetscNonlinearSolver context
384  PetscErrorCode libmesh_petsc_linesearch_shellfunc (SNESLineSearch linesearch, void * ctx)
385  {
386  // No way to safety-check this cast, since we got a void *...
388  static_cast<PetscNonlinearSolver<Number> *> (ctx);
389 
390  solver->linesearch_object->linesearch(linesearch);
391  return 0;
392  }
393 
394 #ifdef LIBMESH_ENABLE_DEPRECATED
395  PetscErrorCode
397 #if PETSC_RELEASE_LESS_THAN(3,5,0)
398  SNES snes, Vec x, Mat * jac, Mat * pc, MatStructure * msflag, void * ctx
399 #else
400  SNES snes, Vec x, Mat jac, Mat pc, void * ctx
401 #endif
402  )
403  {
404  libmesh_deprecated();
406 #if PETSC_RELEASE_LESS_THAN(3,5,0)
407  snes, x, jac, pc, msflag, ctx
408 #else
409  snes, x, jac, pc, ctx
410 #endif
411  );
412  }
413 #endif
414 
415  // This function gets called by PETSc after the SNES linesearch is
416  // complete. We use it to exactly enforce any constraints on the
417  // solution which may have drifted during the linear solve. In the
418  // PETSc nomenclature:
419  // * "x" is the old solution vector,
420  // * "y" is the search direction (Newton step) vector,
421  // * "w" is the candidate solution vector, and
422  // the user is responsible for setting changed_y and changed_w
423  // appropriately, depending on whether or not the search
424  // direction or solution vector was changed, respectively.
426 #if PETSC_VERSION_LESS_THAN(3,3,0)
427  SNES, Vec x, Vec y, Vec w, void * context, PetscBool * changed_y, PetscBool * changed_w
428 #else
429  SNESLineSearch, Vec x, Vec y, Vec w, PetscBool * changed_y, PetscBool * changed_w, void * context
430 #endif
431  )
432  {
433  LOG_SCOPE("postcheck()", "PetscNonlinearSolver");
434 
435  PetscErrorCode ierr = 0;
436 
437  // PETSc almost certainly initializes these to false already, but
438  // it doesn't hurt to be explicit.
439  *changed_w = PETSC_FALSE;
440  *changed_y = PETSC_FALSE;
441 
442  libmesh_assert(context);
443 
444  // Cast the context to a NonlinearSolver object.
446  static_cast<PetscNonlinearSolver<Number> *> (context);
447 
448  // If the user has provided both postcheck function pointer and
449  // object, this is ambiguous, so throw an error.
450  if (solver->postcheck && solver->postcheck_object)
451  libmesh_error_msg("ERROR: cannot specify both a function and object for performing the solve postcheck!");
452 
453  // It's also possible that we don't need to do anything at all, in
454  // that case return early...
455  NonlinearImplicitSystem & sys = solver->system();
456  DofMap & dof_map = sys.get_dof_map();
457 
458  if (!dof_map.n_constrained_dofs() && !solver->postcheck && !solver->postcheck_object)
459  return ierr;
460 
461  // We definitely need to wrap at least "w"
462  PetscVector<Number> petsc_w(w, sys.comm());
463 
464  // The user sets these flags in his/her postcheck function to
465  // indicate whether they changed something.
466  bool
467  changed_search_direction = false,
468  changed_new_soln = false;
469 
470  if (solver->postcheck || solver->postcheck_object)
471  {
472  PetscVector<Number> petsc_x(x, sys.comm());
473  PetscVector<Number> petsc_y(y, sys.comm());
474 
475  if (solver->postcheck)
476  solver->postcheck(petsc_x,
477  petsc_y,
478  petsc_w,
479  changed_search_direction,
480  changed_new_soln,
481  sys);
482 
483  else if (solver->postcheck_object)
484  solver->postcheck_object->postcheck(petsc_x,
485  petsc_y,
486  petsc_w,
487  changed_search_direction,
488  changed_new_soln,
489  sys);
490  }
491 
492  // Record whether the user changed the solution or the search direction.
493  if (changed_search_direction)
494  *changed_y = PETSC_TRUE;
495 
496  if (changed_new_soln)
497  *changed_w = PETSC_TRUE;
498 
499  if (dof_map.n_constrained_dofs())
500  {
501  PetscVector<Number> & system_soln = *cast_ptr<PetscVector<Number> *>(sys.solution.get());
502 
503  // ... and swap it in before enforcing the constraints.
504  petsc_w.swap(system_soln);
505 
506  dof_map.enforce_constraints_exactly(sys);
507 
508  // If we have constraints, we'll assume that we did change the
509  // solution w (hopefully slightly). Enforcing constraints
510  // does not change the search direction, y, but the user may
511  // have, so we leave it alone.
512  *changed_w = PETSC_TRUE;
513 
514  // Swap back
515  petsc_w.swap(system_soln);
516  }
517 
518  return ierr;
519  }
520 
521 #ifdef LIBMESH_ENABLE_DEPRECATED
523 #if PETSC_VERSION_LESS_THAN(3,3,0)
524  SNES, Vec x, Vec y, Vec w, void * context, PetscBool * changed_y, PetscBool * changed_w
525 #else
526  SNESLineSearch, Vec x, Vec y, Vec w, PetscBool * changed_y, PetscBool * changed_w, void * context
527 #endif
528  )
529  {
530  libmesh_deprecated();
532 #if PETSC_VERSION_LESS_THAN(3,3,0)
533  nullptr, x, y, w, context, changed_y, changed_w
534 #else
535  nullptr, x, y, w, changed_y, changed_w, context
536 #endif
537  );
538  }
539 #endif
540 } // end extern "C"
541 
542 
543 
544 //---------------------------------------------------------------------
545 // PetscNonlinearSolver<> methods
546 template <typename T>
548  NonlinearSolver<T>(system_in),
549  linesearch_object(nullptr),
550  _reason(SNES_CONVERGED_ITERATING/*==0*/), // Arbitrary initial value...
551  _n_linear_iterations(0),
552  _current_nonlinear_iteration_number(0),
553  _zero_out_residual(true),
554  _zero_out_jacobian(true),
555  _default_monitor(true),
556  _snesmf_reuse_base(true)
557 {
558 }
559 
560 
561 
562 template <typename T>
564 {
565  this->clear ();
566 }
567 
568 
569 
570 template <typename T>
572 {
573  if (this->initialized())
574  {
575  this->_is_initialized = false;
576 
577  PetscErrorCode ierr=0;
578 
579  ierr = LibMeshSNESDestroy(&_snes);
580  LIBMESH_CHKERR(ierr);
581 
582  // Reset the nonlinear iteration counter. This information is only relevant
583  // *during* the solve(). After the solve is completed it should return to
584  // the default value of 0.
585  _current_nonlinear_iteration_number = 0;
586  }
587 }
588 
589 
590 
591 template <typename T>
593 {
594  // Initialize the data structures if not done so already.
595  if (!this->initialized())
596  {
597  this->_is_initialized = true;
598 
599  PetscErrorCode ierr=0;
600 
601  ierr = SNESCreate(this->comm().get(),&_snes);
602  LIBMESH_CHKERR(ierr);
603 
604  if (name)
605  {
606  ierr = SNESSetOptionsPrefix(_snes, name);
607  LIBMESH_CHKERR(ierr);
608  }
609 
610 #if !PETSC_RELEASE_LESS_THAN(3,3,0)
611  // Attaching a DM to SNES.
612  DM dm;
613  ierr = DMCreate(this->comm().get(), &dm);LIBMESH_CHKERR(ierr);
614  ierr = DMSetType(dm,DMLIBMESH);LIBMESH_CHKERR(ierr);
615  ierr = DMlibMeshSetSystem(dm,this->system());LIBMESH_CHKERR(ierr);
616  if (name)
617  {
618  ierr = DMSetOptionsPrefix(dm,name); LIBMESH_CHKERR(ierr);
619  }
620  ierr = DMSetFromOptions(dm); LIBMESH_CHKERR(ierr);
621  ierr = DMSetUp(dm); LIBMESH_CHKERR(ierr);
622  ierr = SNESSetDM(this->_snes, dm); LIBMESH_CHKERR(ierr);
623  // SNES now owns the reference to dm.
624  ierr = DMDestroy(&dm); LIBMESH_CHKERR(ierr);
625 
626 #endif
627 
628  if (_default_monitor)
629  {
630  ierr = SNESMonitorSet (_snes, libmesh_petsc_snes_monitor,
631  this, PETSC_NULL);
632  LIBMESH_CHKERR(ierr);
633  }
634 
635  // If the SolverConfiguration object is provided, use it to set
636  // options during solver initialization.
637  if (this->_solver_configuration)
638  {
639  this->_solver_configuration->set_options_during_init();
640  }
641 
642 #if PETSC_VERSION_LESS_THAN(3,1,0)
643  // Cannot call SNESSetOptions before SNESSetFunction when using
644  // any matrix free options with PETSc 3.1.0+
645  ierr = SNESSetFromOptions(_snes);
646  LIBMESH_CHKERR(ierr);
647 #endif
648 
649  if (this->_preconditioner)
650  {
651  KSP ksp;
652  ierr = SNESGetKSP (_snes, &ksp);
653  LIBMESH_CHKERR(ierr);
654  PC pc;
655  ierr = KSPGetPC(ksp,&pc);
656  LIBMESH_CHKERR(ierr);
657 
658  this->_preconditioner->init();
659 
660  PCSetType(pc, PCSHELL);
661  PCShellSetContext(pc,(void *)this->_preconditioner);
662 
663  //Re-Use the shell functions from petsc_linear_solver
664  PCShellSetSetUp(pc,libmesh_petsc_preconditioner_setup);
665  PCShellSetApply(pc,libmesh_petsc_preconditioner_apply);
666  }
667  }
668 
669 
670  // Tell PETSc about our linesearch "post-check" function, but only
671  // if the user has provided one. There seem to be extra,
672  // unnecessary residual calculations if a postcheck function is
673  // attached for no reason.
674  if (this->postcheck || this->postcheck_object)
675  {
676 #if PETSC_VERSION_LESS_THAN(3,3,0)
677  PetscErrorCode ierr = SNESLineSearchSetPostCheck(_snes, libmesh_petsc_snes_postcheck, this);
678  LIBMESH_CHKERR(ierr);
679 
680 #else
681  SNESLineSearch linesearch;
682  PetscErrorCode ierr = SNESGETLINESEARCH(_snes, &linesearch);
683  LIBMESH_CHKERR(ierr);
684 
685  ierr = SNESLineSearchSetPostCheck(linesearch, libmesh_petsc_snes_postcheck, this);
686  LIBMESH_CHKERR(ierr);
687 #endif
688  }
689 }
690 
691 #if !PETSC_VERSION_LESS_THAN(3,3,0)
692 template <typename T>
693 void
695  void (*computeSubspace)(std::vector<NumericVector<Number> *> &, sys_type &),
696  MatNullSpace * msp)
697 {
698  PetscErrorCode ierr;
699  std::vector<NumericVector<Number> *> sp;
700  if (computeSubspaceObject)
701  (*computeSubspaceObject)(sp, this->system());
702  else
703  (*computeSubspace)(sp, this->system());
704 
705  *msp = PETSC_NULL;
706  if (sp.size())
707  {
708  Vec * modes;
709  PetscScalar * dots;
710  PetscInt nmodes = cast_int<PetscInt>(sp.size());
711 
712 #if PETSC_RELEASE_LESS_THAN(3,5,0)
713  ierr = PetscMalloc2(nmodes,Vec,&modes,nmodes,PetscScalar,&dots);
714 #else
715  ierr = PetscMalloc2(nmodes,&modes,nmodes,&dots);
716 #endif
717  LIBMESH_CHKERR(ierr);
718 
719  for (PetscInt i=0; i<nmodes; ++i)
720  {
721  PetscVector<T> * pv = cast_ptr<PetscVector<T> *>(sp[i]);
722  Vec v = pv->vec();
723 
724  ierr = VecDuplicate(v, modes+i);
725  LIBMESH_CHKERR(ierr);
726 
727  ierr = VecCopy(v,modes[i]);
728  LIBMESH_CHKERR(ierr);
729  }
730 
731  // Normalize.
732  ierr = VecNormalize(modes[0],PETSC_NULL);
733  LIBMESH_CHKERR(ierr);
734 
735  for (PetscInt i=1; i<nmodes; i++)
736  {
737  // Orthonormalize vec[i] against vec[0:i-1]
738  ierr = VecMDot(modes[i],i,modes,dots);
739  LIBMESH_CHKERR(ierr);
740 
741  for (PetscInt j=0; j<i; j++)
742  dots[j] *= -1.;
743 
744  ierr = VecMAXPY(modes[i],i,dots,modes);
745  LIBMESH_CHKERR(ierr);
746 
747  ierr = VecNormalize(modes[i],PETSC_NULL);
748  LIBMESH_CHKERR(ierr);
749  }
750 
751  ierr = MatNullSpaceCreate(this->comm().get(), PETSC_FALSE, nmodes, modes, msp);
752  LIBMESH_CHKERR(ierr);
753 
754  for (PetscInt i=0; i<nmodes; ++i)
755  {
756  ierr = VecDestroy(modes+i);
757  LIBMESH_CHKERR(ierr);
758  }
759 
760  ierr = PetscFree2(modes,dots);
761  LIBMESH_CHKERR(ierr);
762  }
763 }
764 #endif
765 
766 template <typename T>
767 std::pair<unsigned int, Real>
768 PetscNonlinearSolver<T>::solve (SparseMatrix<T> & pre_in, // System Preconditioning Matrix
769  NumericVector<T> & x_in, // Solution vector
770  NumericVector<T> & r_in, // Residual vector
771  const double, // Stopping tolerance
772  const unsigned int)
773 {
774  LOG_SCOPE("solve()", "PetscNonlinearSolver");
775  this->init ();
776 
777  // Make sure the data passed in are really of Petsc types
778  PetscMatrix<T> * pre = cast_ptr<PetscMatrix<T> *>(&pre_in);
779  PetscVector<T> * x = cast_ptr<PetscVector<T> *>(&x_in);
780  PetscVector<T> * r = cast_ptr<PetscVector<T> *>(&r_in);
781 
782  PetscErrorCode ierr=0;
783  PetscInt n_iterations =0;
784  // Should actually be a PetscReal, but I don't know which version of PETSc first introduced PetscReal
785  Real final_residual_norm=0.;
786 
787  ierr = SNESSetFunction (_snes, r->vec(), libmesh_petsc_snes_residual, this);
788  LIBMESH_CHKERR(ierr);
789 
790  // Only set the jacobian function if we've been provided with something to call.
791  // This allows a user to set their own jacobian function if they want to
792  if (this->jacobian || this->jacobian_object || this->residual_and_jacobian_object)
793  {
794  ierr = SNESSetJacobian (_snes, pre->mat(), pre->mat(), libmesh_petsc_snes_jacobian, this);
795  LIBMESH_CHKERR(ierr);
796  }
797 
798 
799 #if !PETSC_VERSION_LESS_THAN(3,3,0)
800  // Only set the nullspace if we have a way of computing it and the result is non-empty.
801  if (this->nullspace || this->nullspace_object)
802  {
803  MatNullSpace msp;
804  this->build_mat_null_space(this->nullspace_object, this->nullspace, &msp);
805  if (msp)
806  {
807  ierr = MatSetNullSpace(pre->mat(), msp);
808  LIBMESH_CHKERR(ierr);
809  ierr = MatNullSpaceDestroy(&msp);
810  LIBMESH_CHKERR(ierr);
811  }
812  }
813 
814  // Only set the transpose nullspace if we have a way of computing it and the result is non-empty.
815  if (this->transpose_nullspace || this->transpose_nullspace_object)
816  {
817 #if PETSC_VERSION_LESS_THAN(3,6,0)
818  libmesh_warning("MatSetTransposeNullSpace is only supported for PETSc >= 3.6, transpose nullspace will be ignored.");
819 #else
820  MatNullSpace msp = PETSC_NULL;
821  this->build_mat_null_space(this->transpose_nullspace_object, this->transpose_nullspace, &msp);
822  if (msp)
823  {
824  ierr = MatSetTransposeNullSpace(pre->mat(), msp);
825  LIBMESH_CHKERR(ierr);
826  ierr = MatNullSpaceDestroy(&msp);
827  LIBMESH_CHKERR(ierr);
828  }
829 #endif
830  }
831 
832  // Only set the nearnullspace if we have a way of computing it and the result is non-empty.
833  if (this->nearnullspace || this->nearnullspace_object)
834  {
835  MatNullSpace msp = PETSC_NULL;
836  this->build_mat_null_space(this->nearnullspace_object, this->nearnullspace, &msp);
837 
838  if (msp)
839  {
840  ierr = MatSetNearNullSpace(pre->mat(), msp);
841  LIBMESH_CHKERR(ierr);
842  ierr = MatNullSpaceDestroy(&msp);
843  LIBMESH_CHKERR(ierr);
844  }
845  }
846 #endif
847 
848  // Have the Krylov subspace method use our good initial guess rather than 0
849  KSP ksp;
850  ierr = SNESGetKSP (_snes, &ksp);
851  LIBMESH_CHKERR(ierr);
852 
853  // Set the tolerances for the iterative solver. Use the user-supplied
854  // tolerance for the relative residual & leave the others at default values
855  ierr = KSPSetTolerances (ksp, this->initial_linear_tolerance, PETSC_DEFAULT,
856  PETSC_DEFAULT, this->max_linear_iterations);
857  LIBMESH_CHKERR(ierr);
858 
859  // Set the tolerances for the non-linear solver.
860  ierr = SNESSetTolerances(_snes, this->absolute_residual_tolerance, this->relative_residual_tolerance,
861  this->relative_step_tolerance, this->max_nonlinear_iterations, this->max_function_evaluations);
862  LIBMESH_CHKERR(ierr);
863 
864  //Pull in command-line options
865 #if PETSC_VERSION_LESS_THAN(3,7,0)
866  KSPSetFromOptions(ksp);
867 #endif
868  SNESSetFromOptions(_snes);
869 
870  if (this->user_presolve)
871  this->user_presolve(this->system());
872 
873  //Set the preconditioning matrix
874  if (this->_preconditioner)
875  {
876  this->_preconditioner->set_matrix(pre_in);
877  this->_preconditioner->init();
878  }
879 
880  // If the SolverConfiguration object is provided, use it to override
881  // solver options.
882  if (this->_solver_configuration)
883  {
884  this->_solver_configuration->configure_solver();
885  }
886 
887  // In PETSc versions before 3.5.0, it is not possible to call
888  // SNESSetUp() before the solution and rhs vectors are initialized, as
889  // this triggers the
890  //
891  // "Solution vector cannot be right hand side vector!"
892  //
893  // error message. It is also not possible to call SNESSetSolution()
894  // in those versions of PETSc to work around the problem, since that
895  // API was removed in 3.0.0 and only restored in 3.6.0. The
896  // overzealous check was moved out of SNESSetUp in PETSc 3.5.0
897  // (petsc/petsc@154060b), so this code block should be safe to use
898  // in 3.5.0 and later.
899 #if !PETSC_VERSION_LESS_THAN(3,5,0)
900 #if !PETSC_VERSION_LESS_THAN(3,6,0)
901  ierr = SNESSetSolution(_snes, x->vec());
902  LIBMESH_CHKERR(ierr);
903 #endif
904  ierr = SNESSetUp(_snes);
905  LIBMESH_CHKERR(ierr);
906 
907  Mat J;
908  ierr = SNESGetJacobian(_snes,&J,PETSC_NULL,PETSC_NULL,PETSC_NULL);
909  LIBMESH_CHKERR(ierr);
910  ierr = MatMFFDSetFunction(J, libmesh_petsc_snes_mffd_interface, this);
911  LIBMESH_CHKERR(ierr);
912 #if !PETSC_VERSION_LESS_THAN(3, 8, 4)
913  // Resue the residual vector from SNES
914  ierr = MatSNESMFSetReuseBase(J, static_cast<PetscBool>(_snesmf_reuse_base));
915  LIBMESH_CHKERR(ierr);
916 #endif
917 #endif
918 
919 #if PETSC_VERSION_LESS_THAN(3, 3, 0)
920  if (linesearch_object)
921  libmesh_error_msg("Line search setter interface introduced in petsc version 3.3!");
922 #else
923  SNESLineSearch linesearch;
924  if (linesearch_object)
925  {
926  ierr = SNESGETLINESEARCH(_snes, &linesearch);
927  LIBMESH_CHKERR(ierr);
928  ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHSHELL);
929  LIBMESH_CHKERR(ierr);
930  ierr = SNESLineSearchShellSetUserFunc(linesearch, libmesh_petsc_linesearch_shellfunc, this);
931  LIBMESH_CHKERR(ierr);
932  }
933 #endif
934 
935  ierr = SNESSolve (_snes, PETSC_NULL, x->vec());
936  LIBMESH_CHKERR(ierr);
937 
938  ierr = SNESGetIterationNumber(_snes,&n_iterations);
939  LIBMESH_CHKERR(ierr);
940 
941  ierr = SNESGetLinearSolveIterations(_snes, &_n_linear_iterations);
942  LIBMESH_CHKERR(ierr);
943 
944  // Enforce constraints exactly now that the solve is done. We have
945  // been enforcing them on the current_local_solution during the
946  // solve, but now need to be sure they are enforced on the parallel
947  // solution vector as well.
948  this->system().get_dof_map().enforce_constraints_exactly(this->system());
949 
950  // SNESGetFunction has been around forever and should work on all
951  // versions of PETSc. This is also now the recommended approach
952  // according to the documentation for the PETSc 3.5.1 release:
953  // http://www.mcs.anl.gov/petsc/documentation/changes/35.html
954  Vec f;
955  ierr = SNESGetFunction(_snes, &f, 0, 0);
956  LIBMESH_CHKERR(ierr);
957  ierr = VecNorm(f, NORM_2, &final_residual_norm);
958  LIBMESH_CHKERR(ierr);
959 
960  // Get and store the reason for convergence
961  SNESGetConvergedReason(_snes, &_reason);
962 
963  //Based on Petsc 2.3.3 documentation all diverged reasons are negative
964  this->converged = (_reason >= 0);
965 
966  this->clear();
967 
968  // return the # of its. and the final residual norm.
969  return std::make_pair(n_iterations, final_residual_norm);
970 }
971 
972 
973 
974 template <typename T>
976 {
977 
978  libMesh::out << "Nonlinear solver convergence/divergence reason: "
979  << SNESConvergedReasons[this->get_converged_reason()] << std::endl;
980 }
981 
982 
983 
984 template <typename T>
986 {
987  PetscErrorCode ierr=0;
988 
989  if (this->initialized())
990  {
991  ierr = SNESGetConvergedReason(_snes, &_reason);
992  LIBMESH_CHKERR(ierr);
993  }
994 
995  return _reason;
996 }
997 
998 template <typename T>
1000 {
1001  return _n_linear_iterations;
1002 }
1003 
1004 
1005 //------------------------------------------------------------------
1006 // Explicit instantiations
1007 template class PetscNonlinearSolver<Number>;
1008 
1009 } // namespace libMesh
1010 
1011 
1012 
1013 #endif // #ifdef LIBMESH_HAVE_PETSC
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
PetscNonlinearSolver< Number > * solver
NumericVector interface to PETSc Vec.
Definition: petsc_vector.h:64
PetscErrorCode libmesh_petsc_snes_postcheck(#if PETSC_VERSION_LESS_THAN(3, 3, 0) SNES, Vec x, Vec y, Vec w, void *context, PetscBool *changed_y, PetscBool *changed_w #else SNESLineSearch, Vec x, Vec y, Vec w, PetscBool *changed_y, PetscBool *changed_w, void *context #endif)
void(* residual)(const NumericVector< Number > &X, NumericVector< Number > &R, sys_type &S)
NonlinearImplicitSystem::ComputeResidualandJacobian * residual_and_jacobian_object
virtual void residual(const NumericVector< Number > &X, NumericVector< Number > &R, sys_type &S)=0
virtual void zero() override
Definition: petsc_vector.h:862
PetscErrorCode __libmesh_petsc_snes_residual(SNES snes, Vec x, Vec r, void *ctx)
NonlinearImplicitSystem & sys
PetscErrorCode __libmesh_petsc_snes_monitor(SNES, PetscInt its, PetscReal fnorm, void *)
const Parallel::Communicator & comm() const
NonlinearImplicitSystem::ComputeResidual * residual_object
void(* jacobian)(const NumericVector< Number > &X, SparseMatrix< Number > &J, sys_type &S)
PetscErrorCode libmesh_petsc_snes_monitor(SNES, PetscInt its, PetscReal fnorm, void *)
PetscErrorCode __libmesh_petsc_snes_mffd_interface(void *ctx, Vec x, Vec r)
virtual void postcheck(const NumericVector< Number > &old_soln, NumericVector< Number > &search_direction, NumericVector< Number > &new_soln, bool &changed_search_direction, bool &changed_new_soln, sys_type &S)=0
virtual void swap(NumericVector< T > &v) override
PetscErrorCode libmesh_petsc_snes_mffd_residual(SNES snes, Vec x, Vec r, void *ctx)
Manages the degrees of freedom (DOFs) in a simulation.
Definition: dof_map.h:176
PetscErrorCode __libmesh_petsc_snes_jacobian(#if PETSC_RELEASE_LESS_THAN(3, 5, 0) SNES snes, Vec x, Mat *jac, Mat *pc, MatStructure *msflag, void *ctx #else SNES snes, Vec x, Mat jac, Mat pc, void *ctx #endif)
const sys_type & system() const
NonlinearImplicitSystem::ComputeJacobian * jacobian_object
PetscErrorCode libmesh_petsc_snes_residual(SNES snes, Vec x, Vec r, void *ctx)
PetscErrorCode libmesh_petsc_preconditioner_setup(void *ctx)
ResidualContext(PetscNonlinearSolver< Number > *solver_in, NonlinearImplicitSystem &sys_in, PetscErrorCode ierr_in)
virtual void jacobian(const NumericVector< Number > &X, SparseMatrix< Number > &J, sys_type &S)=0
std::unique_ptr< NumericVector< Number > > solution
Definition: system.h:1523
void init(triangulateio &t)
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and non-linear solv...
ResidualContext libmesh_petsc_snes_residual_helper(SNES snes, Vec x, void *ctx)
PetscErrorCode libmesh_petsc_snes_jacobian(#if PETSC_RELEASE_LESS_THAN(3, 5, 0) SNES snes, Vec x, Mat *jac, Mat *pc, MatStructure *msflag, void *ctx #else SNES snes, Vec x, Mat jac, Mat pc, void *ctx #endif)
PetscErrorCode libmesh_petsc_snes_fd_residual(SNES snes, Vec x, Vec r, void *ctx)
PetscErrorCode libmesh_petsc_snes_mffd_interface(void *ctx, Vec x, Vec r)
PetscErrorCode DMlibMeshSetSystem(DM dm, libMesh::NonlinearImplicitSystem &sys)
virtual void update()
Definition: system.C:408
PetscTruth PetscBool
Definition: petsc_macro.h:67
PetscErrorCode libmesh_petsc_linesearch_shellfunc(SNESLineSearch linesearch, void *ctx)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::unique_ptr< ComputeLineSearchObject > linesearch_object
NonlinearImplicitSystem::ComputePostCheck * postcheck_object
dof_id_type n_constrained_dofs() const
PetscErrorCode ierr
PetscErrorCode __libmesh_petsc_snes_postcheck(#if PETSC_VERSION_LESS_THAN(3, 3, 0) SNES, Vec x, Vec y, Vec w, void *context, PetscBool *changed_y, PetscBool *changed_w #else SNESLineSearch, Vec x, Vec y, Vec w, PetscBool *changed_y, PetscBool *changed_w, void *context #endif)
PetscErrorCode libmesh_petsc_preconditioner_apply(void *ctx, Vec x, Vec y)
NonlinearImplicitSystem::ComputeResidual * mffd_residual_object
SparseMatrix interface to PETSc Mat.
std::unique_ptr< NumericVector< Number > > current_local_solution
Definition: system.h:1535
bool initialized()
Definition: libmesh.C:258
void(* postcheck)(const NumericVector< Number > &old_soln, NumericVector< Number > &search_direction, NumericVector< Number > &new_soln, bool &changed_search_direction, bool &changed_new_soln, sys_type &S)
OStreamProxy out(std::cout)
const DofMap & get_dof_map() const
Definition: system.h:2049
virtual void residual_and_jacobian(const NumericVector< Number > &X, NumericVector< Number > *R, SparseMatrix< Number > *J, sys_type &S)=0
PetscErrorCode __libmesh_petsc_snes_fd_residual(SNES snes, Vec x, Vec r, void *ctx)
void(* matvec)(const NumericVector< Number > &X, NumericVector< Number > *R, SparseMatrix< Number > *J, sys_type &S)
NonlinearImplicitSystem::ComputeResidual * fd_residual_object
void enforce_constraints_exactly(const System &system, NumericVector< Number > *v=nullptr, bool homogeneous=false) const