Go to the documentation of this file.
20 #ifndef LIBMESH_DIFF_SYSTEM_H
21 #define LIBMESH_DIFF_SYSTEM_H
24 #include "libmesh/diff_context.h"
25 #include "libmesh/diff_physics.h"
26 #include "libmesh/diff_qoi.h"
27 #include "libmesh/implicit_system.h"
28 #include "libmesh/time_solver.h"
39 template <
typename T>
class NumericVector;
64 const std::string &
name,
65 const unsigned int number);
86 virtual void clear ()
override;
92 virtual void reinit ()
override;
117 virtual std::pair<unsigned int, Real>
136 virtual void assembly (
bool get_residual,
138 bool apply_heterogeneous_constraints =
false,
139 bool apply_no_constraints =
false)
override = 0;
146 virtual void solve ()
override;
152 virtual std::pair<unsigned int, Real>
160 libmesh_not_implemented();
162 return std::unique_ptr<DifferentiablePhysics>(
this);
168 virtual std::unique_ptr<DifferentiableQoI>
clone()
override
170 libmesh_not_implemented();
172 return std::unique_ptr<DifferentiableQoI>(
this);
399 #ifdef LIBMESH_ENABLE_DIRICHLET
419 libmesh_assert_equal_to (&(
time_solver->system()),
this);
427 libmesh_assert_equal_to (&(
time_solver->system()),
this);
434 #endif // LIBMESH_DIFF_SYSTEM_H
void swap_physics(DifferentiablePhysics *&swap_physics)
Swap current physics object with external object.
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time.
virtual void release_linear_solver(LinearSolver< Number > *) const override
Releases a pointer to a linear solver acquired by this->get_linear_solver()
bool print_jacobian_norms
Set print_jacobian_norms to true to print |J| whenever it is assembled.
bool print_element_jacobians
Set print_element_jacobians to true to print each J_elem contribution.
Manages consistently variables, degrees of freedom, coefficient vectors, and matrices for implicit sy...
bool print_residual_norms
Set print_residual_norms to true to print |F| whenever it is assembled.
virtual std::unique_ptr< DifferentiablePhysics > clone_physics() override
We don't allow systems to be attached to each other.
virtual void postprocess()
Executes a postprocessing loop over all elements, and if postprocess_sides is true over all sides.
DifferentiableSystem sys_type
The type of system.
bool print_element_solutions
Set print_element_solutions to true to print each U_elem input.
This class provides a specific system class.
const DifferentiableQoI * get_qoi() const
virtual std::unique_ptr< DifferentiableQoI > clone() override
We don't allow systems to be attached to each other.
bool print_residuals
Set print_residuals to true to print F whenever it is assembled.
virtual void reinit() override
Reinitializes the member data fields associated with the system, so that, e.g., assemble() may be use...
The libMesh namespace provides an interface to certain functionality in the library.
DifferentiableSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
void attach_qoi(DifferentiableQoI *qoi_in)
Attach external QoI object.
bool print_solution_norms
Set print_residual_norms to true to print |U| whenever it is used in an assembly() call.
unsigned int number() const
This class provides a specific system class.
std::unique_ptr< TimeSolver > time_solver
A pointer to the solver object we're going to use.
DifferentiablePhysics * _diff_physics
Pointer to object to use for physics assembly evaluations.
This class provides a specific system class.
std::vector< Number > qoi
Values of the quantities of interest.
virtual void init_qoi(std::vector< Number > &)
Initialize system qoi.
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
virtual void element_postprocess(DiffContext &)
Does any work that needs to be done on elem in a postprocessing loop.
bool print_jacobians
Set print_jacobians to true to print J whenever it is assembled.
virtual void solve() override
Invokes the solver associated with the system.
TimeSolver & get_time_solver()
virtual void assemble() override
Prepares matrix and rhs for matrix assembly.
virtual ~DifferentiableSystem()
Destructor.
virtual LinearSolver< Number > * get_linear_solver() const override
virtual void clear() override
Clear all the data structures associated with the system.
virtual std::pair< unsigned int, Real > adjoint_solve(const QoISet &qoi_indices=QoISet()) override
This function sets the _is_adjoint boolean member of TimeSolver to true and then calls the adjoint_so...
void add_second_order_dot_vars()
Helper function to add "velocity" variables that are cousins to second order-in-time variables in the...
ImplicitSystem Parent
The type of the parent.
This class provides all data required for a physics package (e.g.
This is the EquationSystems class.
bool have_first_order_scalar_vars() const
Check for any first order vars that are also belong to FEFamily::SCALAR.
void attach_physics(DifferentiablePhysics *physics_in)
Attach external Physics object.
DifferentiableQoI * diff_qoi
Pointer to object to use for quantity of interest assembly evaluations.
const std::string & name() const
const DifferentiablePhysics * get_physics() const
DifferentiableQoI * get_qoi()
virtual void assembly(bool get_residual, bool get_jacobian, bool apply_heterogeneous_constraints=false, bool apply_no_constraints=false) override=0
Assembles a residual in rhs and/or a jacobian in matrix, as requested.
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
This is a generic class that defines a solver to handle time integration of DifferentiableSystems.
bool postprocess_sides
If postprocess_sides is true (it is false by default), the postprocessing loop will loop over all sid...
bool print_solutions
Set print_solutions to true to print U whenever it is used in an assembly() call.
unsigned int get_second_order_dot_var(unsigned int var) const
For a given second order (in time) variable var, this method will return the index to the correspondi...
virtual void side_postprocess(DiffContext &)
Does any work that needs to be done on side of elem in a postprocessing loop.
virtual std::pair< unsigned int, Real > get_linear_solve_parameters() const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void set_time_solver(std::unique_ptr< TimeSolver > _time_solver)
Sets the time_solver FIXME: This code is a little dangerous as it transfers ownership from the TimeSo...
virtual void init_physics(const System &sys)
Initialize any data structures associated with the physics.
bool print_element_residuals
Set print_element_residuals to true to print each R_elem contribution.
virtual std::unique_ptr< DifferentiablePhysics > clone_physics()=0
Copy of this object.
DifferentiablePhysics * get_physics()
virtual std::unique_ptr< DifferentiableQoI > clone()=0
Copy of this object.
virtual std::unique_ptr< DiffContext > build_context()
Builds a DiffContext object with enough information to do evaluations on each element.
void add_dot_var_dirichlet_bcs(unsigned int var_idx, unsigned int dot_var_idx)
Helper function to and Dirichlet boundary conditions to "dot" variable cousins of second order variab...
bool have_second_order_scalar_vars() const
Check for any second order vars that are also belong to FEFamily::SCALAR.