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 #ifndef ABSTRACTNONLINEARASSEMBLERSOLVERHYBRID_HPP_
00030 #define ABSTRACTNONLINEARASSEMBLERSOLVERHYBRID_HPP_
00031
00032
00033
00034 #include "BoundaryConditionsContainer.hpp"
00035 #include "AbstractNonlinearEllipticPde.hpp"
00036 #include "AbstractNonlinearSolver.hpp"
00037 #include "PetscTools.hpp"
00038 #include "AbstractFeObjectAssembler.hpp"
00039 #include "LinearSystem.hpp"
00040 #include "SimplePetscNonlinearSolver.hpp"
00041
00042
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
00075
00076
00093 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00094 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian(SNES snes,
00095 Vec currentGuess,
00096 Mat* pGlobalJacobian,
00097 Mat* pPreconditioner,
00098 MatStructure* pMatStructure,
00099 void* pContext);
00100
00101
00102
00103
00104
00105
00106
00107
00117
00118
00119
00124 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00125 class AbstractNonlinearAssemblerSolverHybrid : public AbstractFeObjectAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,true,true,NONLINEAR>
00126 {
00127 protected:
00129 AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* mpMesh;
00130
00132 BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* mpBoundaryConditions;
00133
00135 AbstractNonlinearSolver* mpSolver;
00136
00138 bool mWeAllocatedSolverMemory;
00139
00141 bool mUseAnalyticalJacobian;
00142
00143
00144
00154 void ApplyDirichletConditions(Vec currentGuess, Vec residual, Mat* pMat);
00155
00163 void ComputeJacobianNumerically(const Vec currentGuess, Mat* pJacobian);
00164
00165 public :
00166
00176 void ComputeResidual(const Vec currentGuess, Vec residualVector);
00177
00189 void ComputeJacobian(const Vec currentGuess, Mat* pJacobian);
00190
00191
00192
00200 AbstractNonlinearAssemblerSolverHybrid(AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* pMesh,
00201 BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* pBoundaryConditions,
00202 unsigned numQuadPoints = 2);
00203
00207 virtual ~AbstractNonlinearAssemblerSolverHybrid();
00208
00217 virtual Vec Solve(Vec initialGuess, bool useAnalyticalJacobian=false);
00218
00226 void SetNonlinearSolver(AbstractNonlinearSolver* pNonlinearSolver);
00227
00228
00242 bool VerifyJacobian(double tol);
00243 };
00244
00245
00246 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00247 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::AbstractNonlinearAssemblerSolverHybrid(
00248 AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* pMesh,
00249 BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* pBoundaryConditions,
00250 unsigned numQuadPoints)
00251 : AbstractFeObjectAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,true,true,NONLINEAR>(pMesh,numQuadPoints),
00252 mpMesh(pMesh),
00253 mpBoundaryConditions(pBoundaryConditions)
00254 {
00255
00256
00257 assert(SPACE_DIM==ELEMENT_DIM);
00258 assert(pMesh!=NULL);
00259 assert(pBoundaryConditions!=NULL);
00260
00261
00262 mpSolver = new SimplePetscNonlinearSolver;
00263 mWeAllocatedSolverMemory = true;
00264
00265 assert(mpMesh->GetNumNodes() == mpMesh->GetDistributedVectorFactory()->GetProblemSize());
00266 mpMesh->SetElementOwnerships(mpMesh->GetDistributedVectorFactory()->GetLow(),
00267 mpMesh->GetDistributedVectorFactory()->GetHigh());
00268 }
00269
00270
00271 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00272 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::~AbstractNonlinearAssemblerSolverHybrid()
00273 {
00274 if (mWeAllocatedSolverMemory)
00275 {
00276 delete mpSolver;
00277 }
00278 }
00279
00280 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00281 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ApplyDirichletConditions(Vec currentGuess, Vec residual, Mat* pJacobian)
00282 {
00283 if (residual)
00284 {
00285 mpBoundaryConditions->ApplyDirichletToNonlinearResidual( currentGuess,
00286 residual,
00287 *(mpMesh->GetDistributedVectorFactory()));
00288 }
00289 if (pJacobian)
00290 {
00291 mpBoundaryConditions->ApplyDirichletToNonlinearJacobian(*pJacobian);
00292 }
00293 }
00294
00295
00296 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00297 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ComputeResidual(const Vec currentGuess, Vec residualVector)
00298 {
00299 this->SetVectorToAssemble(residualVector,true);
00300 this->SetCurrentSolution(currentGuess);
00301 this->SetApplyNeummanBoundaryConditionsToVector(mpBoundaryConditions);
00302 this->AssembleVector();
00303
00304 VecAssemblyBegin(residualVector);
00305 VecAssemblyEnd(residualVector);
00306
00307 ApplyDirichletConditions(currentGuess, residualVector, NULL);
00308 }
00309
00310 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00311 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ComputeJacobian(const Vec currentGuess, Mat* pJacobian)
00312 {
00313 if (mUseAnalyticalJacobian)
00314 {
00315 this->SetMatrixToAssemble(*pJacobian);
00316 this->SetCurrentSolution(currentGuess);
00317 this->AssembleMatrix();
00318
00319 MatAssemblyBegin(*pJacobian, MAT_FLUSH_ASSEMBLY);
00320 MatAssemblyEnd(*pJacobian, MAT_FLUSH_ASSEMBLY);
00321
00322 ApplyDirichletConditions(currentGuess, NULL, pJacobian);
00323
00324 MatAssemblyBegin(*pJacobian, MAT_FINAL_ASSEMBLY);
00325 MatAssemblyEnd(*pJacobian, MAT_FINAL_ASSEMBLY);
00326 }
00327 else
00328 {
00329 ComputeJacobianNumerically(currentGuess, pJacobian);
00330 }
00331 }
00332
00333
00334 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00335 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ComputeJacobianNumerically(const Vec currentGuess, Mat* pJacobian)
00336 {
00337 unsigned num_unknowns = PROBLEM_DIM * this->mpMesh->GetNumNodes();
00338
00339
00340 Vec residual;
00341 Vec perturbed_residual;
00342 Vec result;
00343 residual = PetscTools::CreateVec(num_unknowns);
00344 result = PetscTools::CreateVec(num_unknowns);
00345 perturbed_residual = PetscTools::CreateVec(num_unknowns);
00346
00347
00348 Vec current_guess_copy;
00349 PETSCEXCEPT( VecDuplicate(currentGuess, ¤t_guess_copy) );
00350 PETSCEXCEPT( VecCopy(currentGuess, current_guess_copy) );
00351
00352
00353 ComputeResidual(currentGuess, residual);
00354
00355
00356 double h = 0.00001;
00357 PetscScalar subtract = -1;
00358 PetscScalar one_over_h = 1.0/h;
00359
00360 PetscInt ilo, ihi;
00361 VecGetOwnershipRange(current_guess_copy, &ilo, &ihi);
00362 unsigned lo = ilo;
00363 unsigned hi = ihi;
00364
00365
00366 for (unsigned global_index_outer = 0; global_index_outer < num_unknowns; global_index_outer++)
00367 {
00368
00369 if (lo<=global_index_outer && global_index_outer<hi)
00370 {
00371 PETSCEXCEPT( VecSetValue(current_guess_copy, global_index_outer,h, ADD_VALUES) );
00372 }
00373 ComputeResidual(current_guess_copy, perturbed_residual);
00374
00375
00376 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00377 PETSCEXCEPT( VecWAXPY(&subtract, residual, perturbed_residual, result) );
00378 PETSCEXCEPT( VecScale(&one_over_h, result) );
00379 #else
00380 PETSCEXCEPT( VecWAXPY(result, subtract, residual, perturbed_residual) );
00381 PETSCEXCEPT( VecScale(result, one_over_h) );
00382 #endif
00383
00384 double* p_result;
00385 PETSCEXCEPT( VecGetArray(result, &p_result) );
00386 for (unsigned global_index=lo; global_index < hi; global_index++)
00387 {
00388 unsigned local_index = global_index - lo;
00389 PETSCEXCEPT( MatSetValue(*pJacobian, global_index, global_index_outer,
00390 p_result[local_index], INSERT_VALUES) );
00391 }
00392 PETSCEXCEPT( VecRestoreArray(result, &p_result) );
00393
00394 if (lo<=global_index_outer && global_index_outer<hi)
00395 {
00396 PETSCEXCEPT( VecSetValue(current_guess_copy, global_index_outer, -h, ADD_VALUES) );
00397 }
00398 }
00399
00400 VecDestroy(residual);
00401 VecDestroy(perturbed_residual);
00402 VecDestroy(result);
00403 VecDestroy(current_guess_copy);
00404
00405 MatAssemblyBegin(*pJacobian, MAT_FINAL_ASSEMBLY);
00406 MatAssemblyEnd(*pJacobian, MAT_FINAL_ASSEMBLY);
00407 }
00408
00409
00410 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00411 Vec AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::Solve(Vec initialGuess,
00412 bool useAnalyticalJacobian)
00413 {
00414 assert(initialGuess!=NULL);
00415 mUseAnalyticalJacobian = useAnalyticalJacobian;
00416
00417 PetscInt size_of_init_guess;
00418 VecGetSize(initialGuess, &size_of_init_guess);
00419 PetscInt problem_size = PROBLEM_DIM * this->mpMesh->GetNumNodes();
00420 if (size_of_init_guess != problem_size)
00421 {
00422 std::stringstream error_message;
00423 error_message << "Size of initial guess vector, " << size_of_init_guess
00424 << ", does not match size of problem, " << problem_size;
00425 EXCEPTION(error_message.str());
00426 }
00427
00428
00429
00430 return mpSolver->Solve(&AbstractNonlinearAssemblerSolverHybrid_ComputeResidual<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>,
00431 &AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>,
00432 initialGuess,
00433 PROBLEM_DIM * this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00434 this);
00435 }
00436
00437
00438
00439
00440 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00441 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::SetNonlinearSolver(AbstractNonlinearSolver* pNonlinearSolver)
00442 {
00443 if (mWeAllocatedSolverMemory)
00444 {
00445 delete mpSolver;
00446 }
00447 mpSolver = pNonlinearSolver;
00448 mWeAllocatedSolverMemory = false;
00449 }
00450
00451
00452 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00453 bool AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::VerifyJacobian(double tol)
00454 {
00455 unsigned size = PROBLEM_DIM * this->mpMesh->GetNumNodes();
00456
00457 Vec initial_guess = PetscTools::CreateAndSetVec(size, 0.0);
00458
00459 Mat analytic_jacobian;
00460 Mat numerical_jacobian;
00461
00462 PetscTools::SetupMat(analytic_jacobian, size, size, size);
00463 PetscTools::SetupMat(numerical_jacobian, size, size, size);
00464
00465 mUseAnalyticalJacobian = true;
00466 ComputeJacobian(initial_guess, &analytic_jacobian);
00467
00468 mUseAnalyticalJacobian = false;
00469 ComputeJacobian(initial_guess, &numerical_jacobian);
00470
00471 double minus_one = -1.0;
00472 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
00473
00474 MatAYPX(&minus_one, analytic_jacobian, numerical_jacobian);
00475 #elif (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 3 && PETSC_VERSION_SUBMINOR == 1) //PETSc 2.3.1
00476
00477 MatAYPX(numerical_jacobian, minus_one, analytic_jacobian);
00478 #else
00479
00480 MatAYPX(numerical_jacobian, minus_one, analytic_jacobian, DIFFERENT_NONZERO_PATTERN);
00481 #endif
00482 double norm;
00483 MatNorm(numerical_jacobian,NORM_INFINITY,&norm);
00484
00485 MatDestroy(numerical_jacobian);
00486 MatDestroy(analytic_jacobian);
00487 VecDestroy(initial_guess);
00488
00489 return (norm<tol);
00490 }
00491
00492
00493
00494
00504
00505
00506 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00507 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeResidual(SNES snes, Vec currentGuess,
00508 Vec residualVector, 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->ComputeResidual(currentGuess, residualVector);
00515
00516 return 0;
00517 }
00518
00519
00520
00521
00522 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00523 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian(SNES snes, Vec currentGuess,
00524 Mat* pGlobalJacobian, Mat* pPreconditioner,
00525 MatStructure* pMatStructure, void* pContext)
00526 {
00527
00528 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* p_solver =
00529 (AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>*) pContext;
00530
00531 p_solver->ComputeJacobian(currentGuess, pGlobalJacobian);
00532
00533 return 0;
00534 }
00535
00536
00537 #endif