SimplePetscNonlinearSolver.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2011
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
00026 
00027 */
00028 
00029 
00036 #include "SimplePetscNonlinearSolver.hpp"
00037 #include "Exception.hpp"
00038 #include "petscsnes.h"
00039 #include "PetscTools.hpp"
00040 #include <sstream>
00041 
00042 
00043 
00044 Vec SimplePetscNonlinearSolver::Solve(PetscErrorCode (*pComputeResidual)(SNES,Vec,Vec,void*),
00045                                       PetscErrorCode (*pComputeJacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),
00046                                       Vec initialGuess,
00047                                       unsigned fill,
00048                                       void* pContext)
00049 {
00050     SNES snes;
00051 
00052     // create the residual vector by copying the structure of the initial guess
00053     Vec residual;
00054     VecDuplicate(initialGuess, &residual);
00055 
00056     Mat jacobian; //Jacobian Matrix
00057 
00058     PetscInt N; //number of elements
00059     //get the size of the jacobian from the residual
00060     VecGetSize(initialGuess,&N);
00061 
00062     PetscTools::SetupMat(jacobian, N, N, fill);
00063 
00064     SNESCreate(PETSC_COMM_WORLD, &snes);
00065     SNESSetFunction(snes, residual, pComputeResidual, pContext);
00066     SNESSetJacobian(snes, jacobian, jacobian, pComputeJacobian, pContext);
00067     SNESSetType(snes,SNESLS);
00068     SNESSetTolerances(snes,1.0e-5,1.0e-5,1.0e-5,PETSC_DEFAULT,PETSC_DEFAULT);
00069 
00070     // x is the iteration vector SNES uses when solving, set equal to initialGuess to start with
00071     Vec x;
00072     VecDuplicate(initialGuess, &x);
00073     VecCopy(initialGuess, x);
00074 
00075 
00076 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00077     SNESSolve(snes, x);
00078 #else
00079     SNESSolve(snes, PETSC_NULL, x);
00080 #endif
00081 
00082     VecDestroy(residual);
00083     MatDestroy(jacobian); // Free Jacobian
00084 
00085     SNESConvergedReason reason;
00086     SNESGetConvergedReason(snes,&reason);
00087 #define COVERAGE_IGNORE
00088     if (reason<0)
00089     {
00090         std::stringstream reason_stream;
00091         reason_stream << reason;
00092         VecDestroy(x); // Since caller can't free the memory in this case
00093         SNESDestroy(snes);
00094         EXCEPTION("Nonlinear Solver did not converge. Petsc reason code:"
00095                   +reason_stream.str()+" .");
00096     }
00097 #undef COVERAGE_IGNORE
00098     SNESDestroy(snes);
00099 
00100     return x;
00101 }

Generated on Tue May 31 14:31:42 2011 for Chaste by  doxygen 1.5.5