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
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
00052
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 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
00075
00093 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00094
00095 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian(SNES snes,
00096 Vec currentGuess,
00097 Mat globalJacobian,
00098 Mat preconditioner,
00099 void* pContext);
00100 #else
00101
00117 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00118 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian(SNES snes,
00119 Vec currentGuess,
00120 Mat* pGlobalJacobian,
00121 Mat* pPreconditioner,
00122 MatStructure* pMatStructure,
00123 void* pContext);
00124 #endif
00125
00131 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00132 class AbstractNonlinearAssemblerSolverHybrid : public AbstractFeVolumeIntegralAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,true,true,NONLINEAR>
00133 {
00134 protected:
00135
00137 AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* mpMesh;
00138
00140 BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* mpBoundaryConditions;
00141
00143 AbstractNonlinearSolver* mpSolver;
00144
00146 bool mWeAllocatedSolverMemory;
00147
00149 bool mUseAnalyticalJacobian;
00150
00152 NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* mpNeumannSurfaceTermsAssembler;
00153
00163 void ApplyDirichletConditions(Vec currentGuess, Vec residual, Mat* pMat);
00164
00172 void ComputeJacobianNumerically(const Vec currentGuess, Mat* pJacobian);
00173
00174 public:
00175
00185 void ComputeResidual(const Vec currentGuess, Vec residualVector);
00186
00198 void ComputeJacobian(const Vec currentGuess, Mat* pJacobian);
00199
00206 AbstractNonlinearAssemblerSolverHybrid(AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* pMesh,
00207 BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* pBoundaryConditions);
00208
00212 virtual ~AbstractNonlinearAssemblerSolverHybrid();
00213
00223 virtual Vec Solve(Vec initialGuess, bool useAnalyticalJacobian=false);
00224
00232 void SetNonlinearSolver(AbstractNonlinearSolver* pNonlinearSolver);
00233
00247 bool VerifyJacobian(double tol);
00248 };
00249
00250 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00251 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::AbstractNonlinearAssemblerSolverHybrid(
00252 AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* pMesh,
00253 BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* pBoundaryConditions)
00254 : AbstractFeVolumeIntegralAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,true,true,NONLINEAR>(pMesh),
00255 mpMesh(pMesh),
00256 mpBoundaryConditions(pBoundaryConditions)
00257 {
00258
00259
00260
00261
00262 assert(SPACE_DIM==ELEMENT_DIM);
00263 assert(pMesh!=NULL);
00264 assert(pBoundaryConditions!=NULL);
00265
00266 mpSolver = new SimplePetscNonlinearSolver;
00267 mWeAllocatedSolverMemory = true;
00268
00269 assert(mpMesh->GetNumNodes() == mpMesh->GetDistributedVectorFactory()->GetProblemSize());
00270
00271 mpNeumannSurfaceTermsAssembler = new NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>(pMesh,pBoundaryConditions);
00272 }
00273
00274 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00275 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::~AbstractNonlinearAssemblerSolverHybrid()
00276 {
00277 if (mWeAllocatedSolverMemory)
00278 {
00279 delete mpSolver;
00280 }
00281 delete mpNeumannSurfaceTermsAssembler;
00282 }
00283
00284 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00285 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ApplyDirichletConditions(Vec currentGuess, Vec residual, Mat* pJacobian)
00286 {
00287 if (residual)
00288 {
00289 mpBoundaryConditions->ApplyDirichletToNonlinearResidual(currentGuess,
00290 residual,
00291 *(mpMesh->GetDistributedVectorFactory()));
00292 }
00293 if (pJacobian)
00294 {
00295 mpBoundaryConditions->ApplyDirichletToNonlinearJacobian(*pJacobian);
00296 }
00297 }
00298
00299 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00300 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ComputeResidual(const Vec currentGuess, Vec residualVector)
00301 {
00302
00303
00304 this->SetVectorToAssemble(residualVector,true);
00305 this->SetCurrentSolution(currentGuess);
00306 this->AssembleVector();
00307
00308
00309
00310
00311 mpNeumannSurfaceTermsAssembler->SetVectorToAssemble(residualVector, false);
00312 mpNeumannSurfaceTermsAssembler->SetScaleFactor(-1.0);
00313 mpNeumannSurfaceTermsAssembler->Assemble();
00314
00315 PetscVecTools::Finalise(residualVector);
00316
00317 ApplyDirichletConditions(currentGuess, residualVector, NULL);
00318 }
00319
00320 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00321 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ComputeJacobian(const Vec currentGuess, Mat* pJacobian)
00322 {
00323 if (mUseAnalyticalJacobian)
00324 {
00325 this->SetMatrixToAssemble(*pJacobian);
00326 this->SetCurrentSolution(currentGuess);
00327 this->AssembleMatrix();
00328
00329 PetscMatTools::SwitchWriteMode(*pJacobian);
00330
00331 ApplyDirichletConditions(currentGuess, NULL, pJacobian);
00332
00333 PetscMatTools::Finalise(*pJacobian);
00334 }
00335 else
00336 {
00337 ComputeJacobianNumerically(currentGuess, pJacobian);
00338 }
00339 }
00340
00341 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00342 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ComputeJacobianNumerically(const Vec currentGuess, Mat* pJacobian)
00343 {
00344 unsigned num_unknowns = PROBLEM_DIM * this->mpMesh->GetNumNodes();
00345
00346
00347 Vec residual;
00348 Vec perturbed_residual;
00349 Vec result;
00350 residual = PetscTools::CreateVec(num_unknowns);
00351 result = PetscTools::CreateVec(num_unknowns);
00352 perturbed_residual = PetscTools::CreateVec(num_unknowns);
00353
00354
00355 Vec current_guess_copy;
00356 PETSCEXCEPT( VecDuplicate(currentGuess, ¤t_guess_copy) );
00357 PETSCEXCEPT( VecCopy(currentGuess, current_guess_copy) );
00358
00359
00360 ComputeResidual(currentGuess, residual);
00361
00362
00363 double h = 1e-5;
00364 double near_hsquared = 1e-9;
00365
00366 PetscInt ilo, ihi;
00367 VecGetOwnershipRange(current_guess_copy, &ilo, &ihi);
00368 unsigned lo = ilo;
00369 unsigned hi = ihi;
00370
00371
00372 for (unsigned global_index_outer = 0; global_index_outer < num_unknowns; global_index_outer++)
00373 {
00374
00375 PetscVecTools::AddToElement(current_guess_copy, global_index_outer, h);
00376 ComputeResidual(current_guess_copy, perturbed_residual);
00377
00378
00379 PetscVecTools::WAXPY(result, -1.0, residual, perturbed_residual);
00380 PetscVecTools::Scale(result, 1.0/h);
00381
00382 double* p_result;
00385 PETSCEXCEPT( VecGetArray(result, &p_result) );
00386 for (unsigned global_index=lo; global_index < hi; global_index++)
00387 {
00388 double result_entry = p_result[ global_index - lo];
00389 if (!CompareDoubles::IsNearZero(result_entry, near_hsquared))
00390 {
00391 PetscMatTools::SetElement(*pJacobian, global_index, global_index_outer, result_entry);
00392 }
00393 }
00394 PETSCEXCEPT( VecRestoreArray(result, &p_result) );
00395
00396 PetscVecTools::AddToElement(current_guess_copy, global_index_outer, -h);
00397 }
00398
00399 PetscTools::Destroy(residual);
00400 PetscTools::Destroy(perturbed_residual);
00401 PetscTools::Destroy(result);
00402 PetscTools::Destroy(current_guess_copy);
00403
00404 PetscMatTools::Finalise(*pJacobian);
00405 }
00406
00407 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00408 Vec AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::Solve(Vec initialGuess,
00409 bool useAnalyticalJacobian)
00410 {
00411 assert(initialGuess != NULL);
00412 mUseAnalyticalJacobian = useAnalyticalJacobian;
00413
00414 PetscInt size_of_init_guess;
00415 VecGetSize(initialGuess, &size_of_init_guess);
00416 PetscInt problem_size = PROBLEM_DIM * this->mpMesh->GetNumNodes();
00417 if (size_of_init_guess != problem_size)
00418 {
00419 EXCEPTION("Size of initial guess vector, " << size_of_init_guess
00420 << ", does not match size of problem, " << problem_size);
00421 }
00422
00423
00424 return mpSolver->Solve(&AbstractNonlinearAssemblerSolverHybrid_ComputeResidual<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>,
00425 &AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>,
00426 initialGuess,
00427 PROBLEM_DIM * this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00428 this);
00429 }
00430
00431 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00432 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::SetNonlinearSolver(AbstractNonlinearSolver* pNonlinearSolver)
00433 {
00434 if (mWeAllocatedSolverMemory)
00435 {
00436 delete mpSolver;
00437 }
00438 mpSolver = pNonlinearSolver;
00439 mWeAllocatedSolverMemory = false;
00440 }
00441
00442 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00443 bool AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::VerifyJacobian(double tol)
00444 {
00445 unsigned size = PROBLEM_DIM * this->mpMesh->GetNumNodes();
00446
00447 Vec initial_guess = PetscTools::CreateAndSetVec(size, 0.0);
00448
00449 Mat analytic_jacobian;
00450 Mat numerical_jacobian;
00451
00452 PetscTools::SetupMat(analytic_jacobian, size, size, size);
00453 PetscTools::SetupMat(numerical_jacobian, size, size, size);
00454
00455 mUseAnalyticalJacobian = true;
00456 ComputeJacobian(initial_guess, &analytic_jacobian);
00457
00458 mUseAnalyticalJacobian = false;
00459 ComputeJacobian(initial_guess, &numerical_jacobian);
00460
00461 bool ok = PetscMatTools::CheckEquality(numerical_jacobian, analytic_jacobian, tol);
00462 PetscTools::Destroy(numerical_jacobian);
00463 PetscTools::Destroy(analytic_jacobian);
00464 PetscTools::Destroy(initial_guess);
00465
00466 return ok;
00467 }
00468
00469
00470
00471 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00472 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeResidual(SNES snes,
00473 Vec currentGuess,
00474 Vec residualVector,
00475 void* pContext)
00476 {
00477
00478 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* p_solver =
00479 (AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>*) pContext;
00480
00481 p_solver->ComputeResidual(currentGuess, residualVector);
00482
00483 return 0;
00484 }
00485
00486 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00487 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
00488 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian(SNES snes,
00489 Vec currentGuess,
00490 Mat globalJacobian,
00491 Mat preconditioner,
00492 void* pContext)
00493 {
00494
00495 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* p_solver =
00496 (AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>*) pContext;
00497
00498 p_solver->ComputeJacobian(currentGuess, &globalJacobian);
00499
00500 return 0;
00501 }
00502 #else
00503 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian(SNES snes,
00504 Vec currentGuess,
00505 Mat* pGlobalJacobian,
00506 Mat* pPreconditioner,
00507 MatStructure* pMatStructure,
00508 void* pContext)
00509 {
00510
00511 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* p_solver =
00512 (AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>*) pContext;
00513
00514 p_solver->ComputeJacobian(currentGuess, pGlobalJacobian);
00515
00516 return 0;
00517 }
00518 #endif
00519
00520 #endif