Go to the documentation of this file.
20 #ifndef LIBMESH_CONTINUATION_SYSTEM_H
21 #define LIBMESH_CONTINUATION_SYSTEM_H
24 #include "libmesh/fem_system.h"
32 template <
typename T>
class LinearSolver;
63 const std::string &
name,
64 const unsigned int number);
85 virtual void clear ()
override;
90 virtual void solve ()
override;
429 #endif // LIBMESH_CONTINUATION_SYSTEM_H
Real Theta_LOCA
Another normalization parameter, which is described in the LOCA manual.
virtual ~ContinuationSystem()
Destructor.
NumericVector< Number > * du_ds
Extra work vectors used by the continuation algorithm.
Real ds_current
Value of stepsize currently in use.
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
bool tangent_initialized
False until initialize_tangent() is called.
RHS_Mode
There are (so far) two different vectors which may be assembled using the assembly routine: 1....
void set_max_arclength_stepsize(Real maxds)
Sets (initializes) the max-allowable ds value and the current ds value.
Real ds_min
The minimum-allowed steplength, defaults to 1.e-8.
virtual void solve() override
Perform a standard "solve" of the system, without doing continuation.
The libMesh namespace provides an interface to certain functionality in the library.
NumericVector< Number > * z
Temporary vector "z" ...
unsigned int number() const
Real ds
The initial arclength stepsize, selected by the user.
NumericVector< Number > * y
Temporary vector "y" ...
This class provides a specific system class.
Real min_continuation_parameter
The minimum-allowable value of the continuation parameter.
void continuation_solve()
Perform a continuation solve of the system.
Real * continuation_parameter
The continuation parameter must be a member variable of the derived class, and the "continuation_para...
First-order Euler predictor.
std::unique_ptr< LinearSolver< Number > > linear_solver
We maintain our own linear solver interface, for solving custom systems of equations and/or things wh...
Real previous_ds
The previous arcstep length used.
ContinuationSystem sys_type
The type of system.
Real previous_dlambda_ds
The old parameter tangent value.
NumericVector< Number > * delta_u
Temporary vector "delta u" ...
void solve_tangent()
Special solve algorithm for solving the tangent system.
NumericVector< Number > * y_old
Temporary vector "y_old" ...
Real max_continuation_parameter
The maximum-allowable value of the continuation parameter.
NumericVector< Number > * previous_du_ds
The value of du_ds from the previous solution.
void set_Theta_LOCA()
A centralized function for setting the other normalization parameter, i.e.
bool newton_progress_check
True by default, the Newton progress check allows the Newton loop to exit if half the allowed iterati...
unsigned int n_backtrack_steps
Another scaling parameter suggested by the LOCA people.
void apply_predictor()
Applies the predictor to the current solution to get a guess for the next solution.
unsigned int newton_step
Loop counter for nonlinear (Newton) iteration loop.
double solution_tolerance
How tightly should the Newton iterations attempt to converge ||delta_u|| Defaults to 1....
This class defines a solver which uses the default libMesh linear solver in a quasiNewton method to h...
Real old_continuation_parameter
The system also keeps track of the old value of the continuation parameter.
void save_current_solution()
Stores the current solution and continuation parameter (as "previous_u" and "old_continuation_paramet...
virtual void clear() override
Clear all the data structures associated with the system.
This is the EquationSystems class.
FEMSystem Parent
The type of the parent.
bool quiet
If quiet==false, the System prints extra information about what it is doing.
const std::string & name() const
Real newton_stepgrowth_aggressiveness
A measure of how rapidly one should attempt to grow the arclength stepsize based on the number of New...
void advance_arcstep()
Call this function after a continuation solve to compute the tangent and get the next initial guess.
unsigned int n_arclength_reductions
Number of arclength reductions to try when Newton fails to reduce the residual.
Predictor
The code provides the ability to select from different predictor schemes for getting the initial gues...
Second-order explicit Adams-Bashforth predictor.
void initialize_tangent()
Before starting arclength continuation, we need at least 2 prior solutions (both solution and u_previ...
ContinuationSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Real Theta
Arclength normalization parameter.
double initial_newton_tolerance
How much to try to reduce the residual by at the first (inexact) Newton step.
void set_Theta()
A centralized function for setting the normalization parameter theta.
double continuation_parameter_tolerance
How tightly should the Newton iterations attempt to converge delta_lambda.
void update_solution()
This function (which assumes the most recent tangent vector has been computed) updates the solution a...
This class inherits from the FEMSystem.
NumericVector< Number > * previous_u
The solution at the previous value of the continuation variable.
Real dlambda_ds
The most recent value of the derivative of the continuation parameter with respect to s.
NewtonSolver * newton_solver
A pointer to the underlying Newton solver used by the DiffSystem.