Chaste Release::3.1
|
00001 /* 00002 00003 Copyright (c) 2005-2012, University of Oxford. 00004 All rights reserved. 00005 00006 University of Oxford means the Chancellor, Masters and Scholars of the 00007 University of Oxford, having an administrative office at Wellington 00008 Square, Oxford OX1 2JD, UK. 00009 00010 This file is part of Chaste. 00011 00012 Redistribution and use in source and binary forms, with or without 00013 modification, are permitted provided that the following conditions are met: 00014 * Redistributions of source code must retain the above copyright notice, 00015 this list of conditions and the following disclaimer. 00016 * Redistributions in binary form must reproduce the above copyright notice, 00017 this list of conditions and the following disclaimer in the documentation 00018 and/or other materials provided with the distribution. 00019 * Neither the name of the University of Oxford nor the names of its 00020 contributors may be used to endorse or promote products derived from this 00021 software without specific prior written permission. 00022 00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" 00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE 00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE 00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR 00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE 00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) 00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT 00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT 00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00033 00034 */ 00035 00036 #ifndef ABSTRACTNONLINEARASSEMBLERSOLVERHYBRID_HPP_ 00037 #define ABSTRACTNONLINEARASSEMBLERSOLVERHYBRID_HPP_ 00038 00039 #include "BoundaryConditionsContainer.hpp" 00040 #include "AbstractNonlinearEllipticPde.hpp" 00041 #include "AbstractNonlinearSolver.hpp" 00042 #include "PetscTools.hpp" 00043 #include "PetscMatTools.hpp" 00044 #include "PetscVecTools.hpp" 00045 #include "AbstractFeVolumeIntegralAssembler.hpp" 00046 #include "LinearSystem.hpp" 00047 #include "SimplePetscNonlinearSolver.hpp" 00048 #include "NaturalNeumannSurfaceTermAssembler.hpp" 00049 00050 /* 00051 * Definitions of GLOBAL functions used by PETSc nonlinear solver 00052 * (implementations are at the bottom of this file). 00053 */ 00054 00069 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM> 00070 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeResidual(SNES snes, 00071 Vec currentGuess, 00072 Vec residualVector, 00073 void* pContext); 00074 00091 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM> 00092 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian(SNES snes, 00093 Vec currentGuess, 00094 Mat* pGlobalJacobian, 00095 Mat* pPreconditioner, 00096 MatStructure* pMatStructure, 00097 void* pContext); 00098 00104 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM> 00105 class AbstractNonlinearAssemblerSolverHybrid : public AbstractFeVolumeIntegralAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,true,true,NONLINEAR> 00106 { 00107 protected: 00108 00110 AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* mpMesh; 00111 00113 BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* mpBoundaryConditions; 00114 00116 AbstractNonlinearSolver* mpSolver; 00117 00119 bool mWeAllocatedSolverMemory; 00120 00122 bool mUseAnalyticalJacobian; 00123 00125 NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* mpNeumannSurfaceTermsAssembler; 00126 00136 void ApplyDirichletConditions(Vec currentGuess, Vec residual, Mat* pMat); 00137 00145 void ComputeJacobianNumerically(const Vec currentGuess, Mat* pJacobian); 00146 00147 public: 00148 00158 void ComputeResidual(const Vec currentGuess, Vec residualVector); 00159 00171 void ComputeJacobian(const Vec currentGuess, Mat* pJacobian); 00172 00180 AbstractNonlinearAssemblerSolverHybrid(AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* pMesh, 00181 BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* pBoundaryConditions, 00182 unsigned numQuadPoints = 2); 00183 00187 virtual ~AbstractNonlinearAssemblerSolverHybrid(); 00188 00198 virtual Vec Solve(Vec initialGuess, bool useAnalyticalJacobian=false); 00199 00207 void SetNonlinearSolver(AbstractNonlinearSolver* pNonlinearSolver); 00208 00222 bool VerifyJacobian(double tol); 00223 }; 00224 00225 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM> 00226 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::AbstractNonlinearAssemblerSolverHybrid( 00227 AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* pMesh, 00228 BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* pBoundaryConditions, 00229 unsigned numQuadPoints) 00230 : AbstractFeVolumeIntegralAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,true,true,NONLINEAR>(pMesh,numQuadPoints), 00231 mpMesh(pMesh), 00232 mpBoundaryConditions(pBoundaryConditions) 00233 { 00234 /* 00235 * Note: if this is run with SPACE_DIM != ELEMENT_DIM the class has to be checked: 00236 * there may be lots of places where we should be using SPACE_DIM not ELEMENT_DIM. 00237 */ 00238 assert(SPACE_DIM==ELEMENT_DIM); 00239 assert(pMesh!=NULL); 00240 assert(pBoundaryConditions!=NULL); 00241 00242 // mpAssembler = new SimpleNonlinearEllipticPdeAssembler<ELEMENT_DIM,SPACE_DIM>(mpMesh,pPde,numQuadPoints); 00243 mpSolver = new SimplePetscNonlinearSolver; 00244 mWeAllocatedSolverMemory = true; 00245 00246 assert(mpMesh->GetNumNodes() == mpMesh->GetDistributedVectorFactory()->GetProblemSize()); 00247 00248 mpNeumannSurfaceTermsAssembler = new NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>(pMesh,pBoundaryConditions,numQuadPoints); 00249 } 00250 00251 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM> 00252 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::~AbstractNonlinearAssemblerSolverHybrid() 00253 { 00254 if (mWeAllocatedSolverMemory) 00255 { 00256 delete mpSolver; 00257 } 00258 delete mpNeumannSurfaceTermsAssembler; 00259 } 00260 00261 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM> 00262 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ApplyDirichletConditions(Vec currentGuess, Vec residual, Mat* pJacobian) 00263 { 00264 if (residual) 00265 { 00266 mpBoundaryConditions->ApplyDirichletToNonlinearResidual(currentGuess, 00267 residual, 00268 *(mpMesh->GetDistributedVectorFactory())); 00269 } 00270 if (pJacobian) 00271 { 00272 mpBoundaryConditions->ApplyDirichletToNonlinearJacobian(*pJacobian); 00273 } 00274 } 00275 00276 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM> 00277 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ComputeResidual(const Vec currentGuess, Vec residualVector) 00278 { 00279 // Add the volume integral part of the residual. This will call ComputeVectorTerm, which needs to be 00280 // implemented in the concrete class. 00281 this->SetVectorToAssemble(residualVector,true); 00282 this->SetCurrentSolution(currentGuess); 00283 this->AssembleVector(); 00284 00285 // Add the surface integral contribution - note the negative sign (scale factor). 00286 // This contribution is, for example for a 1 unknown problem 00287 // -integral(g\phi_id dS), where g is the specfied neumann boundary condition 00288 mpNeumannSurfaceTermsAssembler->SetVectorToAssemble(residualVector, false); 00289 mpNeumannSurfaceTermsAssembler->SetScaleFactor(-1.0); 00290 mpNeumannSurfaceTermsAssembler->Assemble(); 00291 00292 PetscVecTools::Finalise(residualVector); 00293 00294 ApplyDirichletConditions(currentGuess, residualVector, NULL); 00295 } 00296 00297 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM> 00298 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ComputeJacobian(const Vec currentGuess, Mat* pJacobian) 00299 { 00300 if (mUseAnalyticalJacobian) 00301 { 00302 this->SetMatrixToAssemble(*pJacobian); 00303 this->SetCurrentSolution(currentGuess); 00304 this->AssembleMatrix(); 00305 00306 PetscMatTools::SwitchWriteMode(*pJacobian); 00307 00308 ApplyDirichletConditions(currentGuess, NULL, pJacobian); 00309 00310 PetscMatTools::Finalise(*pJacobian); 00311 } 00312 else 00313 { 00314 ComputeJacobianNumerically(currentGuess, pJacobian); 00315 } 00316 } 00317 00318 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM> 00319 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ComputeJacobianNumerically(const Vec currentGuess, Mat* pJacobian) 00320 { 00321 unsigned num_unknowns = PROBLEM_DIM * this->mpMesh->GetNumNodes(); 00322 00323 // Set up working vectors 00324 Vec residual; 00325 Vec perturbed_residual; 00326 Vec result; 00327 residual = PetscTools::CreateVec(num_unknowns); 00328 result = PetscTools::CreateVec(num_unknowns); 00329 perturbed_residual = PetscTools::CreateVec(num_unknowns); 00330 00331 // Copy the currentGuess vector; we perturb the copy 00332 Vec current_guess_copy; 00333 PETSCEXCEPT( VecDuplicate(currentGuess, ¤t_guess_copy) ); 00334 PETSCEXCEPT( VecCopy(currentGuess, current_guess_copy) ); 00335 00336 // Compute the current residual 00337 ComputeResidual(currentGuess, residual); 00338 00339 // Amount to perturb each input element by 00340 double h = 1e-5; 00341 00342 PetscInt ilo, ihi; 00343 VecGetOwnershipRange(current_guess_copy, &ilo, &ihi); 00344 unsigned lo = ilo; 00345 unsigned hi = ihi; 00346 00347 // Iterate over entries in the input vector 00348 for (unsigned global_index_outer = 0; global_index_outer < num_unknowns; global_index_outer++) 00349 { 00350 // Only perturb if we own it 00351 PetscVecTools::AddToElement(current_guess_copy, global_index_outer, h); 00352 ComputeResidual(current_guess_copy, perturbed_residual); 00353 00354 // result = (perturbed_residual - residual) / h 00355 PetscVecTools::WAXPY(result, -1.0, residual, perturbed_residual); 00356 PetscVecTools::Scale(result, 1.0/h); 00357 00358 double* p_result; 00359 PETSCEXCEPT( VecGetArray(result, &p_result) ); 00360 for (unsigned global_index=lo; global_index < hi; global_index++) 00361 { 00362 unsigned local_index = global_index - lo; 00363 PetscMatTools::SetElement(*pJacobian, global_index, global_index_outer, p_result[local_index]); 00364 } 00365 PETSCEXCEPT( VecRestoreArray(result, &p_result) ); 00366 00367 PetscVecTools::AddToElement(current_guess_copy, global_index_outer, -h); 00368 } 00369 00370 PetscTools::Destroy(residual); 00371 PetscTools::Destroy(perturbed_residual); 00372 PetscTools::Destroy(result); 00373 PetscTools::Destroy(current_guess_copy); 00374 00375 PetscMatTools::Finalise(*pJacobian); 00376 } 00377 00378 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM> 00379 Vec AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::Solve(Vec initialGuess, 00380 bool useAnalyticalJacobian) 00381 { 00382 assert(initialGuess != NULL); 00383 mUseAnalyticalJacobian = useAnalyticalJacobian; 00384 00385 PetscInt size_of_init_guess; 00386 VecGetSize(initialGuess, &size_of_init_guess); 00387 PetscInt problem_size = PROBLEM_DIM * this->mpMesh->GetNumNodes(); 00388 if (size_of_init_guess != problem_size) 00389 { 00390 EXCEPTION("Size of initial guess vector, " << size_of_init_guess 00391 << ", does not match size of problem, " << problem_size); 00392 } 00393 00394 // Run the solver, telling it which global functions to call in order to assemble the residual or jacobian 00395 return mpSolver->Solve(&AbstractNonlinearAssemblerSolverHybrid_ComputeResidual<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>, 00396 &AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>, 00397 initialGuess, 00398 PROBLEM_DIM * this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(), 00399 this); 00400 } 00401 00402 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM> 00403 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::SetNonlinearSolver(AbstractNonlinearSolver* pNonlinearSolver) 00404 { 00405 if (mWeAllocatedSolverMemory) 00406 { 00407 delete mpSolver; 00408 } 00409 mpSolver = pNonlinearSolver; 00410 mWeAllocatedSolverMemory = false; 00411 } 00412 00413 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM> 00414 bool AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::VerifyJacobian(double tol) 00415 { 00416 unsigned size = PROBLEM_DIM * this->mpMesh->GetNumNodes(); 00417 00418 Vec initial_guess = PetscTools::CreateAndSetVec(size, 0.0); 00419 00420 Mat analytic_jacobian; 00421 Mat numerical_jacobian; 00422 00423 PetscTools::SetupMat(analytic_jacobian, size, size, size); 00424 PetscTools::SetupMat(numerical_jacobian, size, size, size); 00425 00426 mUseAnalyticalJacobian = true; 00427 ComputeJacobian(initial_guess, &analytic_jacobian); 00428 00429 mUseAnalyticalJacobian = false; 00430 ComputeJacobian(initial_guess, &numerical_jacobian); 00431 00432 bool ok = PetscMatTools::CheckEquality(numerical_jacobian, analytic_jacobian, tol); 00433 PetscTools::Destroy(numerical_jacobian); 00434 PetscTools::Destroy(analytic_jacobian); 00435 PetscTools::Destroy(initial_guess); 00436 00437 return ok; 00438 } 00439 00440 // Implementations of GLOBAL functions called by PETSc nonlinear solver 00441 00442 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM> 00443 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeResidual(SNES snes, 00444 Vec currentGuess, 00445 Vec residualVector, 00446 void* pContext) 00447 { 00448 // Extract the solver from the void* 00449 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* p_solver = 00450 (AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>*) pContext; 00451 00452 p_solver->ComputeResidual(currentGuess, residualVector); 00453 00454 return 0; 00455 } 00456 00457 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM> 00458 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian(SNES snes, 00459 Vec currentGuess, 00460 Mat* pGlobalJacobian, 00461 Mat* pPreconditioner, 00462 MatStructure* pMatStructure, 00463 void* pContext) 00464 { 00465 // Extract the solver from the void* 00466 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* p_solver = 00467 (AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>*) pContext; 00468 00469 p_solver->ComputeJacobian(currentGuess, pGlobalJacobian); 00470 00471 return 0; 00472 } 00473 00474 #endif /*ABSTRACTNONLINEARASSEMBLERSOLVERHYBRID_HPP_*/