MAST
MAST::PseudoArclengthContinuationSolver Class Reference

The constraint equation is defined along the path $ s $ as

\[ g(X, p, ds) = (Y - \tilde{Y} )^T t_1 = 0 , \]

where, $ Y = \{ X^T p \}^T $, $ X $ is the solution, $ p $ is the load parameter, $\tilde{Y} $ is the predictor based on the search direction $ t_1 $. More...

#include <pseudo_arclength_continuation_solver.h>

Inheritance diagram for MAST::PseudoArclengthContinuationSolver:
Collaboration diagram for MAST::PseudoArclengthContinuationSolver:

Public Member Functions

 PseudoArclengthContinuationSolver ()
 
virtual ~PseudoArclengthContinuationSolver ()
 
virtual void initialize (Real dp)
 initializes the search direction using the specified load step dp. More...
 
- Public Member Functions inherited from MAST::ContinuationSolverBase
 ContinuationSolverBase ()
 
virtual ~ContinuationSolverBase ()
 
void set_assembly_and_load_parameter (MAST::AssemblyElemOperations &elem_ops, MAST::AssemblyBase &assembly, MAST::Parameter &p)
 sets the assembly object for this solver More...
 
void clear_assembly_and_load_parameters ()
 clears the assembly object from this solver More...
 
virtual void solve ()
 solves for the next load step More...
 

Protected Member Functions

virtual void _solve_NR_iterate (libMesh::NumericVector< Real > &X, MAST::Parameter &p)
 
void _update_search_direction (const libMesh::NumericVector< Real > &X, const MAST::Parameter &p, libMesh::SparseMatrix< Real > &jac, libMesh::NumericVector< Real > &t1_X, Real &t1_p)
 updates $ t_1 $ for the current iterate X and p, and stores these values to _t0_X and _t0_p for next computation. More...
 
virtual Real _g (const libMesh::NumericVector< Real > &X, const MAST::Parameter &p)
 

\[ g(X, p, ds) = X_{scale} * (X-X0) (dX/ds)_{scaled} + p_{scale} * (p-p0) (dp/ds)_{scaled} - ds = 0, \]

where, $ t_1^{X} = (dX/ds)_{scaled} $ and $ t_1^{p} = (dp/ds)_{scaled} $. More...

 
void _g (const libMesh::NumericVector< Real > &X, const MAST::Parameter &p, libMesh::NumericVector< Real > &t1_X, Real &t1_p, Real &g)
 

\begin{eqnarray*} g(X, p, ds) & = & (Y-Y_0)^T t_1 - ds = 0 \\ Y & = & \{ X_{scale} X^T, p_{scale} p\}^T \\ dg/dp & = & p_{scale} t_1^p \\ dg/dX & = & X_{scale} t_1^X, \end{eqnarray*}

where $ $, $ t_1^{X} = (dX/ds)_{scaled} $ and $ t_1^{p} = (dp/ds)_{scaled} $. More...

 
virtual void _save_iteration_data ()
 method saves any data for possible resuse if the solution step is restarted More...
 
virtual void _reset_iterations ()
 method resets any data if a solution step is restarted More...
 
- Protected Member Functions inherited from MAST::ContinuationSolverBase
void _solve (const libMesh::NumericVector< Real > &X, const MAST::Parameter &p, libMesh::NumericVector< Real > &f, bool update_f, libMesh::NumericVector< Real > &dfdp, bool update_dfdp, const libMesh::NumericVector< Real > &dgdX, const Real dgdp, const Real g, libMesh::NumericVector< Real > &dX, Real &dp)
 solves for the linear system of equation as a monolithic system

\[ \left[ \begin{array}{cc} df/dx & df/dp \\ dg/dx & dg/dp \end{array}\right] \left\{ \begin{array}{c} dx \\ dp \end{array} \right\} = - \left\{ \begin{array}{c} f \\ g \end{array} \right\} \]

dX and dp are returned from the solution More...

 
void _solve_schur_factorization (const libMesh::NumericVector< Real > &X, const MAST::Parameter &p, libMesh::SparseMatrix< Real > &jac, bool update_jac, libMesh::NumericVector< Real > &f, bool update_f, libMesh::NumericVector< Real > &dfdp, bool update_dfdp, libMesh::NumericVector< Real > &dXdp, bool update_dXdp, const libMesh::NumericVector< Real > &dgdX, const Real dgdp, const Real g, libMesh::NumericVector< Real > &dX, Real &dp)
 solves for the linear system of equation using Schur factorization. More...
 
Real _res_norm (const libMesh::NumericVector< Real > &X, const MAST::Parameter &p)
 

Protected Attributes

std::unique_ptr< libMesh::NumericVector< Real > > _t0_X
 
std::unique_ptr< libMesh::NumericVector< Real > > _t0_X_orig
 
Real _t0_p
 
Real _t0_p_orig
 
- Protected Attributes inherited from MAST::ContinuationSolverBase
bool _initialized
 
MAST::AssemblyElemOperations_elem_ops
 
MAST::AssemblyBase_assembly
 
MAST::Parameter_p
 
Real _p0
 
Real _X_scale
 
Real _p_scale
 
std::unique_ptr< libMesh::NumericVector< Real > > _X0
 

Additional Inherited Members

- Public Attributes inherited from MAST::ContinuationSolverBase
unsigned int max_it
 Maximum number of Newton-Raphson iterations for the solver. More...
 
Real abs_tol
 Absolute tolerance for the solver. More...
 
Real rel_tol
 Relative tolerance for the solver. More...
 
Real arc_length
 arc length that the solver is required to satisfy for the update. More...
 
Real min_step
 minimum step size allowed with adaptivity More...
 
Real max_step
 maximum step size allowed with adaptivity More...
 
Real step_size_change_exponent
 exponent used in step size update. More...
 
unsigned int step_desired_iters
 desired N-R iterations per load-step. More...
 
bool schur_factorization
 flag to use Schur-factorizaiton (default) or monolithic solver More...
 

Detailed Description

The constraint equation is defined along the path $ s $ as

\[ g(X, p, ds) = (Y - \tilde{Y} )^T t_1 = 0 , \]

where, $ Y = \{ X^T p \}^T $, $ X $ is the solution, $ p $ is the load parameter, $\tilde{Y} $ is the predictor based on the search direction $ t_1 $.

Given that the predictor is defined as $ \tilde{Y} = Y_0 - ds t_1 $, with $ ds $ as the step size, the constraint is rewritten as

\[ (Y - Y_0 - ds * t_1)^T t_1 = 0 , \]

or, assuming $ \| t_1 \| = 1 $,

\[ Y^T t_1 - Y_0^T t_1 - ds = 0 . \]

The search direction is evaluated based on:

\[ \left[ \begin{array}{cc} df/dx & df/dp \\ t_0^X & t_0^p \end{array} \right] \left\{ \begin{array}{c} t_1^X \\ t_1^p \end{array} \right\} = \left\{ \begin{array}{c} 0 \\ 1 \end{array} \right\} \]

and, $ t_1 $ is then scaled to unit norm.

Definition at line 49 of file pseudo_arclength_continuation_solver.h.

Constructor & Destructor Documentation

MAST::PseudoArclengthContinuationSolver::PseudoArclengthContinuationSolver ( )

Definition at line 31 of file pseudo_arclength_continuation_solver.cpp.

MAST::PseudoArclengthContinuationSolver::~PseudoArclengthContinuationSolver ( )
virtual

Definition at line 39 of file pseudo_arclength_continuation_solver.cpp.

Member Function Documentation

Real MAST::PseudoArclengthContinuationSolver::_g ( const libMesh::NumericVector< Real > &  X,
const MAST::Parameter p 
)
protectedvirtual

\[ g(X, p, ds) = X_{scale} * (X-X0) (dX/ds)_{scaled} + p_{scale} * (p-p0) (dp/ds)_{scaled} - ds = 0, \]

where, $ t_1^{X} = (dX/ds)_{scaled} $ and $ t_1^{p} = (dp/ds)_{scaled} $.

Implements MAST::ContinuationSolverBase.

Definition at line 205 of file pseudo_arclength_continuation_solver.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::PseudoArclengthContinuationSolver::_g ( const libMesh::NumericVector< Real > &  X,
const MAST::Parameter p,
libMesh::NumericVector< Real > &  t1_X,
Real t1_p,
Real g 
)
protected

\begin{eqnarray*} g(X, p, ds) & = & (Y-Y_0)^T t_1 - ds = 0 \\ Y & = & \{ X_{scale} X^T, p_{scale} p\}^T \\ dg/dp & = & p_{scale} t_1^p \\ dg/dX & = & X_{scale} t_1^X, \end{eqnarray*}

where $ $, $ t_1^{X} = (dX/ds)_{scaled} $ and $ t_1^{p} = (dp/ds)_{scaled} $.

Definition at line 228 of file pseudo_arclength_continuation_solver.cpp.

void MAST::PseudoArclengthContinuationSolver::_reset_iterations ( )
protectedvirtual

method resets any data if a solution step is restarted

Implements MAST::ContinuationSolverBase.

Definition at line 259 of file pseudo_arclength_continuation_solver.cpp.

void MAST::PseudoArclengthContinuationSolver::_save_iteration_data ( )
protectedvirtual

method saves any data for possible resuse if the solution step is restarted

Implements MAST::ContinuationSolverBase.

Definition at line 248 of file pseudo_arclength_continuation_solver.cpp.

Here is the caller graph for this function:

void MAST::PseudoArclengthContinuationSolver::_solve_NR_iterate ( libMesh::NumericVector< Real > &  X,
MAST::Parameter p 
)
protectedvirtual

Implements MAST::ContinuationSolverBase.

Definition at line 86 of file pseudo_arclength_continuation_solver.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::PseudoArclengthContinuationSolver::_update_search_direction ( const libMesh::NumericVector< Real > &  X,
const MAST::Parameter p,
libMesh::SparseMatrix< Real > &  jac,
libMesh::NumericVector< Real > &  t1_X,
Real t1_p 
)
protected

updates $ t_1 $ for the current iterate X and p, and stores these values to _t0_X and _t0_p for next computation.

Definition at line 144 of file pseudo_arclength_continuation_solver.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::PseudoArclengthContinuationSolver::initialize ( Real  dp)
virtual

initializes the search direction using the specified load step dp.

Implements MAST::ContinuationSolverBase.

Definition at line 45 of file pseudo_arclength_continuation_solver.cpp.

Here is the call graph for this function:

Member Data Documentation

Real MAST::PseudoArclengthContinuationSolver::_t0_p
protected

Definition at line 128 of file pseudo_arclength_continuation_solver.h.

Real MAST::PseudoArclengthContinuationSolver::_t0_p_orig
protected

Definition at line 128 of file pseudo_arclength_continuation_solver.h.

std::unique_ptr<libMesh::NumericVector<Real> > MAST::PseudoArclengthContinuationSolver::_t0_X
protected

Definition at line 127 of file pseudo_arclength_continuation_solver.h.

std::unique_ptr<libMesh::NumericVector<Real> > MAST::PseudoArclengthContinuationSolver::_t0_X_orig
protected

Definition at line 127 of file pseudo_arclength_continuation_solver.h.


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