20 #ifndef LIBMESH_TIME_SOLVER_H 21 #define LIBMESH_TIME_SOLVER_H 24 #include "libmesh/libmesh_common.h" 25 #include "libmesh/reference_counted_object.h" 36 class DifferentiablePhysics;
37 class DifferentiableSystem;
38 class ParameterVector;
39 class SensitivityData;
40 class SolutionHistory;
43 class AdjointRefinementEstimator;
114 virtual void solve ();
154 #ifdef LIBMESH_ENABLE_AMR 163 #endif // LIBMESH_ENABLE_AMR 348 #endif // LIBMESH_TIME_SOLVER_H bool quiet
Print extra debugging information if quiet == false.
virtual bool side_residual(bool request_jacobian, DiffContext &)=0
This method uses the DifferentiablePhysics side_time_derivative(), side_constraint(), and side_mass_residual() to build a full residual on an element's side.
virtual Real last_completed_timestep_size()
Returns system.deltat if fixed timestep solver is used, the complete timestep size (sum of all subste...
This class provides all data required for a physics package (e.g.
std::unique_ptr< DiffSolver > _diff_solver
An implicit linear or nonlinear solver to use at each timestep.
virtual void advance_timestep()
This method advances the solution to the next timestep, after a solve() has been performed.
Data structure for specifying which Parameters should be independent variables in a parameter sensiti...
This is a generic class that defines a solver to handle time integration of DifferentiableSystems.
This class implements a "brute force" goal-oriented error estimator which computes an estimate of err...
virtual void init_data()
The data initialization function.
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
std::unique_ptr< SolutionHistory > solution_history
A std::unique_ptr to a SolutionHistory object.
virtual void reinit()
The reinitialization function.
virtual void integrate_adjoint_refinement_error_estimate(AdjointRefinementEstimator &adjoint_refinement_error_estimator, ErrorVector &QoI_elementwise_error)
A method to compute the adjoint refinement error estimate at the current timestep.
This class defines a norm/seminorm to be applied to a NumericVector which contains coefficients in a ...
virtual bool is_steady() const =0
Is this effectively a steady-state solver?
The libMesh namespace provides an interface to certain functionality in the library.
void(DiffContext::* ReinitFuncType)(Real)
virtual bool nonlocal_residual(bool request_jacobian, DiffContext &)=0
This method uses the DifferentiablePhysics nonlocal_time_derivative(), nonlocal_constraint(), and nonlocal_mass_residual() to build a full residual of non-local terms.
virtual std::unique_ptr< DiffSolver > & diff_solver()
An implicit linear or nonlinear solver to use at each timestep.
A SolutionHistory class that enables the storage and retrieval of timesteps and (in the future) adapt...
const sys_type & system() const
This class provides a specific system class.
virtual bool element_residual(bool request_jacobian, DiffContext &)=0
This method uses the DifferentiablePhysics element_time_derivative(), element_constraint(), and mass_residual() to build a full residual on an element.
virtual void integrate_adjoint_sensitivity(const QoISet &qois, const ParameterVector ¶meter_vector, SensitivityData &sensitivities)
A method to integrate the adjoint sensitivity w.r.t a given parameter vector.
sys_type & _system
A reference to the system we are solving.
void set_is_adjoint(bool _is_adjoint_value)
Accessor for setting whether we need to do a primal or adjoint solve.
TimeSolver(sys_type &s)
Constructor.
Data structure for holding completed parameter sensitivity calculations.
virtual void init()
The initialization function.
bool is_adjoint() const
Accessor for querying whether we need to do a primal or adjoint solve.
bool _is_adjoint
This boolean tells the TimeSolver whether we are solving a primal or adjoint problem.
virtual std::pair< unsigned int, Real > adjoint_solve(const QoISet &qoi_indices)
This method solves for the adjoint solution at the next adjoint timestep (or a steady state adjoint s...
bool(DifferentiablePhysics::* ResFuncType)(bool, DiffContext &)
Definitions of argument types for use in refactoring subclasses.
unsigned int reduce_deltat_on_diffsolver_failure
This value (which defaults to zero) is the number of times the TimeSolver is allowed to halve deltat ...
This class implements reference counting.
virtual void retrieve_timestep()
This method retrieves all the stored solutions at the current system.time.
std::unique_ptr< LinearSolver< Number > > _linear_solver
An implicit linear solver to use for adjoint problems.
virtual void before_timestep()
This method is for subclasses or users to override to do arbitrary processing between timesteps...
This class provides a specific system class.
DifferentiableSystem sys_type
The type of system.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void set_solution_history(const SolutionHistory &_solution_history)
A setter function users will employ if they need to do something other than save no solution history...
virtual std::unique_ptr< LinearSolver< Number > > & linear_solver()
An implicit linear solver to use for adjoint and sensitivity problems.
virtual ~TimeSolver()
Destructor.
virtual void integrate_qoi_timestep()
A method to integrate the system::QoI functionals.
Real last_deltat
The deltat for the last completed timestep before the current one.
virtual void adjoint_advance_timestep()
This method advances the adjoint solution to the previous timestep, after an adjoint_solve() has been...
virtual Real du(const SystemNorm &norm) const =0
Computes the size of ||u^{n+1} - u^{n}|| in some norm.
virtual void solve()
This method solves for the solution at the next timestep (or solves for a steady-state solution)...
virtual void init_adjoints()
Initialize any adjoint related data structures, based on the number of qois.
SolutionHistory & get_solution_history()
A getter function that returns a reference to the solution history object owned by TimeSolver...