Loading [MathJax]/extensions/tex2jax.js
libMesh
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
diff_system.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_DIFF_SYSTEM_H
21 #define LIBMESH_DIFF_SYSTEM_H
22 
23 // Local Includes
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"
29 
30 // C++ includes
31 #include <memory>
32 #include <stack>
33 
34 namespace libMesh
35 {
36 
37 // Forward Declarations
38 class TimeSolver;
39 
40 template <typename T> class NumericVector;
41 
55  public virtual DifferentiablePhysics,
56  public virtual DifferentiableQoI
57 {
58 public:
59 
64  const std::string & name,
65  const unsigned int number);
66 
72  DifferentiableSystem (const DifferentiableSystem &) = delete;
76  virtual ~DifferentiableSystem ();
77 
82 
87 
92  virtual void clear () override;
93 
98  virtual void reinit () override;
99 
110  virtual void assemble () override;
111 
116  virtual LinearSolver<Number> * get_linear_solver() const override;
117 
123  virtual std::pair<unsigned int, Real>
124  get_linear_solve_parameters() const override;
125 
136  virtual void assembly (bool get_residual,
137  bool get_jacobian,
138  bool apply_heterogeneous_constraints = false,
139  bool apply_no_constraints = false) override = 0;
140 
146  virtual void solve () override;
147 
152  virtual std::pair<unsigned int, Real>
153  adjoint_solve (const QoISet & qoi_indices = QoISet()) override;
154 
158  virtual std::unique_ptr<DifferentiablePhysics> clone_physics() override
159  {
160  libmesh_not_implemented();
161  // dummy to avoid compiler warnings, not a real implementation
162  return std::unique_ptr<DifferentiablePhysics>(nullptr);
163  }
164 
168  virtual std::unique_ptr<DifferentiableQoI> clone() override
169  {
170  libmesh_not_implemented();
171  // dummy to avoid compiler warnings, not a real implementation
172  return std::unique_ptr<DifferentiableQoI>(nullptr);
173  }
174 
182  { if (this->_diff_physics.empty())
183  return this;
184  return this->_diff_physics.top().get(); }
185 
193  { if (this->_diff_physics.empty())
194  return this;
195  return this->_diff_physics.top().get(); }
196 
201  { this->_diff_physics.push(physics_in->clone_physics());
202  this->_diff_physics.top()->init_physics(*this);}
203 
204 #ifdef LIBMESH_ENABLE_DEPRECATED
205 
210 #endif // LIBMESH_ENABLE_DEPRECATED
211 
217  void push_physics ( DifferentiablePhysics & new_physics );
218 
222  void pop_physics ();
223 
229  const DifferentiableQoI * get_qoi() const
230  { if (this->_diff_qoi.empty())
231  return this;
232  return this->_diff_qoi.top().get(); }
233 
240  { if (this->_diff_qoi.empty())
241  return this;
242  return this->_diff_qoi.top().get(); }
243 
247  void attach_qoi( DifferentiableQoI * qoi_in );
248 
253  std::unique_ptr<TimeSolver> time_solver;
254 
261  void set_time_solver(std::unique_ptr<TimeSolver> _time_solver)
262  {
263  time_solver.reset(_time_solver.release());
264  }
265 
270 
274  const TimeSolver & get_time_solver() const;
275 
281 
290  virtual std::unique_ptr<DiffContext> build_context();
291 
296  virtual void postprocess () {}
297 
301  virtual void element_postprocess (DiffContext &) {}
302 
307  virtual void side_postprocess (DiffContext &) {}
308 
318  unsigned int get_second_order_dot_var( unsigned int var ) const;
319 
323  bool have_first_order_scalar_vars() const;
324 
328  bool have_second_order_scalar_vars() const;
329 
335 
341 
347 
352 
357 
362 
367 
372 
377 
382 
387  virtual void set_constrain_in_solver(bool enable);
388 
389  virtual bool get_constrain_in_solver()
390  {
391  return _constrain_in_solver;
392  }
393 
394 protected:
395 
407  virtual void init_data () override;
408 
415 
416 #ifdef LIBMESH_ENABLE_DIRICHLET
417 
425  void add_dot_var_dirichlet_bcs( unsigned int var_idx, unsigned int dot_var_idx);
426 #endif
427 
433 
434 private:
441  std::stack<std::unique_ptr<DifferentiablePhysics>,
442  std::vector<std::unique_ptr<DifferentiablePhysics>>> _diff_physics;
443 
449  std::stack<std::unique_ptr<DifferentiableQoI>,
450  std::vector<std::unique_ptr<DifferentiableQoI>>> _diff_qoi;
451 };
452 
453 // --------------------------------------------------------------
454 // DifferentiableSystem inline methods
455 inline
457 {
459  libmesh_assert_equal_to (&(time_solver->system()), this);
460  return *time_solver;
461 }
462 
463 inline
465 {
467  libmesh_assert_equal_to (&(time_solver->system()), this);
468  return *time_solver;
469 }
470 
471 } // namespace libMesh
472 
473 
474 #endif // LIBMESH_DIFF_SYSTEM_H
void attach_physics(DifferentiablePhysics *physics_in)
Attach external Physics object.
Definition: diff_system.h:200
const DifferentiableQoI * get_qoi() const
Definition: diff_system.h:229
This is the EquationSystems class.
virtual void assemble() override
Prepares matrix and rhs for matrix assembly.
Definition: diff_system.C:131
bool _constrain_in_solver
_constrain_in_solver defaults to true; if false then we apply constraints only via residual terms in ...
Definition: diff_system.h:432
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
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.
This is a generic class that defines a solver to handle time integration of DifferentiableSystems.
Definition: time_solver.h:63
virtual void set_constrain_in_solver(bool enable)
set_constrain_in_solver to false to apply constraints only via residual terms in the systems to be so...
Definition: diff_system.C:424
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
bool print_jacobian_norms
Set print_jacobian_norms to true to print |J| whenever it is assembled.
Definition: diff_system.h:361
std::unique_ptr< TimeSolver > time_solver
A pointer to the solver object we&#39;re going to use.
Definition: diff_system.h:253
DifferentiableQoI * get_qoi()
Definition: diff_system.h:239
void add_second_order_dot_vars()
Helper function to add "velocity" variables that are cousins to second order-in-time variables in the...
Definition: diff_system.C:182
bool postprocess_sides
If postprocess_sides is true (it is false by default), the postprocessing loop will loop over all sid...
Definition: diff_system.h:334
DifferentiablePhysics * get_physics()
Definition: diff_system.h:192
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...
Definition: diff_system.C:214
The libMesh namespace provides an interface to certain functionality in the library.
virtual std::unique_ptr< DiffContext > build_context()
Builds a DiffContext object with enough information to do evaluations on each element.
Definition: diff_system.C:123
bool print_jacobians
Set print_jacobians to true to print J whenever it is assembled.
Definition: diff_system.h:366
void pop_physics()
Pop a physics object off of our stack.
Definition: diff_system.C:411
bool print_element_residuals
Set print_element_residuals to true to print each R_elem contribution.
Definition: diff_system.h:376
virtual void reinit() override
Reinitializes the member data fields associated with the system, so that, e.g., assemble() may be use...
Definition: diff_system.C:82
bool print_element_solutions
Set print_element_solutions to true to print each U_elem input.
Definition: diff_system.h:371
This class provides a specific system class.
Definition: diff_system.h:54
virtual bool get_constrain_in_solver()
Definition: diff_system.h:389
bool print_solution_norms
Set print_residual_norms to true to print |U| whenever it is used in an assembly() call...
Definition: diff_system.h:340
void swap_physics(DifferentiablePhysics *&swap_physics)
Swap current physics object with external object.
Definition: diff_system.C:350
unsigned int number() const
Definition: system.h:2337
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...
Definition: diff_system.C:306
bool print_residual_norms
Set print_residual_norms to true to print |F| whenever it is assembled.
Definition: diff_system.h:351
virtual void side_postprocess(DiffContext &)
Does any work that needs to be done on side of elem in a postprocessing loop.
Definition: diff_system.h:307
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time...
Definition: diff_system.h:280
std::stack< std::unique_ptr< DifferentiablePhysics >, std::vector< std::unique_ptr< DifferentiablePhysics > > > _diff_physics
Stack of pointers to objects to use for physics assembly evaluations.
Definition: diff_system.h:442
std::stack< std::unique_ptr< DifferentiableQoI >, std::vector< std::unique_ptr< DifferentiableQoI > > > _diff_qoi
Pointer to object to use for quantity of interest assembly evaluations.
Definition: diff_system.h:450
bool print_residuals
Set print_residuals to true to print F whenever it is assembled.
Definition: diff_system.h:356
libmesh_assert(ctx)
virtual std::unique_ptr< DifferentiablePhysics > clone_physics() override
We don&#39;t allow systems to be attached to each other.
Definition: diff_system.h:158
DifferentiableSystem sys_type
The type of system.
Definition: diff_system.h:81
void attach_qoi(DifferentiableQoI *qoi_in)
Attach external QoI object.
Definition: diff_system.C:279
virtual void element_postprocess(DiffContext &)
Does any work that needs to be done on elem in a postprocessing loop.
Definition: diff_system.h:301
bool have_second_order_scalar_vars() const
Check for any second order vars that are also belong to FEFamily::SCALAR.
Definition: diff_system.C:335
bool print_solutions
Set print_solutions to true to print U whenever it is used in an assembly() call. ...
Definition: diff_system.h:346
const DifferentiablePhysics * get_physics() const
Definition: diff_system.h:181
This class provides a specific system class.
Definition: diff_physics.h:76
DifferentiableSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
Definition: diff_system.C:36
This class provides a specific system class.
Definition: diff_qoi.h:52
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used...
Definition: diff_system.C:94
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual LinearSolver< Number > * get_linear_solver() const override
Definition: diff_system.C:163
virtual std::pair< unsigned int, Real > get_linear_solve_parameters() const override
Definition: diff_system.C:172
ImplicitSystem Parent
The type of the parent.
Definition: diff_system.h:86
virtual std::unique_ptr< DifferentiableQoI > clone() override
We don&#39;t allow systems to be attached to each other.
Definition: diff_system.h:168
bool print_element_jacobians
Set print_element_jacobians to true to print each J_elem contribution.
Definition: diff_system.h:381
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...
Definition: diff_system.h:261
virtual std::unique_ptr< DifferentiablePhysics > clone_physics()=0
Copy of this object.
void push_physics(DifferentiablePhysics &new_physics)
Push a clone of a new physics object onto our stack, overriding the current physics until the new phy...
Definition: diff_system.C:399
virtual void solve() override
Invokes the solver associated with the system.
Definition: diff_system.C:138
virtual void clear() override
Clear all the data structures associated with the system.
Definition: diff_system.C:64
const std::string & name() const
Definition: system.h:2329
bool have_first_order_scalar_vars() const
Check for any first order vars that are also belong to FEFamily::SCALAR.
Definition: diff_system.C:323
virtual void postprocess()
Executes a postprocessing loop over all elements, and if postprocess_sides is true over all sides...
Definition: diff_system.h:296
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...
Definition: diff_system.C:150
DifferentiableSystem & operator=(const DifferentiableSystem &)=delete
TimeSolver & get_time_solver()
Definition: diff_system.h:456
Manages consistently variables, degrees of freedom, coefficient vectors, and matrices for implicit sy...