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 
00033 #include "SimplePetscNonlinearSolver.hpp"
00034 #include "Exception.hpp"
00035 #include "petscsnes.h"
00036 #include "PetscTools.hpp"
00037 #include <sstream>
00038 
00039 Vec SimplePetscNonlinearSolver::Solve(PetscErrorCode (*pComputeResidual)(SNES,Vec,Vec,void*),
00040                                       PetscErrorCode (*pComputeJacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),
00041                                       Vec initialGuess,
00042                                       unsigned fill,
00043                                       void* pContext)
00044 {
00045     SNES snes;
00046 
00047     // Create the residual vector by copying the structure of the initial guess
00048     Vec residual;
00049     VecDuplicate(initialGuess, &residual);
00050 
00051     Mat jacobian; // Jacobian Matrix
00052 
00053     PetscInt N; // number of elements
00054 
00055     // Get the size of the Jacobian from the residual
00056     VecGetSize(initialGuess,&N);
00057 
00058     PetscTools::SetupMat(jacobian, N, N, fill);
00059 
00060     SNESCreate(PETSC_COMM_WORLD, &snes);
00061     SNESSetFunction(snes, residual, pComputeResidual, pContext);
00062     SNESSetJacobian(snes, jacobian, jacobian, pComputeJacobian, pContext);
00063     SNESSetType(snes,SNESLS);
00064     SNESSetTolerances(snes,1.0e-5,1.0e-5,1.0e-5,PETSC_DEFAULT,PETSC_DEFAULT);
00065 
00066     // x is the iteration vector SNES uses when solving, set equal to initialGuess to start with
00067     Vec x;
00068     VecDuplicate(initialGuess, &x);
00069     VecCopy(initialGuess, x);
00070 
00071 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00072     SNESSolve(snes, x);
00073 #else
00074     SNESSolve(snes, PETSC_NULL, x);
00075 #endif
00076 
00077     VecDestroy(residual);
00078     MatDestroy(jacobian); // Free Jacobian
00079 
00080     SNESConvergedReason reason;
00081     SNESGetConvergedReason(snes,&reason);
00082 #define COVERAGE_IGNORE
00083     if (reason < 0)
00084     {
00085         std::stringstream reason_stream;
00086         reason_stream << reason;
00087         VecDestroy(x); // Since caller can't free the memory in this case
00088         SNESDestroy(snes);
00089         EXCEPTION("Nonlinear Solver did not converge. PETSc reason code:"
00090                   +reason_stream.str()+" .");
00091     }
00092 #undef COVERAGE_IGNORE
00093     SNESDestroy(snes);
00094 
00095     return x;
00096 }
Generated on Thu Dec 22 13:00:07 2011 for Chaste by  doxygen 1.6.3