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
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
00179 AbstractNonlinearAssemblerSolverHybrid(AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* pMesh,
00180 BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* pBoundaryConditions);
00181
00185 virtual ~AbstractNonlinearAssemblerSolverHybrid();
00186
00196 virtual Vec Solve(Vec initialGuess, bool useAnalyticalJacobian=false);
00197
00205 void SetNonlinearSolver(AbstractNonlinearSolver* pNonlinearSolver);
00206
00220 bool VerifyJacobian(double tol);
00221 };
00222
00223 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00224 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::AbstractNonlinearAssemblerSolverHybrid(
00225 AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>* pMesh,
00226 BoundaryConditionsContainer<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* pBoundaryConditions)
00227 : AbstractFeVolumeIntegralAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM,true,true,NONLINEAR>(pMesh),
00228 mpMesh(pMesh),
00229 mpBoundaryConditions(pBoundaryConditions)
00230 {
00231
00232
00233
00234
00235 assert(SPACE_DIM==ELEMENT_DIM);
00236 assert(pMesh!=NULL);
00237 assert(pBoundaryConditions!=NULL);
00238
00239 mpSolver = new SimplePetscNonlinearSolver;
00240 mWeAllocatedSolverMemory = true;
00241
00242 assert(mpMesh->GetNumNodes() == mpMesh->GetDistributedVectorFactory()->GetProblemSize());
00243
00244 mpNeumannSurfaceTermsAssembler = new NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>(pMesh,pBoundaryConditions);
00245 }
00246
00247 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00248 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::~AbstractNonlinearAssemblerSolverHybrid()
00249 {
00250 if (mWeAllocatedSolverMemory)
00251 {
00252 delete mpSolver;
00253 }
00254 delete mpNeumannSurfaceTermsAssembler;
00255 }
00256
00257 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00258 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ApplyDirichletConditions(Vec currentGuess, Vec residual, Mat* pJacobian)
00259 {
00260 if (residual)
00261 {
00262 mpBoundaryConditions->ApplyDirichletToNonlinearResidual(currentGuess,
00263 residual,
00264 *(mpMesh->GetDistributedVectorFactory()));
00265 }
00266 if (pJacobian)
00267 {
00268 mpBoundaryConditions->ApplyDirichletToNonlinearJacobian(*pJacobian);
00269 }
00270 }
00271
00272 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00273 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ComputeResidual(const Vec currentGuess, Vec residualVector)
00274 {
00275
00276
00277 this->SetVectorToAssemble(residualVector,true);
00278 this->SetCurrentSolution(currentGuess);
00279 this->AssembleVector();
00280
00281
00282
00283
00284 mpNeumannSurfaceTermsAssembler->SetVectorToAssemble(residualVector, false);
00285 mpNeumannSurfaceTermsAssembler->SetScaleFactor(-1.0);
00286 mpNeumannSurfaceTermsAssembler->Assemble();
00287
00288 PetscVecTools::Finalise(residualVector);
00289
00290 ApplyDirichletConditions(currentGuess, residualVector, NULL);
00291 }
00292
00293 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00294 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ComputeJacobian(const Vec currentGuess, Mat* pJacobian)
00295 {
00296 if (mUseAnalyticalJacobian)
00297 {
00298 this->SetMatrixToAssemble(*pJacobian);
00299 this->SetCurrentSolution(currentGuess);
00300 this->AssembleMatrix();
00301
00302 PetscMatTools::SwitchWriteMode(*pJacobian);
00303
00304 ApplyDirichletConditions(currentGuess, NULL, pJacobian);
00305
00306 PetscMatTools::Finalise(*pJacobian);
00307 }
00308 else
00309 {
00310 ComputeJacobianNumerically(currentGuess, pJacobian);
00311 }
00312 }
00313
00314 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00315 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::ComputeJacobianNumerically(const Vec currentGuess, Mat* pJacobian)
00316 {
00317 unsigned num_unknowns = PROBLEM_DIM * this->mpMesh->GetNumNodes();
00318
00319
00320 Vec residual;
00321 Vec perturbed_residual;
00322 Vec result;
00323 residual = PetscTools::CreateVec(num_unknowns);
00324 result = PetscTools::CreateVec(num_unknowns);
00325 perturbed_residual = PetscTools::CreateVec(num_unknowns);
00326
00327
00328 Vec current_guess_copy;
00329 PETSCEXCEPT( VecDuplicate(currentGuess, ¤t_guess_copy) );
00330 PETSCEXCEPT( VecCopy(currentGuess, current_guess_copy) );
00331
00332
00333 ComputeResidual(currentGuess, residual);
00334
00335
00336 double h = 1e-5;
00337 double near_hsquared = 1e-9;
00338
00339 PetscInt ilo, ihi;
00340 VecGetOwnershipRange(current_guess_copy, &ilo, &ihi);
00341 unsigned lo = ilo;
00342 unsigned hi = ihi;
00343
00344
00345 for (unsigned global_index_outer = 0; global_index_outer < num_unknowns; global_index_outer++)
00346 {
00347
00348 PetscVecTools::AddToElement(current_guess_copy, global_index_outer, h);
00349 ComputeResidual(current_guess_copy, perturbed_residual);
00350
00351
00352 PetscVecTools::WAXPY(result, -1.0, residual, perturbed_residual);
00353 PetscVecTools::Scale(result, 1.0/h);
00354
00355 double* p_result;
00358 PETSCEXCEPT( VecGetArray(result, &p_result) );
00359 for (unsigned global_index=lo; global_index < hi; global_index++)
00360 {
00361 double result_entry = p_result[ global_index - lo];
00362 if (!CompareDoubles::IsNearZero(result_entry, near_hsquared))
00363 {
00364 PetscMatTools::SetElement(*pJacobian, global_index, global_index_outer, result_entry);
00365 }
00366 }
00367 PETSCEXCEPT( VecRestoreArray(result, &p_result) );
00368
00369 PetscVecTools::AddToElement(current_guess_copy, global_index_outer, -h);
00370 }
00371
00372 PetscTools::Destroy(residual);
00373 PetscTools::Destroy(perturbed_residual);
00374 PetscTools::Destroy(result);
00375 PetscTools::Destroy(current_guess_copy);
00376
00377 PetscMatTools::Finalise(*pJacobian);
00378 }
00379
00380 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00381 Vec AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::Solve(Vec initialGuess,
00382 bool useAnalyticalJacobian)
00383 {
00384 assert(initialGuess != NULL);
00385 mUseAnalyticalJacobian = useAnalyticalJacobian;
00386
00387 PetscInt size_of_init_guess;
00388 VecGetSize(initialGuess, &size_of_init_guess);
00389 PetscInt problem_size = PROBLEM_DIM * this->mpMesh->GetNumNodes();
00390 if (size_of_init_guess != problem_size)
00391 {
00392 EXCEPTION("Size of initial guess vector, " << size_of_init_guess
00393 << ", does not match size of problem, " << problem_size);
00394 }
00395
00396
00397 return mpSolver->Solve(&AbstractNonlinearAssemblerSolverHybrid_ComputeResidual<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>,
00398 &AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>,
00399 initialGuess,
00400 PROBLEM_DIM * this->mpMesh->CalculateMaximumNodeConnectivityPerProcess(),
00401 this);
00402 }
00403
00404 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00405 void AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::SetNonlinearSolver(AbstractNonlinearSolver* pNonlinearSolver)
00406 {
00407 if (mWeAllocatedSolverMemory)
00408 {
00409 delete mpSolver;
00410 }
00411 mpSolver = pNonlinearSolver;
00412 mWeAllocatedSolverMemory = false;
00413 }
00414
00415 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00416 bool AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>::VerifyJacobian(double tol)
00417 {
00418 unsigned size = PROBLEM_DIM * this->mpMesh->GetNumNodes();
00419
00420 Vec initial_guess = PetscTools::CreateAndSetVec(size, 0.0);
00421
00422 Mat analytic_jacobian;
00423 Mat numerical_jacobian;
00424
00425 PetscTools::SetupMat(analytic_jacobian, size, size, size);
00426 PetscTools::SetupMat(numerical_jacobian, size, size, size);
00427
00428 mUseAnalyticalJacobian = true;
00429 ComputeJacobian(initial_guess, &analytic_jacobian);
00430
00431 mUseAnalyticalJacobian = false;
00432 ComputeJacobian(initial_guess, &numerical_jacobian);
00433
00434 bool ok = PetscMatTools::CheckEquality(numerical_jacobian, analytic_jacobian, tol);
00435 PetscTools::Destroy(numerical_jacobian);
00436 PetscTools::Destroy(analytic_jacobian);
00437 PetscTools::Destroy(initial_guess);
00438
00439 return ok;
00440 }
00441
00442
00443
00444 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00445 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeResidual(SNES snes,
00446 Vec currentGuess,
00447 Vec residualVector,
00448 void* pContext)
00449 {
00450
00451 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* p_solver =
00452 (AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>*) pContext;
00453
00454 p_solver->ComputeResidual(currentGuess, residualVector);
00455
00456 return 0;
00457 }
00458
00459 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
00460 PetscErrorCode AbstractNonlinearAssemblerSolverHybrid_ComputeJacobian(SNES snes,
00461 Vec currentGuess,
00462 Mat* pGlobalJacobian,
00463 Mat* pPreconditioner,
00464 MatStructure* pMatStructure,
00465 void* pContext)
00466 {
00467
00468 AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>* p_solver =
00469 (AbstractNonlinearAssemblerSolverHybrid<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM>*) pContext;
00470
00471 p_solver->ComputeJacobian(currentGuess, pGlobalJacobian);
00472
00473 return 0;
00474 }
00475
00476 #endif