20 #ifndef LIBMESH_DIFF_SYSTEM_H 21 #define LIBMESH_DIFF_SYSTEM_H 40 template <
typename T>
class NumericVector;
65 const std::string &
name,
66 const unsigned int number);
87 virtual void clear ()
override;
93 virtual void reinit ()
override;
118 virtual std::pair<unsigned int, Real>
137 virtual void assembly (
bool get_residual,
139 bool apply_heterogeneous_constraints =
false,
140 bool apply_no_constraints =
false)
override = 0;
147 virtual void solve ()
override;
153 virtual std::pair<unsigned int, Real>
161 libmesh_not_implemented();
163 return std::unique_ptr<DifferentiablePhysics>(
this);
169 virtual std::unique_ptr<DifferentiableQoI>
clone()
override 171 libmesh_not_implemented();
173 return std::unique_ptr<DifferentiableQoI>(
this);
418 libmesh_assert_equal_to (&(
time_solver->system()),
this);
426 libmesh_assert_equal_to (&(
time_solver->system()),
this);
433 #endif // LIBMESH_DIFF_SYSTEM_H void attach_physics(DifferentiablePhysics *physics_in)
const DifferentiableQoI * get_qoi() const
Manages multiples systems of equations.
virtual void assemble() override
virtual void assembly(bool get_residual, bool get_jacobian, bool apply_heterogeneous_constraints=false, bool apply_no_constraints=false) override=0
Used to specify quantities of interest in a simulation.
bool print_jacobian_norms
std::unique_ptr< TimeSolver > time_solver
DifferentiableQoI * get_qoi()
void add_second_order_dot_vars()
DifferentiablePhysics * get_physics()
void add_dot_var_dirichlet_bcs(unsigned int var_idx, unsigned int dot_var_idx)
virtual std::unique_ptr< DiffContext > build_context()
bool print_element_residuals
DifferentiableQoI * diff_qoi
virtual void reinit() override
bool print_element_solutions
bool print_solution_norms
void swap_physics(DifferentiablePhysics *&swap_physics)
DifferentiablePhysics * _diff_physics
unsigned int number() const
unsigned int get_second_order_dot_var(unsigned int var) const
bool print_residual_norms
virtual void side_postprocess(DiffContext &)
std::vector< Number > qoi
virtual std::unique_ptr< DifferentiablePhysics > clone_physics() override
DifferentiableSystem sys_type
void attach_qoi(DifferentiableQoI *qoi_in)
virtual void element_postprocess(DiffContext &)
bool have_second_order_scalar_vars() const
virtual std::unique_ptr< DifferentiableQoI > clone()=0
const DifferentiablePhysics * get_physics() const
DifferentiableSystem(EquationSystems &es, const std::string &name, const unsigned int number)
virtual void init_data() override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual LinearSolver< Number > * get_linear_solver() const override
virtual std::pair< unsigned int, Real > get_linear_solve_parameters() const override
virtual std::unique_ptr< DifferentiableQoI > clone() override
bool print_element_jacobians
void set_time_solver(std::unique_ptr< TimeSolver > _time_solver)
virtual std::unique_ptr< DifferentiablePhysics > clone_physics()=0
virtual void init_physics(const System &sys)
virtual void solve() override
virtual void clear() override
const std::string & name() const
virtual void release_linear_solver(LinearSolver< Number > *) const override
bool have_first_order_scalar_vars() const
virtual void postprocess()
virtual std::pair< unsigned int, Real > adjoint_solve(const QoISet &qoi_indices=QoISet()) override
virtual ~DifferentiableSystem()
TimeSolver & get_time_solver()
virtual void init_qoi(std::vector< Number > &)
Manages consistently variables, degrees of freedom, coefficient vectors, and matrices for implicit sy...