diff_system.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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/auto_ptr.h" // deprecated
25 #include "libmesh/diff_context.h"
26 #include "libmesh/diff_physics.h"
27 #include "libmesh/diff_qoi.h"
29 #include "libmesh/time_solver.h"
30 
31 // C++ includes
32 #include <memory>
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 
65  const std::string & name,
66  const unsigned int number);
67 
71  virtual ~DifferentiableSystem ();
72 
77 
82 
87  virtual void clear () override;
88 
93  virtual void reinit () override;
94 
105  virtual void assemble () override;
106 
111  virtual LinearSolver<Number> * get_linear_solver() const override;
112 
118  virtual std::pair<unsigned int, Real>
119  get_linear_solve_parameters() const override;
120 
125  virtual void release_linear_solver(LinearSolver<Number> *) const override;
126 
137  virtual void assembly (bool get_residual,
138  bool get_jacobian,
139  bool apply_heterogeneous_constraints = false,
140  bool apply_no_constraints = false) override = 0;
141 
147  virtual void solve () override;
148 
153  virtual std::pair<unsigned int, Real>
154  adjoint_solve (const QoISet & qoi_indices = QoISet()) override;
155 
159  virtual std::unique_ptr<DifferentiablePhysics> clone_physics() override
160  {
161  libmesh_not_implemented();
162  // dummy
163  return std::unique_ptr<DifferentiablePhysics>(this);
164  }
165 
169  virtual std::unique_ptr<DifferentiableQoI> clone() override
170  {
171  libmesh_not_implemented();
172  // dummy
173  return std::unique_ptr<DifferentiableQoI>(this);
174  }
175 
183  { return this->_diff_physics; }
184 
192  { return this->_diff_physics; }
193 
198  { this->_diff_physics = (physics_in->clone_physics()).release();
199  this->_diff_physics->init_physics(*this);}
200 
205 
211  const DifferentiableQoI * get_qoi() const
212  { return this->diff_qoi; }
213 
220  { return this->diff_qoi; }
221 
225  void attach_qoi( DifferentiableQoI * qoi_in )
226  { this->diff_qoi = (qoi_in->clone()).release();
227  // User needs to resize qoi system qoi accordingly
228  this->diff_qoi->init_qoi( this->qoi );}
229 
234  std::unique_ptr<TimeSolver> time_solver;
235 
242  void set_time_solver(std::unique_ptr<TimeSolver> _time_solver)
243  {
244  time_solver.reset(_time_solver.release());
245  }
246 
251 
255  const TimeSolver & get_time_solver() const;
256 
262 
271  virtual std::unique_ptr<DiffContext> build_context();
272 
277  virtual void postprocess () {}
278 
282  virtual void element_postprocess (DiffContext &) {}
283 
288  virtual void side_postprocess (DiffContext &) {}
289 
299  unsigned int get_second_order_dot_var( unsigned int var ) const;
300 
304  bool have_first_order_scalar_vars() const;
305 
309  bool have_second_order_scalar_vars() const;
310 
316 
322 
328 
333 
338 
343 
348 
353 
358 
363 
364 protected:
365 
372 
379 
391  virtual void init_data () override;
392 
399 
408  void add_dot_var_dirichlet_bcs( unsigned int var_idx, unsigned int dot_var_idx);
409 
410 };
411 
412 // --------------------------------------------------------------
413 // DifferentiableSystem inline methods
414 inline
416 {
417  libmesh_assert(time_solver.get());
418  libmesh_assert_equal_to (&(time_solver->system()), this);
419  return *time_solver;
420 }
421 
422 inline
424 {
425  libmesh_assert(time_solver.get());
426  libmesh_assert_equal_to (&(time_solver->system()), this);
427  return *time_solver;
428 }
429 
430 } // namespace libMesh
431 
432 
433 #endif // LIBMESH_DIFF_SYSTEM_H
void attach_physics(DifferentiablePhysics *physics_in)
Definition: diff_system.h:197
const DifferentiableQoI * get_qoi() const
Definition: diff_system.h:211
Manages multiples systems of equations.
virtual void assemble() override
Definition: diff_system.C:145
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.
Definition: qoi_set.h:45
std::unique_ptr< TimeSolver > time_solver
Definition: diff_system.h:234
DifferentiableQoI * get_qoi()
Definition: diff_system.h:219
DifferentiablePhysics * get_physics()
Definition: diff_system.h:191
void add_dot_var_dirichlet_bcs(unsigned int var_idx, unsigned int dot_var_idx)
Definition: diff_system.C:227
virtual std::unique_ptr< DiffContext > build_context()
Definition: diff_system.C:137
DifferentiableQoI * diff_qoi
Definition: diff_system.h:378
virtual void reinit() override
Definition: diff_system.C:96
void swap_physics(DifferentiablePhysics *&swap_physics)
Definition: diff_system.C:349
DifferentiablePhysics * _diff_physics
Definition: diff_system.h:371
unsigned int number() const
Definition: system.h:2025
unsigned int get_second_order_dot_var(unsigned int var) const
Definition: diff_system.C:306
virtual void side_postprocess(DiffContext &)
Definition: diff_system.h:288
std::vector< Number > qoi
Definition: system.h:1558
virtual std::unique_ptr< DifferentiablePhysics > clone_physics() override
Definition: diff_system.h:159
DifferentiableSystem sys_type
Definition: diff_system.h:76
void attach_qoi(DifferentiableQoI *qoi_in)
Definition: diff_system.h:225
virtual void element_postprocess(DiffContext &)
Definition: diff_system.h:282
bool have_second_order_scalar_vars() const
Definition: diff_system.C:335
virtual std::unique_ptr< DifferentiableQoI > clone()=0
const DifferentiablePhysics * get_physics() const
Definition: diff_system.h:182
DifferentiableSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Definition: diff_system.C:34
virtual void init_data() override
Definition: diff_system.C:108
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual LinearSolver< Number > * get_linear_solver() const override
Definition: diff_system.C:175
virtual std::pair< unsigned int, Real > get_linear_solve_parameters() const override
Definition: diff_system.C:184
virtual std::unique_ptr< DifferentiableQoI > clone() override
Definition: diff_system.h:169
void set_time_solver(std::unique_ptr< TimeSolver > _time_solver)
Definition: diff_system.h:242
virtual std::unique_ptr< DifferentiablePhysics > clone_physics()=0
virtual void init_physics(const System &sys)
virtual void solve() override
Definition: diff_system.C:152
virtual void clear() override
Definition: diff_system.C:69
const std::string & name() const
Definition: system.h:2017
virtual void release_linear_solver(LinearSolver< Number > *) const override
Definition: diff_system.C:194
bool have_first_order_scalar_vars() const
Definition: diff_system.C:323
virtual std::pair< unsigned int, Real > adjoint_solve(const QoISet &qoi_indices=QoISet()) override
Definition: diff_system.C:164
TimeSolver & get_time_solver()
Definition: diff_system.h:415
virtual void init_qoi(std::vector< Number > &)
Definition: diff_qoi.h:69
Manages consistently variables, degrees of freedom, coefficient vectors, and matrices for implicit sy...