AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE > Class Template Reference

#include <AbstractNonlinearAssembler.hpp>

Inheritance diagram for AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >:

Inheritance graph
[legend]
Collaboration diagram for AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >:

Collaboration graph
[legend]

List of all members.

Public Member Functions

PetscErrorCode AssembleResidual (const Vec currentGuess, Vec residualVector)
PetscErrorCode AssembleJacobian (const Vec currentGuess, Mat *pGlobalJacobian)
 AbstractNonlinearAssembler (unsigned numQuadPoints=2)
 ~AbstractNonlinearAssembler ()
void SetUseAnalyticalJacobian (bool useAnalyticalJacobian)
virtual Vec Solve (Vec initialGuess, bool useAnalyticalJacobian=false)
void SetNonlinearSolver (AbstractNonlinearSolver *pNonlinearSolver)
Vec CreateConstantInitialGuess (double value)
bool VerifyJacobian (double tol=1e-4, bool print=false)

Protected Member Functions

void ApplyDirichletConditions (Vec currentGuess, bool applyToMatrix)
PetscErrorCode AssembleJacobianNumerically (const Vec currentGuess, Mat *pJacobian)
virtual void AssembleSystem (bool assembleVector, bool assembleMatrix, Vec currentGuess=NULL, double currentTime=0.0)
bool ProblemIsNonlinear ()
void InitialiseForSolve (Vec initialGuess)
Vec StaticSolve (Vec currentSolutionOrGuess=NULL, double currentTime=0.0, bool assembleMatrix=true)

Protected Attributes

AbstractNonlinearSolvermpSolver
bool mWeAllocatedSolverMemory
bool mUseAnalyticalJacobian

Private Types

typedef
AbstractStaticAssembler
< ELEMENT_DIM, SPACE_DIM,
PROBLEM_DIM, true, CONCRETE > 
BaseClassType

Private Attributes

Vec mInitialGuess


Detailed Description

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
class AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >

AbstractNonlinearAssembler

Definition at line 79 of file AbstractNonlinearAssembler.hpp.


Member Typedef Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
typedef AbstractStaticAssembler<ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE> AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::BaseClassType [private]

Base class type (to save typing).

Reimplemented in SimpleNonlinearEllipticAssembler< ELEMENT_DIM, SPACE_DIM >.

Definition at line 81 of file AbstractNonlinearAssembler.hpp.


Constructor & Destructor Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::AbstractNonlinearAssembler ( unsigned  numQuadPoints = 2  )  [inline]

Constructors just call the base class versions.

Parameters:
numQuadPoints number of quadrature points (defaults to 2)

Definition at line 472 of file AbstractNonlinearAssembler.hpp.

References AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::mpSolver.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::~AbstractNonlinearAssembler (  )  [inline]


Member Function Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
void AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::ApplyDirichletConditions ( Vec  currentGuess,
bool  applyToMatrix 
) [inline, protected, virtual]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
PetscErrorCode AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::AssembleResidual ( const Vec  currentGuess,
Vec  residualVector 
) [inline]

Compute the residual vector given the current solution guess.

Parameters:
currentGuess The solution guess for the current iteration.
residualVector We fill this with the residual vector.
Returns:
An error code if any PETSc routines fail.
NOTE: this method is called indirectly by the PETSc iterative solvers, so must be public.

Definition at line 299 of file AbstractNonlinearAssembler.hpp.

References AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::AssembleSystem(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE >::mpLinearSystem, and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE >::PrepareForSolve().

Referenced by AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::AssembleJacobianNumerically().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
PetscErrorCode AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::AssembleJacobian ( const Vec  currentGuess,
Mat *  pGlobalJacobian 
) [inline]

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:
currentGuess The solution guess for the current iteration.
pGlobalJacobian Pointer to object to fill with the jacobian matrix.
Returns:
An error code if any PETSc routines fail.
NOTE: this method is called indirectly by the PETSc iterative solvers, so must be public.

Definition at line 312 of file AbstractNonlinearAssembler.hpp.

References AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::AssembleJacobianNumerically(), AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::AssembleSystem(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE >::mpLinearSystem, and AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::mUseAnalyticalJacobian.

Referenced by AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::VerifyJacobian().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
PetscErrorCode AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::AssembleJacobianNumerically ( const Vec  currentGuess,
Mat *  pJacobian 
) [inline, protected]

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

Parameters:
currentGuess Independent variable, u in f(u), for example
pJacobian A pointer to the Jacobian matrix
Returns:
An error code if any PETSc routines fail.

Definition at line 331 of file AbstractNonlinearAssembler.hpp.

References AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::AssembleResidual(), and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE >::mpMesh.

Referenced by AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::AssembleJacobian().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
void AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::AssembleSystem ( bool  assembleVector,
bool  assembleMatrix,
Vec  currentGuess = NULL,
double  currentTime = 0.0 
) [inline, protected, virtual]

AssembleSystem - the major method for all assemblers

Assemble the linear system for a linear PDE, or the residual or Jacobian for nonlinear PDEs. Loops over each element (and each each surface element if there are non-zero Neumann boundary conditions), calls AssembleOnElement() and adds the contribution to the linear system.

Takes in current solution and time if necessary but only used if the problem is a dynamic one. This method uses PROBLEM_DIM and can assemble linear systems for any number of unknown variables.

Parameters:
assembleVector Whether to assemble the RHS vector of the linear system (i.e. the residual vector for nonlinear problems).
assembleMatrix Whether to assemble the LHS matrix of the linear system (i.e. the jacobian matrix for nonlinear problems).
currentGuess The current solution in a linear dynamic problem, or the current guess in a nonlinear problem. Should be NULL for linear static problems. Defaults to NULL.
currentTime The current time for dynamic problems. Not used in static problems. Defaults to 0.0.

Reimplemented from AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE >.

Definition at line 418 of file AbstractNonlinearAssembler.hpp.

References AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE >::AssembleSystem().

Referenced by AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::AssembleJacobian(), and AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::AssembleResidual().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
bool AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::ProblemIsNonlinear (  )  [inline, protected, virtual]

Whether grad_u should be calculated

Implements AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 428 of file AbstractNonlinearAssembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
void AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::InitialiseForSolve ( Vec  initialGuess  )  [inline, protected, virtual]

No separate initialisation is needed in the nonlinear case; PrepareForSolve does everything. We just check the size of the initial guess.

Parameters:
initialGuess an initial guess

Implements AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 435 of file AbstractNonlinearAssembler.hpp.

References AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE >::mpMesh.

Referenced by AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::Solve().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
Vec AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::StaticSolve ( Vec  currentSolutionOrGuess = NULL,
double  currentTime = 0.0,
bool  assembleMatrix = true 
) [inline, protected, virtual]

Perform the work of a single solve, but without any initialisation.

Parameters:
currentSolutionOrGuess either the current solution (dynamic problem) or initial guess (static problem). MUST be provided.
currentTime for a dynamic problem, the current time
assembleMatrix Whether to assemble the LHS matrix of the linear system (i.e. the jacobian matrix for nonlinear problems).
Returns:
the solution vector

Todo:
do something sensible if assembleMatrix is false.

Implements AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >.

Definition at line 452 of file AbstractNonlinearAssembler.hpp.

References AbstractAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpBoundaryConditions, AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::mpSolver, and AbstractNonlinearSolver::Solve().

Referenced by AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::Solve().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
void AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::SetUseAnalyticalJacobian ( bool  useAnalyticalJacobian  )  [inline]

Set whether to use an analytical jacobian. This is provided for use when solving dynamic nonlinear problems; when solving static problems there is an argument to the Solve method which specifies this property, and overrides any user call to this method.

If this method is not called the default is false i.e. numerical jacobian.

Parameters:
useAnalyticalJacobian Set to true to use an analytically calculated jacobian matrix rather than a numerically approximated one.

Definition at line 493 of file AbstractNonlinearAssembler.hpp.

References AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::mUseAnalyticalJacobian.

Referenced by AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::Solve().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
Vec AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::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 500 of file AbstractNonlinearAssembler.hpp.

References AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::InitialiseForSolve(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE >::PrepareForSolve(), AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::SetUseAnalyticalJacobian(), and AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::StaticSolve().

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
void AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::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 513 of file AbstractNonlinearAssembler.hpp.

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

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
Vec AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::CreateConstantInitialGuess ( double  value  )  [inline]

A helpful method for creating an initial guess vector.

Parameters:
value constant value of the initial guess

Definition at line 525 of file AbstractNonlinearAssembler.hpp.

References PetscTools::CreateVec(), and AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE >::mpMesh.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
bool AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::VerifyJacobian ( double  tol = 1e-4,
bool  print = false 
) [inline]

VerifyJacobian

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 (which only calls ComputeVectorTerm and ComputeVectorSurfaceTerm) to verify the correctness of the code.

Parameters:
tol A tolerance which defaults to 1e-5
print Whether to print the different matrix J_numerical-J_analytical
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.

Definition at line 534 of file AbstractNonlinearAssembler.hpp.

References AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::AssembleJacobian(), AbstractStaticAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, CONCRETE >::mpMesh, AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::mUseAnalyticalJacobian, and PetscTools::SetupMat().


Member Data Documentation

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
Vec AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::mInitialGuess [private]

PETSc vector storing an initial guess for the solution.

Definition at line 86 of file AbstractNonlinearAssembler.hpp.

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
AbstractNonlinearSolver* AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::mpSolver [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
bool AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::mWeAllocatedSolverMemory [protected]

template<unsigned ELEMENT_DIM, unsigned SPACE_DIM, unsigned PROBLEM_DIM, class CONCRETE>
bool AbstractNonlinearAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CONCRETE >::mUseAnalyticalJacobian [protected]


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

Generated on Tue Aug 4 16:10:41 2009 for Chaste by  doxygen 1.5.5