AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > Class Template Reference

#include <AbstractNonlinearAssemblerSolverHybrid.hpp>

Inherits AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR >.

Collaboration diagram for AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >:
Collaboration graph
[legend]

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 98 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 
) [inline]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::~AbstractNonlinearAssemblerSolverHybrid (  )  [inline, virtual]

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 
) [inline, protected]

Apply Dirichlet boundary conditions to either the residual or Jacobian.

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

Definition at line 255 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

References AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpBoundaryConditions, and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh.

Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobian(), and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeResidual().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobian ( const Vec  currentGuess,
Mat *  pJacobian 
) [inline]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobianNumerically ( const Vec  currentGuess,
Mat *  pJacobian 
) [inline, protected]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeResidual ( const Vec  currentGuess,
Vec  residualVector 
) [inline]
template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM>
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::SetNonlinearSolver ( AbstractNonlinearSolver pNonlinearSolver  )  [inline]

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

Parameters:
pNonlinearSolver a nonlinear solver

Definition at line 396 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

References AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpSolver, and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mWeAllocatedSolverMemory.

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

Assemble and solve the system for a nonlinear elliptic PDE.

Parameters:
initialGuess An initial guess for the iterative solver
useAnalyticalJacobian Set 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 372 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

References EXCEPTION, AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh, AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpSolver, AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mUseAnalyticalJacobian, and AbstractNonlinearSolver::Solve().

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

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:
tol A 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 407 of file AbstractNonlinearAssemblerSolverHybrid.hpp.

References PetscMatTools::CheckEquality(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeJacobian(), PetscTools::CreateAndSetVec(), AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh, AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mUseAnalyticalJacobian, 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]
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]
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:
Generated on Thu Dec 22 13:01:18 2011 for Chaste by  doxygen 1.6.3