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 #ifndef _ABSTRACTNONLINEARASSEMBLER_HPP_
00029 #define _ABSTRACTNONLINEARASSEMBLER_HPP_
00030
00036 #include <vector>
00037
00038 #include "AbstractStaticAssembler.hpp"
00039 #include "BoundaryConditionsContainer.hpp"
00040 #include "AbstractNonlinearSolver.hpp"
00041 #include "AbstractNonlinearSolver.hpp"
00042 #include "SimplePetscNonlinearSolver.hpp"
00043 #include "RandomNumberGenerator.hpp"
00044 #include "PetscTools.hpp"
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00061 PetscErrorCode AbstractNonlinearAssembler_AssembleResidual(SNES snes,
00062 Vec currentGuess,
00063 Vec residualVector,
00064 void* pContext);
00065
00066 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00067 PetscErrorCode AbstractNonlinearAssembler_AssembleJacobian(SNES snes,
00068 Vec currentGuess,
00069 Mat* pGlobalJacobian,
00070 Mat* pPreconditioner,
00071 MatStructure* pMatStructure,
00072 void* pContext);
00073
00074
00078 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00079 class AbstractNonlinearAssembler : public AbstractStaticAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE>
00080 {
00081 typedef AbstractStaticAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE> BaseClassType;
00083 private:
00084
00086 Vec mInitialGuess;
00087
00088 protected:
00089
00091 AbstractNonlinearSolver* mpSolver;
00092
00094 bool mWeAllocatedSolverMemory;
00095
00097 bool mUseAnalyticalJacobian;
00098
00105 void ApplyDirichletConditions(Vec currentGuess, bool applyToMatrix);
00106
00107 public :
00108
00119 PetscErrorCode AssembleResidual(const Vec currentGuess, Vec residualVector);
00120
00133 PetscErrorCode AssembleJacobian(const Vec currentGuess, Mat* pGlobalJacobian);
00134
00135 protected:
00136
00146 PetscErrorCode AssembleJacobianNumerically(const Vec currentGuess, Mat* pJacobian);
00147
00170 virtual void AssembleSystem(bool assembleVector, bool assembleMatrix,
00171 Vec currentGuess=NULL, double currentTime=0.0);
00172
00176 bool ProblemIsNonlinear();
00177
00184 void InitialiseForSolve(Vec initialGuess);
00185
00196 Vec StaticSolve(Vec currentSolutionOrGuess=NULL,
00197 double currentTime=0.0,
00198 bool assembleMatrix=true);
00199
00200 public:
00201
00207 AbstractNonlinearAssembler(unsigned numQuadPoints = 2);
00208
00212 ~AbstractNonlinearAssembler();
00213
00225 void SetUseAnalyticalJacobian(bool useAnalyticalJacobian);
00226
00235 virtual Vec Solve(Vec initialGuess, bool useAnalyticalJacobian=false);
00236
00244 void SetNonlinearSolver(AbstractNonlinearSolver* pNonlinearSolver);
00245
00251 Vec CreateConstantInitialGuess(double value);
00252
00268 bool VerifyJacobian(double tol=1e-4);
00269
00270 };
00271
00272
00274
00276
00277
00278 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00279 void AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>::ApplyDirichletConditions(Vec currentGuess, bool applyToMatrix)
00280 {
00281 Vec& residual = this->mpLinearSystem->rGetRhsVector();
00282 Mat& jacobian = this->mpLinearSystem->rGetLhsMatrix();
00283 assert((jacobian && applyToMatrix) || (!jacobian && !applyToMatrix));
00284
00285 if (residual)
00286 {
00287 this->mpBoundaryConditions->ApplyDirichletToNonlinearResidual(
00288 currentGuess, residual, *(this->mpMesh->GetDistributedVectorFactory()));
00289 }
00290 if (jacobian)
00291 {
00292 this->mpBoundaryConditions->ApplyDirichletToNonlinearJacobian(jacobian);
00293 }
00294 }
00295
00296
00297 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00298 PetscErrorCode AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>::AssembleResidual(const Vec currentGuess, Vec residualVector)
00299 {
00300
00301
00302 this->PrepareForSolve();
00303 delete this->mpLinearSystem;
00304 this->mpLinearSystem = new LinearSystem(residualVector, NULL);
00305 this->AssembleSystem(true, false, currentGuess, 0.0);
00306 return 0;
00307 }
00308
00309
00310 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00311 PetscErrorCode AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>::AssembleJacobian(const Vec currentGuess, Mat* pGlobalJacobian)
00312 {
00313 if (mUseAnalyticalJacobian)
00314 {
00315
00316
00317 delete this->mpLinearSystem;
00318 this->mpLinearSystem = new LinearSystem(NULL,* pGlobalJacobian);
00319 this->AssembleSystem(false, true, currentGuess, 0.0);
00320 return 0;
00321 }
00322 else
00323 {
00324 return AssembleJacobianNumerically(currentGuess, pGlobalJacobian);
00325 }
00326 }
00327
00328
00329 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00330 PetscErrorCode AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>::AssembleJacobianNumerically(const Vec currentGuess, Mat* pJacobian)
00331 {
00332 unsigned num_nodes = PROBLEM_DIM*this->mpMesh->GetNumNodes();
00333
00334
00335 Vec residual;
00336 Vec perturbed_residual;
00337 Vec result;
00338
00339 VecCreate(PETSC_COMM_WORLD, &residual);
00340 VecCreate(PETSC_COMM_WORLD, &result);
00341 VecCreate(PETSC_COMM_WORLD, &perturbed_residual);
00342
00343 VecSetSizes(residual, PETSC_DECIDE, num_nodes);
00344 VecSetSizes(result, PETSC_DECIDE, num_nodes);
00345 VecSetSizes(perturbed_residual, PETSC_DECIDE, num_nodes);
00346
00347 VecSetFromOptions(residual);
00348 VecSetFromOptions(result);
00349 VecSetFromOptions(perturbed_residual);
00350
00351
00352 Vec current_guess_copy;
00353 PETSCEXCEPT( VecDuplicate(currentGuess, ¤t_guess_copy) );
00354 PETSCEXCEPT( VecCopy(currentGuess, current_guess_copy) );
00355
00356
00357 AssembleResidual(currentGuess, residual);
00358
00359
00360 double h = 0.00001;
00361 PetscScalar subtract = -1;
00362 PetscScalar one_over_h = 1.0/h;
00363
00364 PetscInt ilo, ihi;
00365 VecGetOwnershipRange(current_guess_copy, &ilo, &ihi);
00366 unsigned lo = ilo;
00367 unsigned hi = ihi;
00368
00369
00370 for (unsigned global_index_outer = 0; global_index_outer < num_nodes; global_index_outer++)
00371 {
00372
00373 if (lo<=global_index_outer && global_index_outer<hi)
00374 {
00375 PETSCEXCEPT( VecSetValue(current_guess_copy, global_index_outer,h, ADD_VALUES) );
00376 }
00377 AssembleResidual(current_guess_copy, perturbed_residual);
00378
00379
00380 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00381 PETSCEXCEPT( VecWAXPY(&subtract, residual, perturbed_residual, result) );
00382 PETSCEXCEPT( VecScale(&one_over_h, result) );
00383 #else
00384 PETSCEXCEPT( VecWAXPY(result, subtract, residual, perturbed_residual) );
00385 PETSCEXCEPT( VecScale(result, one_over_h) );
00386 #endif
00387
00388 double* p_result;
00389 PETSCEXCEPT( VecGetArray(result, &p_result) );
00390 for (unsigned global_index=lo; global_index < hi; global_index++)
00391 {
00392 unsigned local_index = global_index - lo;
00393 PETSCEXCEPT( MatSetValue(*pJacobian, global_index, global_index_outer,
00394 p_result[local_index], INSERT_VALUES) );
00395 }
00396 PETSCEXCEPT( VecRestoreArray(result, &p_result) );
00397
00398 if (lo<=global_index_outer && global_index_outer<hi)
00399 {
00400 PETSCEXCEPT( VecSetValue(current_guess_copy, global_index_outer, -h, ADD_VALUES) );
00401 }
00402 }
00403
00404 VecDestroy(residual);
00405 VecDestroy(perturbed_residual);
00406 VecDestroy(result);
00407 VecDestroy(current_guess_copy);
00408
00409 MatAssemblyBegin(*pJacobian, MAT_FINAL_ASSEMBLY);
00410 MatAssemblyEnd(*pJacobian, MAT_FINAL_ASSEMBLY);
00411
00412 return 0;
00413 }
00414
00415
00416 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00417 void AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>::AssembleSystem(bool assembleVector, bool assembleMatrix,
00418 Vec currentGuess, double currentTime)
00419 {
00420
00421 assert( currentGuess );
00422 BaseClassType::AssembleSystem(assembleVector, assembleMatrix, currentGuess, currentTime);
00423 }
00424
00425
00426 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00427 bool AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>::ProblemIsNonlinear()
00428 {
00429 return true;
00430 }
00431
00432
00433 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00434 void AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>::InitialiseForSolve(Vec initialGuess)
00435 {
00436
00437 PetscInt size_of_init_guess;
00438 VecGetSize(initialGuess, &size_of_init_guess);
00439 PetscInt problem_size = PROBLEM_DIM * this->mpMesh->GetNumNodes();
00440 if (size_of_init_guess != problem_size)
00441 {
00442 std::stringstream error_message;
00443 error_message << "Size of initial guess vector, " << size_of_init_guess
00444 << ", does not match size of problem, " << problem_size;
00445 EXCEPTION(error_message.str());
00446 }
00447 }
00448
00449
00450 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00451 Vec AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>::StaticSolve(Vec currentSolutionOrGuess,
00452 double currentTime,
00453 bool assembleMatrix)
00454 {
00455 assert(this->mpBoundaryConditions!=NULL);
00456 assert(currentSolutionOrGuess!=NULL);
00457 assert(assembleMatrix);
00458
00459
00460
00461 Vec answer = this->mpSolver->Solve(&AbstractNonlinearAssembler_AssembleResidual<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>,
00462 &AbstractNonlinearAssembler_AssembleJacobian<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>,
00463 currentSolutionOrGuess,
00464 this);
00465 return answer;
00466 }
00467
00468
00469
00470 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00471 AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>::AbstractNonlinearAssembler(unsigned numQuadPoints)
00472 : BaseClassType(numQuadPoints),
00473 mInitialGuess(NULL),
00474 mWeAllocatedSolverMemory(true),
00475 mUseAnalyticalJacobian(false)
00476 {
00477 mpSolver = new SimplePetscNonlinearSolver;
00478 }
00479
00480
00481 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00482 AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>::~AbstractNonlinearAssembler()
00483 {
00484 if (mWeAllocatedSolverMemory)
00485 {
00486 delete mpSolver;
00487 }
00488 }
00489
00490
00491 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00492 void AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>::SetUseAnalyticalJacobian(bool useAnalyticalJacobian)
00493 {
00494 mUseAnalyticalJacobian = useAnalyticalJacobian;
00495 }
00496
00497
00498 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00499 Vec AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>::Solve(Vec initialGuess, bool useAnalyticalJacobian)
00500 {
00501 assert(initialGuess!=NULL);
00502 SetUseAnalyticalJacobian(useAnalyticalJacobian);
00503
00504 this->PrepareForSolve();
00505 InitialiseForSolve(initialGuess);
00506
00507 return StaticSolve(initialGuess);
00508 }
00509
00510
00511 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00512 void AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>::SetNonlinearSolver(AbstractNonlinearSolver* pNonlinearSolver)
00513 {
00514 if (mWeAllocatedSolverMemory)
00515 {
00516 delete mpSolver;
00517 }
00518 mpSolver = pNonlinearSolver;
00519 mWeAllocatedSolverMemory = false;
00520 }
00521
00522
00523 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00524 Vec AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>::CreateConstantInitialGuess(double value)
00525 {
00526 assert(this->mpMesh!=NULL);
00527 unsigned size = PROBLEM_DIM * this->mpMesh->GetNumNodes();
00528 return PetscTools::CreateVec(size, value);
00529 }
00530
00531
00532 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00533 bool AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>::VerifyJacobian(double tol)
00534 {
00535 unsigned size = PROBLEM_DIM * this->mpMesh->GetNumNodes();
00536
00537 Vec initial_guess;
00538 VecCreate(PETSC_COMM_WORLD, &initial_guess);
00539 VecSetSizes(initial_guess, PETSC_DECIDE, size);
00540 VecSetFromOptions(initial_guess);
00541
00542 for (unsigned i=0; i<size; i++)
00543 {
00544 VecSetValue(initial_guess, i, 0.0, INSERT_VALUES);
00545 }
00546
00547 VecAssemblyBegin(initial_guess);
00548 VecAssemblyEnd(initial_guess);
00549
00550
00551 Mat analytic_jacobian;
00552 Mat numerical_jacobian;
00553
00554 PetscTools::SetupMat(analytic_jacobian, size, size);
00555 PetscTools::SetupMat(numerical_jacobian, size, size);
00556
00557 mUseAnalyticalJacobian = true;
00558 AssembleJacobian(initial_guess, &analytic_jacobian);
00559
00560 mUseAnalyticalJacobian = false;
00561 AssembleJacobian(initial_guess, &numerical_jacobian);
00562
00563 bool all_less_than_tol = true;
00564
00565 for (unsigned i=0; i<size; i++)
00566 {
00567 for (unsigned j=0; j<size; j++)
00568 {
00569 double val_a[1];
00570 double val_n[1];
00571 PetscInt row[1];
00572 PetscInt col[1];
00573 row[0] = i;
00574 col[0] = j;
00576 MatGetValues(numerical_jacobian,1,row,1,col,val_n);
00577 MatGetValues(analytic_jacobian,1,row,1,col,val_a);
00578
00579 if(fabs(val_n[0]-val_a[0]) > tol)
00580 {
00581 #define COVERAGE_IGNORE // would have to write a bad concrete assembler class just to cover this line
00582 all_less_than_tol = false;
00583 #undef COVERAGE_IGNORE
00584 }
00585 }
00586 }
00587
00588 MatDestroy(numerical_jacobian);
00589 MatDestroy(analytic_jacobian);
00590 VecDestroy(initial_guess);
00591
00592 return all_less_than_tol;
00593 }
00594
00595
00596
00597
00598
00599
00600
00601
00602
00603
00604
00618 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00619 PetscErrorCode AbstractNonlinearAssembler_AssembleResidual(SNES snes, Vec currentGuess,
00620 Vec residualVector, void* pContext)
00621 {
00622
00623 AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>* pAssembler =
00624 (AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>*) pContext;
00625
00626 PetscErrorCode ierr = pAssembler->AssembleResidual(currentGuess, residualVector);
00627
00628
00629
00630
00631
00632 return ierr;
00633 }
00634
00635
00636
00652 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
00653 PetscErrorCode AbstractNonlinearAssembler_AssembleJacobian(SNES snes, Vec currentGuess,
00654 Mat* pGlobalJacobian, Mat* pPreconditioner,
00655 MatStructure* pMatStructure, void* pContext)
00656 {
00657
00658
00659
00660 AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>* pAssembler =
00661 (AbstractNonlinearAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE>*) pContext;
00662
00663 PetscErrorCode ierr = pAssembler->AssembleJacobian(currentGuess, pGlobalJacobian);
00664
00665
00666 return ierr;
00667 }
00668
00669
00670
00671 #endif //_ABSTRACTNONLINEARASSEMBLER_HPP_