libMesh
heatsystem.C
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/getpot.h"
21 
22 #include "heatsystem.h"
23 
24 #include "libmesh/fe_base.h"
25 #include "libmesh/fe_interface.h"
26 #include "libmesh/fem_context.h"
27 #include "libmesh/mesh.h"
28 #include "libmesh/quadrature.h"
29 #include "libmesh/string_to_enum.h"
30 #include "libmesh/system.h"
31 #include "libmesh/equation_systems.h"
32 #include "libmesh/zero_function.h"
33 #include "libmesh/dirichlet_boundaries.h"
34 #include "libmesh/dof_map.h"
35 #include "libmesh/numeric_vector.h"
36 
38 {
39  T_var = this->add_variable ("T", static_cast<Order>(_fe_order),
40  Utility::string_to_enum<FEFamily>(_fe_family));
41 
42  // Make sure the input file heat.in exists, and parse it.
43  {
44  std::ifstream i("heat.in");
45  libmesh_error_msg_if(!i, '[' << this->processor_id() << "] Can't find heat.in; exiting early.");
46  }
47  GetPot infile("heat.in");
48  _k = infile("k", 1.0);
49  _analytic_jacobians = infile("analytic_jacobians", true);
50 
51  parameters.push_back(_k);
52 
53  // Set equation system parameters _k and theta, so they can be read by the exact solution
54  this->get_equation_systems().parameters.set<Real>("_k") = _k;
55 
56  // The temperature is evolving, with a first order time derivative
57  this->time_evolving(T_var, 1);
58 
59  ZeroFunction<Number> zero;
60 
61 #ifdef LIBMESH_ENABLE_DIRICHLET
62  // Most DirichletBoundary users will want to supply a "locally
63  // indexed" functor
65  (DirichletBoundary ({0,1,2,3}, {T_var}, zero,
67 #endif
68 
69  FEMSystem::init_data();
70 }
71 
72 
73 
74 void HeatSystem::init_context(DiffContext & context)
75 {
76  FEMContext & c = cast_ref<FEMContext &>(context);
77 
78  FEBase * elem_fe = nullptr;
79  c.get_element_fe(0, elem_fe);
80 
81  // Now make sure we have requested all the data
82  // we need to build the linear system.
83  elem_fe->get_JxW();
84  elem_fe->get_dphi();
85 
86  // Don't waste time on side computations for T
87  FEBase * side_fe = nullptr;
88  c.get_side_fe(0, side_fe);
89  side_fe->get_nothing();
90 
91  // We'll have a more automatic solution to preparing adjoint
92  // solutions for time integration, eventually...
93  if (c.is_adjoint())
94  {
95  // A reference to the system context is built with
96  const System & sys = c.get_system();
97 
98  // Get a pointer to the adjoint solution vector
99  NumericVector<Number> & adjoint_solution =
100  const_cast<System &>(sys).get_adjoint_solution(0);
101 
102  // Add this adjoint solution to the vectors that diff context should localize
103  c.add_localized_vector(adjoint_solution, sys);
104  }
105 
106  FEMSystem::init_context(context);
107 }
108 
109 #define optassert(X) {if (!(X)) libmesh_error_msg("Assertion " #X " failed.");}
110 //#define optassert(X) libmesh_assert(X);
111 
112 bool HeatSystem::element_time_derivative (bool request_jacobian,
113  DiffContext & context)
114 {
115  bool compute_jacobian = request_jacobian && _analytic_jacobians;
116 
117  FEMContext & c = cast_ref<FEMContext &>(context);
118 
119  // First we get some references to cell-specific data that
120  // will be used to assemble the linear system.
121  FEBase * elem_fe = nullptr;
122  c.get_element_fe(0, elem_fe);
123 
124  // Element Jacobian * quadrature weights for interior integration
125  const std::vector<Real> & JxW = elem_fe->get_JxW();
126 
127  // Element basis functions
128  const std::vector<std::vector<RealGradient>> & dphi = elem_fe->get_dphi();
129 
130  // The number of local degrees of freedom in each variable
131  const unsigned int n_u_dofs = c.n_dof_indices(0);
132 
133  // The subvectors and submatrices we need to fill:
134  DenseSubMatrix<Number> & K = c.get_elem_jacobian(0, 0);
135  DenseSubVector<Number> & F = c.get_elem_residual(0);
136 
137  // Now we will build the element Jacobian and residual.
138  // Constructing the residual requires the solution and its
139  // gradient from the previous timestep. This must be
140  // calculated at each quadrature point by summing the
141  // solution degree-of-freedom values by the appropriate
142  // weight functions.
143  unsigned int n_qpoints = c.get_element_qrule().n_points();
144 
145  for (unsigned int qp=0; qp != n_qpoints; qp++)
146  {
147  // Compute the solution gradient at the Newton iterate
148  Gradient grad_T = c.interior_gradient(0, qp);
149 
150  for (unsigned int i=0; i != n_u_dofs; i++)
151  F(i) += JxW[qp] * -parameters[0] * (grad_T * dphi[i][qp]);
152  if (compute_jacobian)
153  for (unsigned int i=0; i != n_u_dofs; i++)
154  for (unsigned int j=0; j != n_u_dofs; ++j)
155  K(i,j) += JxW[qp] * -parameters[0] * (dphi[i][qp] * dphi[j][qp]);
156  } // end of the quadrature point qp-loop
157 
158  return compute_jacobian;
159 }
160 
161 // Perturb and accumulate dual weighted residuals
162 void HeatSystem::perturb_accumulate_residuals(ParameterVector & parameters_in)
163 {
164  this->update();
165 
166  // Get the timestep size
167  deltat_vector.push_back(this->deltat);
168 
169  for (std::size_t j=0, Np = parameters_in.size(); j != Np; ++j)
170  {
171  Number old_parameter = *parameters_in[j];
172 
173  *parameters_in[j] = old_parameter - dp;
174 
175  this->assembly(true, false);
176 
177  this->rhs->close();
178 
179  std::unique_ptr<NumericVector<Number>> R_minus = this->rhs->clone();
180 
181  // The contribution at a single time step would be [f(z;p+dp) - <partialu/partialt, z>(p+dp) - <g(u),z>(p+dp)] * dt
182  // dt can vary per timestep, so its stored and multiplied separately during the integration
183  R_minus_dp.push_back(-R_minus->dot(this->get_adjoint_solution(0)));
184 
185  *parameters_in[j] = old_parameter + dp;
186 
187  this->assembly(true, false);
188 
189  this->rhs->close();
190 
191  std::unique_ptr<NumericVector<Number>> R_plus = this->rhs->clone();
192 
193  R_plus_dp.push_back(-R_plus->dot(this->get_adjoint_solution(0)));
194 
195  *parameters_in[j] = old_parameter;
196  }
197 }
std::vector< Number > deltat_vector
Definition: heatsystem.h:143
virtual bool element_time_derivative(bool request_jacobian, DiffContext &context)
Adds the time derivative contribution on elem to elem_residual.
Definition: heatsystem.C:78
const EquationSystems & get_equation_systems() const
Definition: system.h:721
NumericVector< Number > * rhs
The system matrix.
virtual std::unique_ptr< NumericVector< T > > clone() const =0
virtual void time_evolving(unsigned int var, unsigned int order)
Tells the DiffSystem that variable var is evolving with respect to time.
Definition: diff_physics.C:41
const Number zero
.
Definition: libmesh.h:304
unsigned int T_var
Definition: heatsystem.h:133
void compute_jacobian(const NumericVector< Number > &, SparseMatrix< Number > &J, NonlinearImplicitSystem &system)
Definition: assembly.C:315
std::vector< Number > R_plus_dp
Definition: heatsystem.h:139
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time...
Definition: diff_system.h:280
bool _analytic_jacobians
Definition: heatsystem.h:136
Number dp
Definition: heatsystem.h:146
virtual void init_context(DiffContext &context)
Definition: heatsystem.C:48
NumberVectorValue Gradient
unsigned int add_variable(std::string_view var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
Definition: system.C:1357
std::vector< Number > parameters
Definition: heatsystem.h:113
FEGenericBase< Real > FEBase
void perturb_accumulate_residuals(ParameterVector &parameters)
Definition: heatsystem.C:162
virtual void close()=0
Calls the NumericVector&#39;s internal assembly routines, ensuring that the values are consistent across ...
virtual void update()
Update the local values to reflect the solution on neighboring processors.
Definition: system.C:510
std::vector< Number > R_minus_dp
Definition: heatsystem.h:140
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
T & set(const std::string &)
Definition: parameters.h:469
std::string _fe_family
Definition: heatsystem.h:129
unsigned int _fe_order
Definition: heatsystem.h:130
virtual void init_data()
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used...
Definition: heatsystem.C:37
NumericVector< Number > & get_adjoint_solution(unsigned int i=0)
Definition: system.C:1245
Parameters parameters
Data structure holding arbitrary parameters.
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
virtual void assembly(bool get_residual, bool get_jacobian, bool apply_heterogeneous_constraints=false, bool apply_no_constraints=false) override
Prepares matrix or rhs for matrix assembly.
Definition: fem_system.C:880
processor_id_type processor_id() const
System(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
Definition: system.C:64
const DofMap & get_dof_map() const
Definition: system.h:2374
template class LIBMESH_EXPORT NumericVector< Number >