libMesh
heatsystem.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 #include "libmesh/enum_fe_family.h"
21 #include "libmesh/fem_system.h"
22 #include "libmesh/parameter_pointer.h"
23 #include "libmesh/parameter_vector.h"
24 
25 using namespace libMesh;
26 
27 // FEMSystem, TimeSolver and NewtonSolver will handle most tasks,
28 // but we must specify element residuals
29 class HeatSystem : public FEMSystem
30 {
31 public:
32  // Constructor
34  const std::string & name_in,
35  const unsigned int number_in)
36  : FEMSystem(es, name_in, number_in),
37  _k(1.0),
38  _fe_family("LAGRANGE"),
39  _fe_order(1),
40  _analytic_jacobians(true)
41  { this->init_qois(2); QoI_time_instant.resize(2);}
42 
43  Real & k() { return _k; }
44  bool & analytic_jacobians() { return _analytic_jacobians; }
45 
46  // A vector to specify which QoIs are instantaneous and what time instant they are evaluated at
47  std::vector<Real> QoI_time_instant;
48 
49 protected:
50  // System initialization
51  virtual void init_data ();
52 
53  // Context initialization
54  virtual void init_context (DiffContext & context);
55 
56  // Element residual and jacobian calculations
57  // Time dependent parts
58  virtual bool element_time_derivative (bool request_jacobian,
59  DiffContext & context);
60 
61  // Library evaluation of QoIs
62  virtual void element_qoi (DiffContext & context, const QoISet & qois);
63 
64  // RHS for adjoint problem
65  virtual void element_qoi_derivative (DiffContext & context,
66  const QoISet & /* qois */);
67 
68  // Parameters associated with the system
69  std::vector<Number> parameters;
70 
71  // The ParameterVector object that will contain pointers to
72  // the system parameters
73  ParameterVector parameter_vector;
74 
75  // The parameters to solve for
76  Real _k;
77 
78  // The FE type to use
79  std::string _fe_family;
80  unsigned int _fe_order;
81 
82  // Calculate Jacobians analytically or not?
83  bool _analytic_jacobians;
84 
85  // Use averaged model
87 
88 };
This is the EquationSystems class.
Real & k()
Definition: heatsystem.h:43
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
std::vector< Real > QoI_time_instant
Definition: heatsystem.h:47
Data structure for specifying which Parameters should be independent variables in a parameter sensiti...
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
The libMesh namespace provides an interface to certain functionality in the library.
This class provides a specific system class.
Definition: fem_system.h:53
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool _averaged_model
Definition: heatsystem.h:86
bool & analytic_jacobians()
Definition: heatsystem.h:44
HeatSystem(EquationSystems &es, const std::string &name_in, const unsigned int number_in)
Definition: heatsystem.h:33