AbstractNonlinearAssemblerSolverHybrid.hpp

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 #ifndef ABSTRACTNONLINEARASSEMBLERSOLVERHYBRID_HPP_
00030 #define ABSTRACTNONLINEARASSEMBLERSOLVERHYBRID_HPP_
00031 
00032 #include "BoundaryConditionsContainer.hpp"
00033 #include "AbstractNonlinearEllipticPde.hpp"
00034 #include "AbstractNonlinearSolver.hpp"
00035 #include "PetscTools.hpp"
00036 #include "PetscMatTools.hpp"
00037 #include "PetscVecTools.hpp"
00038 #include "AbstractFeVolumeIntegralAssembler.hpp"
00039 #include "LinearSystem.hpp"
00040 #include "SimplePetscNonlinearSolver.hpp"
00041 #include "NaturalNeumannSurfaceTermAssembler.hpp"
00042 
00043 /*
00044  * Definitions of GLOBAL functions used by PETSc nonlinear solver
00045  * (implementations are at the bottom of this file).
00046  */
00047 
00062 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00063 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeResidual(SNES snes,
00064                                                                       Vec currentGuess,
00065                                                                       Vec residualVector,
00066                                                                       void* pContext);
00067 
00084 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00085 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian(SNES snes,
00086                                                                       Vec currentGuess,
00087                                                                       Mat* pGlobalJacobian,
00088                                                                       Mat* pPreconditioner,
00089                                                                       MatStructure* pMatStructure,
00090                                                                       void* pContext);
00091 
00097 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00098 class AbstractNonlinearAssemblerSolverHybrid : public AbstractFeVolumeIntegralAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,true,true,NONLINEAR>
00099 {
00100 protected:
00101 
00103     AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* mpMesh;
00104 
00106     BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* mpBoundaryConditions;
00107 
00109     AbstractNonlinearSolver* mpSolver;
00110 
00112     bool mWeAllocatedSolverMemory;
00113 
00115     bool mUseAnalyticalJacobian;
00116 
00118     NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* mpNeumannSurfaceTermsAssembler;
00119 
00129     void ApplyDirichletConditions(Vec currentGuess, Vec residual, Mat* pMat);
00130 
00138     void ComputeJacobianNumerically(const Vec currentGuess, Mat* pJacobian);
00139 
00140 public:
00141 
00151     void ComputeResidual(const Vec currentGuess, Vec residualVector);
00152 
00164     void ComputeJacobian(const Vec currentGuess, Mat* pJacobian);
00165 
00173     AbstractNonlinearAssemblerSolverHybrid(AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* pMesh,
00174                                            BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* pBoundaryConditions,
00175                                            unsigned numQuadPoints = 2);
00176 
00180     virtual ~AbstractNonlinearAssemblerSolverHybrid();
00181 
00191     virtual Vec Solve(Vec initialGuess, bool useAnalyticalJacobian=false);
00192 
00200     void SetNonlinearSolver(AbstractNonlinearSolver* pNonlinearSolver);
00201 
00215     bool VerifyJacobian(double tol);
00216 };
00217 
00218 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00219 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::AbstractNonlinearAssemblerSolverHybrid(
00220             AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* pMesh,
00221             BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* pBoundaryConditions,
00222             unsigned numQuadPoints)
00223     :  AbstractFeVolumeIntegralAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,true,true,NONLINEAR>(pMesh,numQuadPoints),
00224        mpMesh(pMesh),
00225        mpBoundaryConditions(pBoundaryConditions)
00226 {
00227     /*
00228      * Note: if this is run with SPACE_DIM != ELEMENT_DIM the class has to be checked:
00229      * there may be lots of places where we should be using SPACE_DIM not ELEMENT_DIM.
00230      */
00231     assert(SPACE_DIM==ELEMENT_DIM);
00232     assert(pMesh!=NULL);
00233     assert(pBoundaryConditions!=NULL);
00234 
00235    // mpAssembler = new SimpleNonlinearEllipticPdeAssembler<ELEMENT_DIM,SPACE_DIM>(mpMesh,pPde,numQuadPoints);
00236     mpSolver = new SimplePetscNonlinearSolver;
00237     mWeAllocatedSolverMemory = true;
00238 
00239     assert(mpMesh->GetNumNodes() == mpMesh->GetDistributedVectorFactory()->GetProblemSize());
00240 
00241     mpNeumannSurfaceTermsAssembler = new NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>(pMesh,pBoundaryConditions,numQuadPoints);
00242 }
00243 
00244 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00245 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::~AbstractNonlinearAssemblerSolverHybrid()
00246 {
00247     if (mWeAllocatedSolverMemory)
00248     {
00249         delete mpSolver;
00250     }
00251     delete mpNeumannSurfaceTermsAssembler;
00252 }
00253 
00254 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00255 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ApplyDirichletConditions(Vec currentGuess, Vec residual, Mat* pJacobian)
00256 {
00257     if (residual)
00258     {
00259         mpBoundaryConditions->ApplyDirichletToNonlinearResidual(currentGuess,
00260                                                                 residual,
00261                                                                 *(mpMesh->GetDistributedVectorFactory()));
00262     }
00263     if (pJacobian)
00264     {
00265         mpBoundaryConditions->ApplyDirichletToNonlinearJacobian(*pJacobian);
00266     }
00267 }
00268 
00269 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00270 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ComputeResidual(const Vec currentGuess, Vec residualVector)
00271 {
00272     // Add the volume integral part of the residual. This will call ComputeVectorTerm, which needs to be
00273     // implemented in the concrete class.
00274     this->SetVectorToAssemble(residualVector,true);
00275     this->SetCurrentSolution(currentGuess);
00276     this->AssembleVector();
00277 
00278     // Add the surface integral contribution - note the negative sign (scale factor).
00279     // This contribution is, for example for a 1 unknown problem
00280     // -integral(g\phi_id dS), where g is the specfied neumann boundary condition
00281     mpNeumannSurfaceTermsAssembler->SetVectorToAssemble(residualVector, false);
00282     mpNeumannSurfaceTermsAssembler->SetScaleFactor(-1.0);
00283     mpNeumannSurfaceTermsAssembler->Assemble();
00284 
00285     PetscVecTools::Finalise(residualVector);
00286 
00287     ApplyDirichletConditions(currentGuess, residualVector, NULL);
00288 }
00289 
00290 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00291 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ComputeJacobian(const Vec currentGuess, Mat* pJacobian)
00292 {
00293     if (mUseAnalyticalJacobian)
00294     {
00295         this->SetMatrixToAssemble(*pJacobian);
00296         this->SetCurrentSolution(currentGuess);
00297         this->AssembleMatrix();
00298 
00299         PetscMatTools::SwitchWriteMode(*pJacobian);
00300 
00301         ApplyDirichletConditions(currentGuess, NULL, pJacobian);
00302 
00303         PetscMatTools::Finalise(*pJacobian);
00304     }
00305     else
00306     {
00307         ComputeJacobianNumerically(currentGuess, pJacobian);
00308     }
00309 }
00310 
00311 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00312 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ComputeJacobianNumerically(const Vec currentGuess, Mat* pJacobian)
00313 {
00314     unsigned num_unknowns = PROBLEM_DIM * this->mpMesh->GetNumNodes();
00315 
00316     // Set up working vectors
00317     Vec residual;
00318     Vec perturbed_residual;
00319     Vec result;
00320     residual = PetscTools::CreateVec(num_unknowns);
00321     result = PetscTools::CreateVec(num_unknowns);
00322     perturbed_residual = PetscTools::CreateVec(num_unknowns);
00323 
00324     // Copy the currentGuess vector; we perturb the copy
00325     Vec current_guess_copy;
00326     PETSCEXCEPT( VecDuplicate(currentGuess, &current_guess_copy) );
00327     PETSCEXCEPT( VecCopy(currentGuess, current_guess_copy) );
00328 
00329     // Compute the current residual
00330     ComputeResidual(currentGuess, residual);
00331 
00332     // Amount to perturb each input element by
00333     double h = 1e-5;
00334 
00335     PetscInt ilo, ihi;
00336     VecGetOwnershipRange(current_guess_copy, &ilo, &ihi);
00337     unsigned lo = ilo;
00338     unsigned hi = ihi;
00339 
00340     // Iterate over entries in the input vector
00341     for (unsigned global_index_outer = 0; global_index_outer < num_unknowns; global_index_outer++)
00342     {
00343         // Only perturb if we own it
00344         PetscVecTools::AddToElement(current_guess_copy, global_index_outer, h);
00345         ComputeResidual(current_guess_copy, perturbed_residual);
00346 
00347         // result = (perturbed_residual - residual) / h
00348         PetscVecTools::WAXPY(result, -1.0, residual, perturbed_residual);
00349         PetscVecTools::Scale(result, 1.0/h);
00350 
00351         double* p_result;
00352         PETSCEXCEPT( VecGetArray(result, &p_result) );
00353         for (unsigned global_index=lo; global_index < hi; global_index++)
00354         {
00355             unsigned local_index = global_index - lo;
00356             PetscMatTools::SetElement(*pJacobian, global_index, global_index_outer, p_result[local_index]);
00357         }
00358         PETSCEXCEPT( VecRestoreArray(result, &p_result) );
00359 
00360         PetscVecTools::AddToElement(current_guess_copy, global_index_outer, -h);
00361     }
00362 
00363     VecDestroy(residual);
00364     VecDestroy(perturbed_residual);
00365     VecDestroy(result);
00366     VecDestroy(current_guess_copy);
00367 
00368     PetscMatTools::Finalise(*pJacobian);
00369 }
00370 
00371 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00372 Vec AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::Solve(Vec initialGuess,
00373                                                                                        bool useAnalyticalJacobian)
00374 {
00375     assert(initialGuess != NULL);
00376     mUseAnalyticalJacobian = useAnalyticalJacobian;
00377 
00378     PetscInt size_of_init_guess;
00379     VecGetSize(initialGuess, &size_of_init_guess);
00380     PetscInt problem_size = PROBLEM_DIM * this->mpMesh->GetNumNodes();
00381     if (size_of_init_guess != problem_size)
00382     {
00383         EXCEPTION("Size of initial guess vector, " << size_of_init_guess
00384                       << ", does not match size of problem, " << problem_size);
00385     }
00386 
00387     // Run the solver, telling it which global functions to call in order to assemble the residual or jacobian
00388     return mpSolver->Solve(&AbstractNonlinearAssemblerSolverHybrid_ComputeResidual<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>,
00389                            &AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>,
00390                            initialGuess,
00391                            PROBLEM_DIM * this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00392                            this);
00393 }
00394 
00395 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00396 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::SetNonlinearSolver(AbstractNonlinearSolver* pNonlinearSolver)
00397 {
00398     if (mWeAllocatedSolverMemory)
00399     {
00400         delete mpSolver;
00401     }
00402     mpSolver = pNonlinearSolver;
00403     mWeAllocatedSolverMemory = false;
00404 }
00405 
00406 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00407 bool AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::VerifyJacobian(double tol)
00408 {
00409     unsigned size = PROBLEM_DIM * this->mpMesh->GetNumNodes();
00410 
00411     Vec initial_guess = PetscTools::CreateAndSetVec(size, 0.0);
00412 
00413     Mat analytic_jacobian;
00414     Mat numerical_jacobian;
00415 
00416     PetscTools::SetupMat(analytic_jacobian, size, size, size);
00417     PetscTools::SetupMat(numerical_jacobian, size, size, size);
00418 
00419     mUseAnalyticalJacobian = true;
00420     ComputeJacobian(initial_guess, &analytic_jacobian);
00421 
00422     mUseAnalyticalJacobian = false;
00423     ComputeJacobian(initial_guess, &numerical_jacobian);
00424 
00425     bool ok = PetscMatTools::CheckEquality(numerical_jacobian, analytic_jacobian, tol);
00426     MatDestroy(numerical_jacobian);
00427     MatDestroy(analytic_jacobian);
00428     VecDestroy(initial_guess);
00429 
00430     return ok;
00431 }
00432 
00433 // Implementations of GLOBAL functions called by PETSc nonlinear solver
00434 
00435 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00436 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeResidual(SNES snes,
00437                                                                       Vec currentGuess,
00438                                                                       Vec residualVector,
00439                                                                       void* pContext)
00440 {
00441     // Extract the solver from the void*
00442     AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* p_solver =
00443         (AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>*) pContext;
00444 
00445     p_solver->ComputeResidual(currentGuess, residualVector);
00446 
00447     return 0;
00448 }
00449 
00450 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00451 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian(SNES snes,
00452                                                                       Vec currentGuess,
00453                                                                       Mat* pGlobalJacobian,
00454                                                                       Mat* pPreconditioner,
00455                                                                       MatStructure* pMatStructure,
00456                                                                       void* pContext)
00457 {
00458     // Extract the solver from the void*
00459     AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* p_solver =
00460         (AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>*) pContext;
00461 
00462     p_solver->ComputeJacobian(currentGuess, pGlobalJacobian);
00463 
00464     return 0;
00465 }
00466 
00467 #endif /*ABSTRACTNONLINEARASSEMBLERSOLVERHYBRID_HPP_*/
Generated on Thu Dec 22 13:00:17 2011 for Chaste by  doxygen 1.6.3