00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00040 #include <sstream>
00041 #include "petscsnes.h"
00042 #include "SimplePetscNonlinearSolver.hpp"
00043 #include "Exception.hpp"
00044 #include "PetscTools.hpp"
00045
00046 Vec SimplePetscNonlinearSolver::Solve(PetscErrorCode (*pComputeResidual)(SNES,Vec,Vec,void*),
00047 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
00048 PetscErrorCode (*pComputeJacobian)(SNES,Vec,Mat,Mat,void*),
00049 #else
00050 PetscErrorCode (*pComputeJacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),
00051 #endif
00052 Vec initialGuess,
00053 unsigned fill,
00054 void* pContext)
00055 {
00056 SNES snes;
00057
00058
00059 Vec residual;
00060 VecDuplicate(initialGuess, &residual);
00061
00062 Mat jacobian;
00063 PetscInt N;
00064
00065
00066 VecGetSize(initialGuess, &N);
00067
00068
00069 PetscTools::SetupMat(jacobian, N, N, fill, PETSC_DECIDE, PETSC_DECIDE, true, false );
00070
00071 SNESCreate(PETSC_COMM_WORLD, &snes);
00072
00073 SNESSetFunction(snes, residual, pComputeResidual, pContext);
00074 SNESSetJacobian(snes, jacobian, jacobian, pComputeJacobian, pContext);
00075
00076 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 4) //PETSc 3.4 or later
00077 SNESSetType(snes, SNESNEWTONLS);
00078 #else
00079 SNESSetType(snes, SNESLS);
00080 #endif
00081
00082 SNESSetTolerances(snes,1.0e-5,1.0e-5,1.0e-5,PETSC_DEFAULT,PETSC_DEFAULT);
00083 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR == 3) //PETSc 3.3
00084 SNESLineSearch linesearch;
00085 SNESGetSNESLineSearch(snes, &linesearch);
00086 SNESLineSearchSetType(linesearch, "bt");
00087 #elif (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 4) //PETSc 3.4 or later
00088 SNESLineSearch linesearch;
00089 SNESGetLineSearch(snes, &linesearch);
00090 SNESLineSearchSetType(linesearch, "bt");
00091 #endif
00092
00093
00094 Vec x;
00095
00096 VecDuplicate(initialGuess, &x);
00097 VecCopy(initialGuess, x);
00098
00099 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 5)
00100
00101
00102
00103
00104 KSP ksp;
00105 SNESGetKSP(snes,&ksp);
00106 PC pc;
00107 KSPGetPC(ksp,&pc);
00108 PCSetType(pc,PCNONE);
00109 #endif
00110
00111 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00112 SNESSolve(snes, x);
00113 #else
00114 SNESSolve(snes, PETSC_NULL, x);
00115 #endif
00116 PetscTools::Destroy(residual);
00117 PetscTools::Destroy(jacobian);
00118
00119 SNESConvergedReason reason;
00120 SNESGetConvergedReason(snes,&reason);
00121 #define COVERAGE_IGNORE
00122 if (reason < 0)
00123 {
00124 std::stringstream reason_stream;
00125 reason_stream << reason;
00126 PetscTools::Destroy(x);
00127 SNESDestroy(PETSC_DESTROY_PARAM(snes));
00128 EXCEPTION("Nonlinear Solver did not converge. PETSc reason code:"
00129 +reason_stream.str()+" .");
00130 }
00131 #undef COVERAGE_IGNORE
00132 SNESDestroy(PETSC_DESTROY_PARAM(snes));
00133
00134 return x;
00135 }