$extrastylesheet
petsc_nonlinear_solver.C
Go to the documentation of this file.
00001 // The libMesh Finite Element Library.
00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
00003 
00004 // This library is free software; you can redistribute it and/or
00005 // modify it under the terms of the GNU Lesser General Public
00006 // License as published by the Free Software Foundation; either
00007 // version 2.1 of the License, or (at your option) any later version.
00008 
00009 // This library is distributed in the hope that it will be useful,
00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00012 // Lesser General Public License for more details.
00013 
00014 // You should have received a copy of the GNU Lesser General Public
00015 // License along with this library; if not, write to the Free Software
00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00017 
00018 
00019 
00020 #include "libmesh/libmesh_common.h"
00021 
00022 #ifdef LIBMESH_HAVE_PETSC
00023 
00024 
00025 // C++ includes
00026 
00027 // Local Includes
00028 #include "libmesh/nonlinear_implicit_system.h"
00029 #include "libmesh/petsc_nonlinear_solver.h"
00030 #include "libmesh/petsc_linear_solver.h"
00031 #include "libmesh/petsc_vector.h"
00032 #include "libmesh/petsc_matrix.h"
00033 #include "libmesh/dof_map.h"
00034 #include "libmesh/preconditioner.h"
00035 
00036 /* DMlibMesh include. */
00037 #include "libmesh/petscdmlibmesh.h"
00038 
00039 namespace libMesh
00040 {
00041 
00042 //--------------------------------------------------------------------
00043 // Functions with C linkage to pass to PETSc.  PETSc will call these
00044 // methods as needed.
00045 //
00046 // Since they must have C linkage they have no knowledge of a namespace.
00047 // Give them an obscure name to avoid namespace pollution.
00048 extern "C"
00049 {
00050   // Older versions of PETSc do not have the different int typedefs.
00051   // On 64-bit machines, PetscInt may actually be a long long int.
00052   // This change occurred in Petsc-2.2.1.
00053 #if PETSC_VERSION_LESS_THAN(2,2,1)
00054   typedef int PetscErrorCode;
00055   typedef int PetscInt;
00056 #endif
00057 
00058   //-------------------------------------------------------------------
00059   // this function is called by PETSc at the end of each nonlinear step
00060   PetscErrorCode
00061   __libmesh_petsc_snes_monitor (SNES, PetscInt its, PetscReal fnorm, void *)
00062   {
00063     //PetscErrorCode ierr=0;
00064 
00065     //if (its > 0)
00066     libMesh::out << "  NL step "
00067                  << std::setw(2) << its
00068                  << std::scientific
00069                  << ", |residual|_2 = " << fnorm
00070                  << std::endl;
00071 
00072     //return ierr;
00073     return 0;
00074   }
00075 
00076 
00077 
00078   //---------------------------------------------------------------
00079   // this function is called by PETSc to evaluate the residual at X
00080   PetscErrorCode
00081   __libmesh_petsc_snes_residual (SNES snes, Vec x, Vec r, void *ctx)
00082   {
00083     START_LOG("residual()", "PetscNonlinearSolver");
00084 
00085     PetscErrorCode ierr=0;
00086 
00087     libmesh_assert(x);
00088     libmesh_assert(r);
00089     libmesh_assert(ctx);
00090 
00091     // No way to safety-check this cast, since we got a void*...
00092     PetscNonlinearSolver<Number>* solver =
00093       static_cast<PetscNonlinearSolver<Number>*> (ctx);
00094 
00095     // Get the current iteration number from the snes object,
00096     // store it in the PetscNonlinearSolver object for possible use
00097     // by the user's residual function.
00098     {
00099       PetscInt n_iterations = 0;
00100       ierr = SNESGetIterationNumber(snes, &n_iterations);
00101       CHKERRABORT(solver->comm().get(),ierr);
00102       solver->_current_nonlinear_iteration_number = cast_int<unsigned>(n_iterations);
00103     }
00104 
00105     NonlinearImplicitSystem &sys = solver->system();
00106 
00107     PetscVector<Number>& X_sys = *cast_ptr<PetscVector<Number>*>(sys.solution.get());
00108     PetscVector<Number>& R_sys = *cast_ptr<PetscVector<Number>*>(sys.rhs);
00109     PetscVector<Number> X_global(x, sys.comm()), R(r, sys.comm());
00110 
00111     // Use the systems update() to get a good local version of the parallel solution
00112     X_global.swap(X_sys);
00113     R.swap(R_sys);
00114 
00115     sys.get_dof_map().enforce_constraints_exactly(sys);
00116 
00117     sys.update();
00118 
00119     //Swap back
00120     X_global.swap(X_sys);
00121     R.swap(R_sys);
00122 
00123     if (solver->_zero_out_residual)
00124       R.zero();
00125 
00126     //-----------------------------------------------------------------------------
00127     // if the user has provided both function pointers and objects only the pointer
00128     // will be used, so catch that as an error
00129     if (solver->residual && solver->residual_object)
00130       libmesh_error_msg("ERROR: cannot specifiy both a function and object to compute the Residual!");
00131 
00132     if (solver->matvec && solver->residual_and_jacobian_object)
00133       libmesh_error_msg("ERROR: cannot specifiy both a function and object to compute the combined Residual & Jacobian!");
00134 
00135     if (solver->residual != NULL)
00136       solver->residual(*sys.current_local_solution.get(), R, sys);
00137 
00138     else if (solver->residual_object != NULL)
00139       solver->residual_object->residual(*sys.current_local_solution.get(), R, sys);
00140 
00141     else if (solver->matvec != NULL)
00142       solver->matvec (*sys.current_local_solution.get(), &R, NULL, sys);
00143 
00144     else if (solver->residual_and_jacobian_object != NULL)
00145       solver->residual_and_jacobian_object->residual_and_jacobian (*sys.current_local_solution.get(), &R, NULL, sys);
00146 
00147     else
00148       libmesh_error_msg("Error! Unable to compute residual and/or Jacobian!");
00149 
00150     R.close();
00151     X_global.close();
00152 
00153     STOP_LOG("residual()", "PetscNonlinearSolver");
00154 
00155     return ierr;
00156   }
00157 
00158 
00159 
00160   //---------------------------------------------------------------
00161   // this function is called by PETSc to evaluate the Jacobian at X
00162 #if PETSC_RELEASE_LESS_THAN(3,5,0)
00163   PetscErrorCode
00164   __libmesh_petsc_snes_jacobian (SNES snes, Vec x, Mat *jac, Mat *pc, MatStructure *msflag, void *ctx)
00165 #else
00166   PetscErrorCode
00167   __libmesh_petsc_snes_jacobian (SNES snes, Vec x, Mat jac, Mat pc, void *ctx)
00168 #endif
00169   {
00170     START_LOG("jacobian()", "PetscNonlinearSolver");
00171 
00172     PetscErrorCode ierr=0;
00173 
00174     libmesh_assert(ctx);
00175 
00176     // No way to safety-check this cast, since we got a void*...
00177     PetscNonlinearSolver<Number>* solver =
00178       static_cast<PetscNonlinearSolver<Number>*> (ctx);
00179 
00180     // Get the current iteration number from the snes object,
00181     // store it in the PetscNonlinearSolver object for possible use
00182     // by the user's Jacobian function.
00183     {
00184       PetscInt n_iterations = 0;
00185       ierr = SNESGetIterationNumber(snes, &n_iterations);
00186       CHKERRABORT(solver->comm().get(),ierr);
00187       solver->_current_nonlinear_iteration_number = cast_int<unsigned>(n_iterations);
00188     }
00189 
00190     NonlinearImplicitSystem &sys = solver->system();
00191 #if PETSC_RELEASE_LESS_THAN(3,5,0)
00192     PetscMatrix<Number> PC(*pc, sys.comm());
00193     PetscMatrix<Number> Jac(*jac, sys.comm());
00194 #else
00195     PetscMatrix<Number> PC(pc, sys.comm());
00196     PetscMatrix<Number> Jac(jac, sys.comm());
00197 #endif
00198     PetscVector<Number>& X_sys = *cast_ptr<PetscVector<Number>*>(sys.solution.get());
00199     PetscMatrix<Number>& Jac_sys = *cast_ptr<PetscMatrix<Number>*>(sys.matrix);
00200     PetscVector<Number> X_global(x, sys.comm());
00201 
00202     // Set the dof maps
00203     PC.attach_dof_map(sys.get_dof_map());
00204     Jac.attach_dof_map(sys.get_dof_map());
00205 
00206     // Use the systems update() to get a good local version of the parallel solution
00207     X_global.swap(X_sys);
00208     Jac.swap(Jac_sys);
00209 
00210     sys.get_dof_map().enforce_constraints_exactly(sys);
00211     sys.update();
00212 
00213     X_global.swap(X_sys);
00214     Jac.swap(Jac_sys);
00215 
00216     if (solver->_zero_out_jacobian)
00217       PC.zero();
00218 
00219     //-----------------------------------------------------------------------------
00220     // if the user has provided both function pointers and objects only the pointer
00221     // will be used, so catch that as an error
00222     if (solver->jacobian && solver->jacobian_object)
00223       libmesh_error_msg("ERROR: cannot specify both a function and object to compute the Jacobian!");
00224 
00225     if (solver->matvec && solver->residual_and_jacobian_object)
00226       libmesh_error_msg("ERROR: cannot specify both a function and object to compute the combined Residual & Jacobian!");
00227 
00228     if (solver->jacobian != NULL)
00229       solver->jacobian(*sys.current_local_solution.get(), PC, sys);
00230 
00231     else if (solver->jacobian_object != NULL)
00232       solver->jacobian_object->jacobian(*sys.current_local_solution.get(), PC, sys);
00233 
00234     else if (solver->matvec != NULL)
00235       solver->matvec(*sys.current_local_solution.get(), NULL, &PC, sys);
00236 
00237     else if (solver->residual_and_jacobian_object != NULL)
00238       solver->residual_and_jacobian_object->residual_and_jacobian (*sys.current_local_solution.get(), NULL, &PC, sys);
00239 
00240     else
00241       libmesh_error_msg("Error! Unable to compute residual and/or Jacobian!");
00242 
00243     PC.close();
00244     Jac.close();
00245     X_global.close();
00246 #if PETSC_RELEASE_LESS_THAN(3,5,0)
00247     *msflag = SAME_NONZERO_PATTERN;
00248 #endif
00249     STOP_LOG("jacobian()", "PetscNonlinearSolver");
00250 
00251     return ierr;
00252   }
00253 
00254 } // end extern "C"
00255 //---------------------------------------------------------------------
00256 
00257 
00258 
00259 //---------------------------------------------------------------------
00260 // PetscNonlinearSolver<> methods
00261 template <typename T>
00262 PetscNonlinearSolver<T>::PetscNonlinearSolver (sys_type& system_in) :
00263   NonlinearSolver<T>(system_in),
00264   _reason(SNES_CONVERGED_ITERATING/*==0*/), // Arbitrary initial value...
00265   _n_linear_iterations(0),
00266   _current_nonlinear_iteration_number(0),
00267   _zero_out_residual(true),
00268   _zero_out_jacobian(true),
00269   _default_monitor(true)
00270 {
00271 }
00272 
00273 
00274 
00275 template <typename T>
00276 PetscNonlinearSolver<T>::~PetscNonlinearSolver ()
00277 {
00278   this->clear ();
00279 }
00280 
00281 
00282 
00283 template <typename T>
00284 void PetscNonlinearSolver<T>::clear ()
00285 {
00286   if (this->initialized())
00287     {
00288       this->_is_initialized = false;
00289 
00290       PetscErrorCode ierr=0;
00291 
00292       ierr = LibMeshSNESDestroy(&_snes);
00293       LIBMESH_CHKERRABORT(ierr);
00294 
00295       // Reset the nonlinear iteration counter.  This information is only relevant
00296       // *during* the solve().  After the solve is completed it should return to
00297       // the default value of 0.
00298       _current_nonlinear_iteration_number = 0;
00299     }
00300 }
00301 
00302 
00303 
00304 template <typename T>
00305 void PetscNonlinearSolver<T>::init (const char* name)
00306 {
00307   // Initialize the data structures if not done so already.
00308   if (!this->initialized())
00309     {
00310       this->_is_initialized = true;
00311 
00312       PetscErrorCode ierr=0;
00313 
00314 #if PETSC_VERSION_LESS_THAN(2,1,2)
00315       // At least until Petsc 2.1.1, the SNESCreate had a different calling syntax.
00316       // The second argument was of type SNESProblemType, and could have a value of
00317       // either SNES_NONLINEAR_EQUATIONS or SNES_UNCONSTRAINED_MINIMIZATION.
00318       ierr = SNESCreate(this->comm().get(), SNES_NONLINEAR_EQUATIONS, &_snes);
00319       LIBMESH_CHKERRABORT(ierr);
00320 
00321 #else
00322 
00323       ierr = SNESCreate(this->comm().get(),&_snes);
00324       LIBMESH_CHKERRABORT(ierr);
00325 
00326 #endif
00327 
00328       if (name)
00329         {
00330           ierr = SNESSetOptionsPrefix(_snes, name);
00331           LIBMESH_CHKERRABORT(ierr);
00332         }
00333 
00334 #if !PETSC_RELEASE_LESS_THAN(3,3,0)
00335       // Attaching a DM to SNES.
00336       DM dm;
00337       ierr = DMCreate(this->comm().get(), &dm);LIBMESH_CHKERRABORT(ierr);
00338       ierr = DMSetType(dm,DMLIBMESH);LIBMESH_CHKERRABORT(ierr);
00339       ierr = DMlibMeshSetSystem(dm,this->system());LIBMESH_CHKERRABORT(ierr);
00340       if (name)
00341         {
00342           ierr = DMSetOptionsPrefix(dm,name);    LIBMESH_CHKERRABORT(ierr);
00343         }
00344       ierr = DMSetFromOptions(dm);               LIBMESH_CHKERRABORT(ierr);
00345       ierr = DMSetUp(dm);                        LIBMESH_CHKERRABORT(ierr);
00346       ierr = SNESSetDM(this->_snes, dm);         LIBMESH_CHKERRABORT(ierr);
00347       // SNES now owns the reference to dm.
00348       ierr = DMDestroy(&dm);                     LIBMESH_CHKERRABORT(ierr);
00349 
00350 #endif
00351 
00352       if (_default_monitor)
00353         {
00354 #if PETSC_VERSION_LESS_THAN(2,3,3)
00355           ierr = SNESSetMonitor (_snes, __libmesh_petsc_snes_monitor,
00356                                  this, PETSC_NULL);
00357 #else
00358           // API name change in PETSc 2.3.3
00359           ierr = SNESMonitorSet (_snes, __libmesh_petsc_snes_monitor,
00360                                  this, PETSC_NULL);
00361 #endif
00362           LIBMESH_CHKERRABORT(ierr);
00363         }
00364 
00365 #if PETSC_VERSION_LESS_THAN(3,1,0)
00366       // Cannot call SNESSetOptions before SNESSetFunction when using
00367       // any matrix free options with PETSc 3.1.0+
00368       ierr = SNESSetFromOptions(_snes);
00369       LIBMESH_CHKERRABORT(ierr);
00370 #endif
00371 
00372       if(this->_preconditioner)
00373         {
00374           KSP ksp;
00375           ierr = SNESGetKSP (_snes, &ksp);
00376           LIBMESH_CHKERRABORT(ierr);
00377           PC pc;
00378           ierr = KSPGetPC(ksp,&pc);
00379           LIBMESH_CHKERRABORT(ierr);
00380 
00381           this->_preconditioner->init();
00382 
00383           PCSetType(pc, PCSHELL);
00384           PCShellSetContext(pc,(void*)this->_preconditioner);
00385 
00386           //Re-Use the shell functions from petsc_linear_solver
00387           PCShellSetSetUp(pc,__libmesh_petsc_preconditioner_setup);
00388           PCShellSetApply(pc,__libmesh_petsc_preconditioner_apply);
00389         }
00390     }
00391 }
00392 
00393 #if !PETSC_VERSION_LESS_THAN(3,3,0)
00394 template <typename T>
00395 void
00396 PetscNonlinearSolver<T>::build_mat_null_space(NonlinearImplicitSystem::ComputeVectorSubspace* computeSubspaceObject,
00397                                               void (*computeSubspace)(std::vector<NumericVector<Number>*>&, sys_type&),
00398                                               MatNullSpace *msp)
00399 {
00400   PetscErrorCode ierr;
00401   std::vector<NumericVector<Number>* > sp;
00402   if (computeSubspaceObject)
00403     (*computeSubspaceObject)(sp, this->system());
00404   else
00405     (*computeSubspace)(sp, this->system());
00406 
00407   *msp = PETSC_NULL;
00408   if (sp.size())
00409     {
00410       Vec *modes;
00411       PetscScalar *dots;
00412       PetscInt nmodes = cast_int<PetscInt>(sp.size());
00413 
00414 #if PETSC_RELEASE_LESS_THAN(3,5,0)
00415       ierr = PetscMalloc2(nmodes,Vec,&modes,nmodes,PetscScalar,&dots);
00416 #else
00417       ierr = PetscMalloc2(nmodes,&modes,nmodes,&dots);
00418 #endif
00419       LIBMESH_CHKERRABORT(ierr);
00420 
00421       for (PetscInt i=0; i<nmodes; ++i)
00422         {
00423           PetscVector<T>* pv = cast_ptr<PetscVector<T>*>(sp[i]);
00424           Vec v = pv->vec();
00425 
00426           ierr = VecDuplicate(v, modes+i);
00427           LIBMESH_CHKERRABORT(ierr);
00428 
00429           ierr = VecCopy(v,modes[i]);
00430           LIBMESH_CHKERRABORT(ierr);
00431         }
00432 
00433       // Normalize.
00434       ierr = VecNormalize(modes[0],PETSC_NULL);
00435       LIBMESH_CHKERRABORT(ierr);
00436 
00437       for (PetscInt i=1; i<nmodes; i++)
00438         {
00439           // Orthonormalize vec[i] against vec[0:i-1]
00440           ierr = VecMDot(modes[i],i,modes,dots);
00441           LIBMESH_CHKERRABORT(ierr);
00442 
00443           for (PetscInt j=0; j<i; j++)
00444             dots[j] *= -1.;
00445 
00446           ierr = VecMAXPY(modes[i],i,dots,modes);
00447           LIBMESH_CHKERRABORT(ierr);
00448 
00449           ierr = VecNormalize(modes[i],PETSC_NULL);
00450           LIBMESH_CHKERRABORT(ierr);
00451         }
00452 
00453       ierr = MatNullSpaceCreate(this->comm().get(), PETSC_FALSE, nmodes, modes, msp);
00454       LIBMESH_CHKERRABORT(ierr);
00455 
00456       for (PetscInt i=0; i<nmodes; ++i)
00457         {
00458           ierr = VecDestroy(modes+i);
00459           LIBMESH_CHKERRABORT(ierr);
00460         }
00461 
00462       ierr = PetscFree2(modes,dots);
00463       LIBMESH_CHKERRABORT(ierr);
00464     }
00465 }
00466 #endif
00467 
00468 template <typename T>
00469 std::pair<unsigned int, Real>
00470 PetscNonlinearSolver<T>::solve (SparseMatrix<T>&  jac_in,  // System Jacobian Matrix
00471                                 NumericVector<T>& x_in,    // Solution vector
00472                                 NumericVector<T>& r_in,    // Residual vector
00473                                 const double,              // Stopping tolerance
00474                                 const unsigned int)
00475 {
00476   START_LOG("solve()", "PetscNonlinearSolver");
00477   this->init ();
00478 
00479   // Make sure the data passed in are really of Petsc types
00480   PetscMatrix<T>* jac = cast_ptr<PetscMatrix<T>*>(&jac_in);
00481   PetscVector<T>* x   = cast_ptr<PetscVector<T>*>(&x_in);
00482   PetscVector<T>* r   = cast_ptr<PetscVector<T>*>(&r_in);
00483 
00484   PetscErrorCode ierr=0;
00485   PetscInt n_iterations =0;
00486   // Should actually be a PetscReal, but I don't know which version of PETSc first introduced PetscReal
00487   Real final_residual_norm=0.;
00488 
00489   ierr = SNESSetFunction (_snes, r->vec(), __libmesh_petsc_snes_residual, this);
00490   LIBMESH_CHKERRABORT(ierr);
00491 
00492   // Only set the jacobian function if we've been provided with something to call.
00493   // This allows a user to set their own jacobian function if they want to
00494   if (this->jacobian || this->jacobian_object || this->residual_and_jacobian_object)
00495     {
00496       ierr = SNESSetJacobian (_snes, jac->mat(), jac->mat(), __libmesh_petsc_snes_jacobian, this);
00497       LIBMESH_CHKERRABORT(ierr);
00498     }
00499 #if !PETSC_VERSION_LESS_THAN(3,3,0)
00500   // Only set the nullspace if we have a way of computing it and the result is non-empty.
00501   if (this->nullspace || this->nullspace_object)
00502     {
00503       MatNullSpace msp;
00504       this->build_mat_null_space(this->nullspace_object, this->nullspace, &msp);
00505       if (msp)
00506         {
00507           ierr = MatSetNullSpace(jac->mat(), msp);
00508           LIBMESH_CHKERRABORT(ierr);
00509 
00510           ierr = MatNullSpaceDestroy(&msp);
00511           LIBMESH_CHKERRABORT(ierr);
00512         }
00513     }
00514 
00515   // Only set the nearnullspace if we have a way of computing it and the result is non-empty.
00516   if (this->nearnullspace || this->nearnullspace_object)
00517     {
00518       MatNullSpace msp = PETSC_NULL;
00519       this->build_mat_null_space(this->nearnullspace_object, this->nearnullspace, &msp);
00520 
00521       if(msp) {
00522         ierr = MatSetNearNullSpace(jac->mat(), msp);
00523         LIBMESH_CHKERRABORT(ierr);
00524 
00525         ierr = MatNullSpaceDestroy(&msp);
00526         LIBMESH_CHKERRABORT(ierr);
00527       }
00528     }
00529 #endif
00530   // Have the Krylov subspace method use our good initial guess rather than 0
00531   KSP ksp;
00532   ierr = SNESGetKSP (_snes, &ksp);
00533   LIBMESH_CHKERRABORT(ierr);
00534 
00535   // Set the tolerances for the iterative solver.  Use the user-supplied
00536   // tolerance for the relative residual & leave the others at default values
00537   ierr = KSPSetTolerances (ksp, this->initial_linear_tolerance, PETSC_DEFAULT,
00538                            PETSC_DEFAULT, this->max_linear_iterations);
00539   LIBMESH_CHKERRABORT(ierr);
00540 
00541   // Set the tolerances for the non-linear solver.
00542   ierr = SNESSetTolerances(_snes, this->absolute_residual_tolerance, this->relative_residual_tolerance,
00543                            this->relative_step_tolerance, this->max_nonlinear_iterations, this->max_function_evaluations);
00544   LIBMESH_CHKERRABORT(ierr);
00545 
00546   //Pull in command-line options
00547   KSPSetFromOptions(ksp);
00548   SNESSetFromOptions(_snes);
00549 
00550   if (this->user_presolve)
00551     this->user_presolve(this->system());
00552 
00553   //Set the preconditioning matrix
00554   if(this->_preconditioner)
00555     {
00556       this->_preconditioner->set_matrix(jac_in);
00557       this->_preconditioner->init();
00558     }
00559 
00560   //    ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
00561   //           LIBMESH_CHKERRABORT(ierr);
00562 
00563   // Older versions (at least up to 2.1.5) of SNESSolve took 3 arguments,
00564   // the last one being a pointer to an int to hold the number of iterations required.
00565 # if PETSC_VERSION_LESS_THAN(2,2,0)
00566 
00567   ierr = SNESSolve (_snes, x->vec(), &n_iterations);
00568   LIBMESH_CHKERRABORT(ierr);
00569 
00570   // 2.2.x style
00571 #elif PETSC_VERSION_LESS_THAN(2,3,0)
00572 
00573   ierr = SNESSolve (_snes, x->vec());
00574   LIBMESH_CHKERRABORT(ierr);
00575 
00576   // 2.3.x & newer style
00577 #else
00578 
00579   ierr = SNESSolve (_snes, PETSC_NULL, x->vec());
00580   LIBMESH_CHKERRABORT(ierr);
00581 
00582   ierr = SNESGetIterationNumber(_snes,&n_iterations);
00583   LIBMESH_CHKERRABORT(ierr);
00584 
00585   ierr = SNESGetLinearSolveIterations(_snes, &_n_linear_iterations);
00586   LIBMESH_CHKERRABORT(ierr);
00587 
00588 #endif
00589 
00590   // SNESGetFunction has been around forever and should work on all
00591   // versions of PETSc.  This is also now the recommended approach
00592   // according to the documentation for the PETSc 3.5.1 release:
00593   // http://www.mcs.anl.gov/petsc/documentation/changes/35.html
00594   Vec f;
00595   ierr = SNESGetFunction(_snes, &f, 0, 0);
00596   LIBMESH_CHKERRABORT(ierr);
00597   ierr = VecNorm(f, NORM_2, &final_residual_norm);
00598   LIBMESH_CHKERRABORT(ierr);
00599 
00600   // Get and store the reason for convergence
00601   SNESGetConvergedReason(_snes, &_reason);
00602 
00603   //Based on Petsc 2.3.3 documentation all diverged reasons are negative
00604   this->converged = (_reason >= 0);
00605 
00606   this->clear();
00607 
00608   STOP_LOG("solve()", "PetscNonlinearSolver");
00609 
00610   // return the # of its. and the final residual norm.
00611   return std::make_pair(n_iterations, final_residual_norm);
00612 }
00613 
00614 
00615 
00616 template <typename T>
00617 void PetscNonlinearSolver<T>::print_converged_reason()
00618 {
00619 
00620   libMesh::out << "Nonlinear solver convergence/divergence reason: "
00621                << SNESConvergedReasons[this->get_converged_reason()] << std::endl;
00622 }
00623 
00624 
00625 
00626 template <typename T>
00627 SNESConvergedReason PetscNonlinearSolver<T>::get_converged_reason()
00628 {
00629   PetscErrorCode ierr=0;
00630 
00631   if (this->initialized())
00632     {
00633       ierr = SNESGetConvergedReason(_snes, &_reason);
00634       LIBMESH_CHKERRABORT(ierr);
00635     }
00636 
00637   return _reason;
00638 }
00639 
00640 template <typename T>
00641 int PetscNonlinearSolver<T>::get_total_linear_iterations()
00642 {
00643   return _n_linear_iterations;
00644 }
00645 
00646 
00647 //------------------------------------------------------------------
00648 // Explicit instantiations
00649 template class PetscNonlinearSolver<Number>;
00650 
00651 } // namespace libMesh
00652 
00653 
00654 
00655 #endif // #ifdef LIBMESH_HAVE_PETSC