Chaste Release::3.1
AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > Class Template Reference

#include <AbstractNonlinearAssemblerSolverHybrid.hpp>

Inheritance diagram for AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >:
Collaboration diagram for AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >:

List of all members.

Public Member Functions

void ComputeResidual (const Vec currentGuess, Vec residualVector)
void ComputeJacobian (const Vec currentGuess, Mat *pJacobian)
 AbstractNonlinearAssemblerSolverHybrid (AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *pMesh, BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > *pBoundaryConditions, unsigned numQuadPoints=2)
virtual ~AbstractNonlinearAssemblerSolverHybrid ()
virtual Vec Solve (Vec initialGuess, bool useAnalyticalJacobian=false)
void SetNonlinearSolver (AbstractNonlinearSolver *pNonlinearSolver)
bool VerifyJacobian (double tol)

Protected Member Functions

void ApplyDirichletConditions (Vec currentGuess, Vec residual, Mat *pMat)
void ComputeJacobianNumerically (const Vec currentGuess, Mat *pJacobian)

Protected Attributes

AbstractTetrahedralMesh
< ELEMENT_DIM, SPACE_DIM > * 
mpMesh
BoundaryConditionsContainer
< ELEMENT_DIM, SPACE_DIM,
PROBLEM_DIM > * 
mpBoundaryConditions
AbstractNonlinearSolvermpSolver
bool mWeAllocatedSolverMemory
bool mUseAnalyticalJacobian
NaturalNeumannSurfaceTermAssembler
< ELEMENT_DIM, SPACE_DIM,
PROBLEM_DIM > * 
mpNeumannSurfaceTermsAssembler

Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
class AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >

The ASSEMBLER-SOLVER class.

An abstract solver for solving nonlinear elliptic PDEs.

Definition at line 105 of file AbstractNonlinearAssemblerSolverHybrid.hpp.


Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AbstractNonlinearAssemblerSolverHybrid ( AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > *  pMesh,
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > *  pBoundaryConditions,
unsigned  numQuadPoints = 2 
)
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::~AbstractNonlinearAssemblerSolverHybrid ( ) [virtual]

Destructor.

Definition at line 252 of file AbstractNonlinearAssemblerSolverHybrid.hpp.


Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ApplyDirichletConditions ( Vec  currentGuess,
Vec  residual,
Mat pMat 
) [protected]

Apply Dirichlet boundary conditions to either the residual or Jacobian.

Parameters:
currentGuessthe solution guess for the current nonlinear solver iteration
residualthe residual to apply boundary conditions to (can be NULL if bcs to be applied to Jacobian only)
pMatthe Jacobian to apply boundary conditions to (can be NULL if bcs to be applied to residual only)

Definition at line 262 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobian ( const Vec  currentGuess,
Mat pJacobian 
)

Compute the Jacobian matrix given a current guess at the solution. Choose whether to use a numerical or analytical method based on a flag provided by the user (in Solve()).

Parameters:
currentGuessThe solution guess for the current iteration.
pJacobianPointer to object to fill with the Jacobian matrix.

NOTE: this method is called indirectly by the PETSc iterative solvers, so must be public.

Definition at line 298 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

References PetscMatTools::Finalise(), and PetscMatTools::SwitchWriteMode().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobianNumerically ( const Vec  currentGuess,
Mat pJacobian 
) [protected]

Computes the Jacobian numerically i.e. an approximation, using numerical partial derivatives.

Parameters:
currentGuessIndependent variable, u in f(u), for example
pJacobianA pointer to the Jacobian matrix

Definition at line 319 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

References PetscVecTools::AddToElement(), PetscTools::CreateVec(), PetscTools::Destroy(), PetscMatTools::Finalise(), PetscVecTools::Scale(), PetscMatTools::SetElement(), and PetscVecTools::WAXPY().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeResidual ( const Vec  currentGuess,
Vec  residualVector 
)

Compute the residual vector given the current solution guess.

Parameters:
currentGuessThe solution guess for the current nonlinear solve iteration.
residualVectorWe fill this with the residual vector.

NOTE: this method is called indirectly by the PETSc iterative solvers, so must be public.

Definition at line 277 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

References PetscVecTools::Finalise().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetNonlinearSolver ( AbstractNonlinearSolver pNonlinearSolver)

SetNonlinearSolver - by default a SimplePetscNonlinearSolver is created and used in this class, this method can be called to use a different AbstractNonlinearSolver.

Parameters:
pNonlinearSolvera nonlinear solver

Definition at line 403 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
Vec AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::Solve ( Vec  initialGuess,
bool  useAnalyticalJacobian = false 
) [virtual]

Assemble and solve the system for a nonlinear elliptic PDE.

Parameters:
initialGuessAn initial guess for the iterative solver
useAnalyticalJacobianSet to true to use an analytically calculated Jacobian matrix rather than a numerically approximated one.
Returns:
A PETSc vector giving the solution at each mesh node.

Definition at line 379 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

References EXCEPTION.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
bool AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::VerifyJacobian ( double  tol)

A helper method for use when writing concrete assemblers. Once the user has calculated (on paper) the weak form and the form of the ComputeMatrixTerm method, they can check whether the analytic Jacobian matches the numerical Jacobian to verify the correctness of the code.

Parameters:
tolA tolerance which defaults to 1e-5
Returns:
true if the componentwise difference between the matrices is less than the tolerance, false otherwise.

This method should NOT be run in simulations - it is only to verify the correctness of the concrete assembler code. Allocates dense matrices.

Definition at line 414 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

References PetscMatTools::CheckEquality(), PetscTools::CreateAndSetVec(), PetscTools::Destroy(), and PetscTools::SetupMat().


Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
BoundaryConditionsContainer<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpBoundaryConditions [protected]

Boundary conditions container.

Definition at line 113 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh [protected]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
NaturalNeumannSurfaceTermAssembler<ELEMENT_DIM,SPACE_DIM,PROBLEM_DIM>* AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpNeumannSurfaceTermsAssembler [protected]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
AbstractNonlinearSolver* AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpSolver [protected]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
bool AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mUseAnalyticalJacobian [protected]

Whether to use an analytical expression for the Jacobian.

Definition at line 122 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
bool AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mWeAllocatedSolverMemory [protected]

The documentation for this class was generated from the following file: