libMesh
diff_system.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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 
33 namespace libMesh
34 {
35 
36 // Forward Declarations
37 class TimeSolver;
38 
39 template <typename T> class NumericVector;
40 
54  public virtual DifferentiablePhysics,
55  public virtual DifferentiableQoI
56 {
57 public:
58 
64  const std::string & name,
65  const unsigned int number);
66 
70  virtual ~DifferentiableSystem ();
71 
76 
81 
86  virtual void clear () override;
87 
92  virtual void reinit () override;
93 
104  virtual void assemble () override;
105 
110  virtual LinearSolver<Number> * get_linear_solver() const override;
111 
117  virtual std::pair<unsigned int, Real>
118  get_linear_solve_parameters() const override;
119 
124  virtual void release_linear_solver(LinearSolver<Number> *) 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
162  return std::unique_ptr<DifferentiablePhysics>(this);
163  }
164 
168  virtual std::unique_ptr<DifferentiableQoI> clone() override
169  {
170  libmesh_not_implemented();
171  // dummy
172  return std::unique_ptr<DifferentiableQoI>(this);
173  }
174 
182  { return this->_diff_physics; }
183 
191  { return this->_diff_physics; }
192 
197  { this->_diff_physics = (physics_in->clone_physics()).release();
198  this->_diff_physics->init_physics(*this);}
199 
204 
210  const DifferentiableQoI * get_qoi() const
211  { return this->diff_qoi; }
212 
219  { return this->diff_qoi; }
220 
224  void attach_qoi( DifferentiableQoI * qoi_in )
225  { this->diff_qoi = (qoi_in->clone()).release();
226  // User needs to resize qoi system qoi accordingly
227  this->diff_qoi->init_qoi( this->qoi );}
228 
233  std::unique_ptr<TimeSolver> time_solver;
234 
241  void set_time_solver(std::unique_ptr<TimeSolver> _time_solver)
242  {
243  time_solver.reset(_time_solver.release());
244  }
245 
250 
254  const TimeSolver & get_time_solver() const;
255 
261 
270  virtual std::unique_ptr<DiffContext> build_context();
271 
276  virtual void postprocess () {}
277 
281  virtual void element_postprocess (DiffContext &) {}
282 
287  virtual void side_postprocess (DiffContext &) {}
288 
298  unsigned int get_second_order_dot_var( unsigned int var ) const;
299 
303  bool have_first_order_scalar_vars() const;
304 
308  bool have_second_order_scalar_vars() const;
309 
315 
321 
327 
332 
337 
342 
347 
352 
357 
362 
363 protected:
364 
371 
378 
390  virtual void init_data () override;
391 
398 
399 #ifdef LIBMESH_ENABLE_DIRICHLET
400 
408  void add_dot_var_dirichlet_bcs( unsigned int var_idx, unsigned int dot_var_idx);
409 #endif
410 
411 };
412 
413 // --------------------------------------------------------------
414 // DifferentiableSystem inline methods
415 inline
417 {
419  libmesh_assert_equal_to (&(time_solver->system()), this);
420  return *time_solver;
421 }
422 
423 inline
425 {
427  libmesh_assert_equal_to (&(time_solver->system()), this);
428  return *time_solver;
429 }
430 
431 } // namespace libMesh
432 
433 
434 #endif // LIBMESH_DIFF_SYSTEM_H
libMesh::DifferentiableSystem::swap_physics
void swap_physics(DifferentiablePhysics *&swap_physics)
Swap current physics object with external object.
Definition: diff_system.C:353
libMesh::DifferentiableSystem::deltat
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time.
Definition: diff_system.h:260
libMesh::DifferentiableSystem::release_linear_solver
virtual void release_linear_solver(LinearSolver< Number > *) const override
Releases a pointer to a linear solver acquired by this->get_linear_solver()
Definition: diff_system.C:194
libMesh::DifferentiableSystem::print_jacobian_norms
bool print_jacobian_norms
Set print_jacobian_norms to true to print |J| whenever it is assembled.
Definition: diff_system.h:341
libMesh::DifferentiableSystem::print_element_jacobians
bool print_element_jacobians
Set print_element_jacobians to true to print each J_elem contribution.
Definition: diff_system.h:361
libMesh::ImplicitSystem
Manages consistently variables, degrees of freedom, coefficient vectors, and matrices for implicit sy...
Definition: implicit_system.h:57
libMesh::DifferentiableSystem::print_residual_norms
bool print_residual_norms
Set print_residual_norms to true to print |F| whenever it is assembled.
Definition: diff_system.h:331
libMesh::DifferentiableSystem::clone_physics
virtual std::unique_ptr< DifferentiablePhysics > clone_physics() override
We don't allow systems to be attached to each other.
Definition: diff_system.h:158
libMesh::DifferentiableSystem::postprocess
virtual void postprocess()
Executes a postprocessing loop over all elements, and if postprocess_sides is true over all sides.
Definition: diff_system.h:276
libMesh::DifferentiableSystem::sys_type
DifferentiableSystem sys_type
The type of system.
Definition: diff_system.h:75
libMesh::DifferentiableSystem::print_element_solutions
bool print_element_solutions
Set print_element_solutions to true to print each U_elem input.
Definition: diff_system.h:351
libMesh::DifferentiablePhysics
This class provides a specific system class.
Definition: diff_physics.h:76
libMesh::DifferentiableSystem::get_qoi
const DifferentiableQoI * get_qoi() const
Definition: diff_system.h:210
libMesh::DifferentiableSystem::clone
virtual std::unique_ptr< DifferentiableQoI > clone() override
We don't allow systems to be attached to each other.
Definition: diff_system.h:168
libMesh::DifferentiableSystem::print_residuals
bool print_residuals
Set print_residuals to true to print F whenever it is assembled.
Definition: diff_system.h:336
libMesh::DifferentiableSystem::reinit
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:96
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::DifferentiableSystem::DifferentiableSystem
DifferentiableSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
Definition: diff_system.C:34
libMesh::DifferentiableSystem::attach_qoi
void attach_qoi(DifferentiableQoI *qoi_in)
Attach external QoI object.
Definition: diff_system.h:224
libMesh::DifferentiableSystem::print_solution_norms
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:320
libMesh::System::number
unsigned int number() const
Definition: system.h:2075
libMesh::DifferentiableSystem
This class provides a specific system class.
Definition: diff_system.h:53
libMesh::DifferentiableSystem::time_solver
std::unique_ptr< TimeSolver > time_solver
A pointer to the solver object we're going to use.
Definition: diff_system.h:233
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::DifferentiableSystem::_diff_physics
DifferentiablePhysics * _diff_physics
Pointer to object to use for physics assembly evaluations.
Definition: diff_system.h:370
libMesh::DifferentiableQoI
This class provides a specific system class.
Definition: diff_qoi.h:52
libMesh::System::qoi
std::vector< Number > qoi
Values of the quantities of interest.
Definition: system.h:1574
libMesh::DifferentiableQoI::init_qoi
virtual void init_qoi(std::vector< Number > &)
Initialize system qoi.
Definition: diff_qoi.h:71
libMesh::QoISet
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
libMesh::DifferentiableSystem::element_postprocess
virtual void element_postprocess(DiffContext &)
Does any work that needs to be done on elem in a postprocessing loop.
Definition: diff_system.h:281
libMesh::DifferentiableSystem::print_jacobians
bool print_jacobians
Set print_jacobians to true to print J whenever it is assembled.
Definition: diff_system.h:346
libMesh::DifferentiableSystem::solve
virtual void solve() override
Invokes the solver associated with the system.
Definition: diff_system.C:152
libMesh::DifferentiableSystem::get_time_solver
TimeSolver & get_time_solver()
Definition: diff_system.h:416
libMesh::DifferentiableSystem::assemble
virtual void assemble() override
Prepares matrix and rhs for matrix assembly.
Definition: diff_system.C:145
libMesh::DifferentiableSystem::~DifferentiableSystem
virtual ~DifferentiableSystem()
Destructor.
Definition: diff_system.C:56
libMesh::DifferentiableSystem::get_linear_solver
virtual LinearSolver< Number > * get_linear_solver() const override
Definition: diff_system.C:175
libMesh::DifferentiableSystem::clear
virtual void clear() override
Clear all the data structures associated with the system.
Definition: diff_system.C:69
libMesh::DifferentiableSystem::adjoint_solve
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:164
libMesh::DifferentiableSystem::add_second_order_dot_vars
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:198
libMesh::DifferentiableSystem::Parent
ImplicitSystem Parent
The type of the parent.
Definition: diff_system.h:80
libMesh::DiffContext
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::DifferentiableSystem::have_first_order_scalar_vars
bool have_first_order_scalar_vars() const
Check for any first order vars that are also belong to FEFamily::SCALAR.
Definition: diff_system.C:327
libMesh::DifferentiableSystem::attach_physics
void attach_physics(DifferentiablePhysics *physics_in)
Attach external Physics object.
Definition: diff_system.h:196
libMesh::DifferentiableSystem::diff_qoi
DifferentiableQoI * diff_qoi
Pointer to object to use for quantity of interest assembly evaluations.
Definition: diff_system.h:377
libMesh::System::name
const std::string & name() const
Definition: system.h:2067
libMesh::LinearSolver< Number >
libMesh::DifferentiableSystem::get_physics
const DifferentiablePhysics * get_physics() const
Definition: diff_system.h:181
libMesh::DifferentiableSystem::get_qoi
DifferentiableQoI * get_qoi()
Definition: diff_system.h:218
libMesh::DifferentiableSystem::assembly
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.
libMesh::DifferentiableSystem::init_data
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:108
libMesh::TimeSolver
This is a generic class that defines a solver to handle time integration of DifferentiableSystems.
Definition: time_solver.h:58
libMesh::DifferentiableSystem::postprocess_sides
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:314
libMesh::DifferentiableSystem::print_solutions
bool print_solutions
Set print_solutions to true to print U whenever it is used in an assembly() call.
Definition: diff_system.h:326
libMesh::DifferentiableSystem::get_second_order_dot_var
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:310
libMesh::DifferentiableSystem::side_postprocess
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:287
libMesh::DifferentiableSystem::get_linear_solve_parameters
virtual std::pair< unsigned int, Real > get_linear_solve_parameters() const override
Definition: diff_system.C:184
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::DifferentiableSystem::set_time_solver
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:241
libMesh::DifferentiablePhysics::init_physics
virtual void init_physics(const System &sys)
Initialize any data structures associated with the physics.
Definition: diff_physics.C:40
libMesh::DifferentiableSystem::print_element_residuals
bool print_element_residuals
Set print_element_residuals to true to print each R_elem contribution.
Definition: diff_system.h:356
libMesh::DifferentiablePhysics::clone_physics
virtual std::unique_ptr< DifferentiablePhysics > clone_physics()=0
Copy of this object.
libMesh::DifferentiableSystem::get_physics
DifferentiablePhysics * get_physics()
Definition: diff_system.h:190
libMesh::DifferentiableQoI::clone
virtual std::unique_ptr< DifferentiableQoI > clone()=0
Copy of this object.
libMesh::DifferentiableSystem::build_context
virtual std::unique_ptr< DiffContext > build_context()
Builds a DiffContext object with enough information to do evaluations on each element.
Definition: diff_system.C:137
libMesh::DifferentiableSystem::add_dot_var_dirichlet_bcs
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:230
libMesh::DifferentiableSystem::have_second_order_scalar_vars
bool have_second_order_scalar_vars() const
Check for any second order vars that are also belong to FEFamily::SCALAR.
Definition: diff_system.C:339