Chaste Commit::baa90ac2819b962188b7562f2326be23c47859a7
|
#include <AbstractNonlinearAssemblerSolverHybrid.hpp>
Protected Member Functions | |
void | ApplyDirichletConditions (Vec currentGuess, Vec residual, Mat *pMat) |
void | ComputeJacobianNumerically (const Vec currentGuess, Mat *pJacobian) |
Protected Member Functions inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR > | |
void | ComputeTransformedBasisFunctionDerivatives (const ChastePoint< ELEMENT_DIM > &rPoint, const c_matrix< double, ELEMENT_DIM, SPACE_DIM > &rInverseJacobian, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rReturnValue) |
void | DoAssemble () |
virtual c_matrix< double, PROBLEM_DIM *(ELEMENT_DIM+1), PROBLEM_DIM *(ELEMENT_DIM+1)> | ComputeMatrixTerm (c_vector< double, ELEMENT_DIM+1 > &rPhi, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, ChastePoint< SPACE_DIM > &rX, c_vector< double, PROBLEM_DIM > &rU, c_matrix< double, PROBLEM_DIM, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement) |
virtual c_vector< double, PROBLEM_DIM *(ELEMENT_DIM+1)> | ComputeVectorTerm (c_vector< double, ELEMENT_DIM+1 > &rPhi, c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, ChastePoint< SPACE_DIM > &rX, c_vector< double, PROBLEM_DIM > &rU, c_matrix< double, PROBLEM_DIM, SPACE_DIM > &rGradU, Element< ELEMENT_DIM, SPACE_DIM > *pElement) |
virtual void | AssembleOnElement (Element< ELEMENT_DIM, SPACE_DIM > &rElement, c_matrix< double, PROBLEM_DIM *(ELEMENT_DIM+1), PROBLEM_DIM *(ELEMENT_DIM+1) > &rAElem, c_vector< double, PROBLEM_DIM *(ELEMENT_DIM+1)> &rBElem) |
virtual bool | ElementAssemblyCriterion (Element< ELEMENT_DIM, SPACE_DIM > &rElement) |
Protected Member Functions inherited from AbstractFeAssemblerCommon< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, CAN_ASSEMBLE_VECTOR, CAN_ASSEMBLE_MATRIX, INTERPOLATION_LEVEL > | |
virtual double | GetCurrentSolutionOrGuessValue (unsigned nodeIndex, unsigned indexOfUnknown) |
virtual void | ResetInterpolatedQuantities () |
virtual void | IncrementInterpolatedQuantities (double phiI, const Node< SPACE_DIM > *pNode) |
virtual void | IncrementInterpolatedGradientQuantities (const c_matrix< double, SPACE_DIM, ELEMENT_DIM+1 > &rGradPhi, unsigned phiIndex, const Node< SPACE_DIM > *pNode) |
Additional Inherited Members | |
Protected Types inherited from AbstractFeVolumeIntegralAssembler< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM, true, true, NONLINEAR > | |
typedef LinearBasisFunction< ELEMENT_DIM > | BasisFunction |
The ASSEMBLER-SOLVER class.
An abstract solver for solving nonlinear elliptic PDEs.
Definition at line 132 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AbstractNonlinearAssemblerSolverHybrid | ( | AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * | pMesh, |
BoundaryConditionsContainer< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM > * | pBoundaryConditions | ||
) |
Constructor.
pMesh | The mesh |
pBoundaryConditions | The boundary conditions to apply |
Definition at line 251 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpMesh, AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpNeumannSurfaceTermsAssembler, AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mpSolver, and AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::mWeAllocatedSolverMemory.
|
virtual |
Destructor.
Definition at line 275 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
|
protected |
Apply Dirichlet boundary conditions to either the residual or Jacobian.
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 285 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
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()).
currentGuess | The solution guess for the current iteration. |
pJacobian | Pointer 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 321 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References PetscMatTools::Finalise(), and PetscMatTools::SwitchWriteMode().
|
protected |
Computes the Jacobian numerically i.e. an approximation, using numerical partial derivatives.
currentGuess | Independent variable, u in f(u), for example |
pJacobian | A pointer to the Jacobian matrix |
Definition at line 342 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References PetscVecTools::AddToElement(), PetscTools::CreateVec(), PetscTools::Destroy(), PetscMatTools::Finalise(), CompareDoubles::IsNearZero(), PetscVecTools::Scale(), PetscMatTools::SetElement(), and PetscVecTools::WAXPY().
void AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::ComputeResidual | ( | const Vec | currentGuess, |
Vec | residualVector | ||
) |
Compute the residual vector given the current solution guess.
currentGuess | The solution guess for the current nonlinear solve iteration. |
residualVector | We fill this with the residual vector. |
NOTE: this method is called indirectly by the PETSc iterative solvers, so must be public.
Definition at line 300 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References PetscVecTools::Finalise().
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.
pNonlinearSolver | a nonlinear solver |
Definition at line 432 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
|
virtual |
Assemble and solve the system for a nonlinear elliptic PDE.
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. |
Definition at line 408 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References EXCEPTION.
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.
tol | A tolerance which defaults to 1e-5 |
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 443 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
References PetscMatTools::CheckEquality(), PetscTools::CreateAndSetVec(), PetscTools::Destroy(), and PetscTools::SetupMat().
|
protected |
Boundary conditions container.
Definition at line 140 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
|
protected |
The mesh.
Definition at line 137 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AbstractNonlinearAssemblerSolverHybrid().
|
protected |
An assembler for the surface integral terms in the residual
Definition at line 152 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AbstractNonlinearAssemblerSolverHybrid().
|
protected |
The nonlinear solver.
Definition at line 143 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AbstractNonlinearAssemblerSolverHybrid().
|
protected |
Whether to use an analytical expression for the Jacobian.
Definition at line 149 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
|
protected |
Whether memory has been allocated for the solver.
Definition at line 146 of file AbstractNonlinearAssemblerSolverHybrid.hpp.
Referenced by AbstractNonlinearAssemblerSolverHybrid< ELEMENT_DIM, SPACE_DIM, PROBLEM_DIM >::AbstractNonlinearAssemblerSolverHybrid().