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/const_function.h"
34 #include "libmesh/dirichlet_boundaries.h"
35 #include "libmesh/dof_map.h"
36 #include "libmesh/numeric_vector.h"
37 
39 {
40  // Make sure the input file heat.in exists, and parse it.
41  {
42  std::ifstream i("heat.in");
43  libmesh_error_msg_if(!i, '[' << this->processor_id() << "] Can't find heat.in; exiting early.");
44  }
45  GetPot infile("heat.in");
46  _k = infile("k", 1.0);
47  _analytic_jacobians = infile("analytic_jacobians", true);
48  _averaged_model = infile("averaged_model", false);
49 
50  // The temperature is evolving, with a first order time derivative
51  this->time_evolving(0, 1);
52 
53  ZeroFunction<Number> zero;
54 
55  ConstFunction<Number> one(1.0);
56 
57 #ifdef LIBMESH_ENABLE_DIRICHLET
58  // Most DirichletBoundary users will want to supply a "locally
59  // indexed" functor
61  (DirichletBoundary ({0,1,2,3}, { /* T_var = */ 0 }, one,
63 #endif
64 
65  this->get_dof_map().add_adjoint_dirichlet_boundary(DirichletBoundary ({0,1,2,3}, {0}, &zero), 0);
66  this->get_dof_map().add_adjoint_dirichlet_boundary(DirichletBoundary ({0,1,2,3}, {0}, &zero), 1);
67 
68  FEMSystem::init_data();
69 }
70 
71 void HeatSystem::init_context(DiffContext & context)
72 {
73  FEMContext & c = cast_ref<FEMContext &>(context);
74 
75  FEBase * elem_fe = nullptr;
76  c.get_element_fe(0, elem_fe);
77 
78  // Now make sure we have requested all the data
79  // we need to build the linear system.
80  elem_fe->get_JxW();
81  elem_fe->get_dphi();
82  elem_fe->get_phi();
83  elem_fe->get_xyz();
84 
85  // Don't waste time on side computations for T
86  FEBase * side_fe = nullptr;
87  c.get_side_fe(0, side_fe);
88  side_fe->get_nothing();
89 
90  // We'll have a more automatic solution to preparing adjoint
91  // solutions for time integration, eventually...
92  if (c.is_adjoint())
93  {
94  // A reference to the system context is built with
95  const System & sys = c.get_system();
96 
97  // Get a pointer to the adjoint solution vector
98  NumericVector<Number> & adjoint_solution0 =
99  const_cast<System &>(sys).get_adjoint_solution(0);
100 
101  // Add this adjoint solution to the vectors that diff context should localize
102  c.add_localized_vector(adjoint_solution0, sys);
103 
104  // Get a pointer to the adjoint solution vector
105  NumericVector<Number> & adjoint_solution1 =
106  const_cast<System &>(sys).get_adjoint_solution(1);
107 
108  // Add this adjoint solution to the vectors that diff context should localize
109  c.add_localized_vector(adjoint_solution1, sys);
110  }
111 
112  FEMSystem::init_context(context);
113 }
114 
115 bool HeatSystem::element_time_derivative (bool request_jacobian,
116  DiffContext & context)
117 {
118  bool compute_jacobian = request_jacobian && _analytic_jacobians;
119 
120  FEMContext & c = cast_ref<FEMContext &>(context);
121 
122  // First we get some references to cell-specific data that
123  // will be used to assemble the linear system.
124  FEBase * elem_fe = nullptr;
125  c.get_element_fe(0, elem_fe);
126 
127  // Element Jacobian * quadrature weights for interior integration
128  const std::vector<Real> & JxW = elem_fe->get_JxW();
129 
130  // Element basis functions
131  const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
132  const std::vector<std::vector<RealGradient>> & dphi = elem_fe->get_dphi();
133 
134  // The number of local degrees of freedom in each variable
135  const unsigned int n_u_dofs = c.n_dof_indices(0);
136 
137  // The subvectors and submatrices we need to fill:
138  DenseSubMatrix<Number> & K = c.get_elem_jacobian(0, 0);
139  DenseSubVector<Number> & F = c.get_elem_residual(0);
140 
141  // Quadrature point locations
142  const std::vector<Point > & q_point = elem_fe->get_xyz();
143 
144  // Now we will build the element Jacobian and residual.
145  // Constructing the residual requires the solution and its
146  // gradient from the previous timestep. This must be
147  // calculated at each quadrature point by summing the
148  // solution degree-of-freedom values by the appropriate
149  // weight functions.
150  unsigned int n_qpoints = c.get_element_qrule().n_points();
151 
152  // Conductivity
153  Real sigma = 1.0;
154 
155  // Forcing function
156  Real f = 1.0;
157 
158  for (unsigned int qp=0; qp != n_qpoints; qp++)
159  {
160  // Compute the solution gradient at the Newton iterate
161  Gradient grad_T = c.interior_gradient(0, qp);
162 
163  // Location of the current qp
164  const Real x = q_point[qp](0);
165 
166  // Spatially varying conductivity
167  if(_averaged_model)
168  {
169 sigma = 0.01;
170  }
171  else
172  {
173 sigma = 0.001 + x;
174  }
175 
176  for (unsigned int i=0; i != n_u_dofs; i++)
177  {
178 F(i) += JxW[qp] * ( ( -sigma * (grad_T * dphi[i][qp]) ) + (f * phi[i][qp]) );
179  }
180 
181  if (compute_jacobian)
182  for (unsigned int i=0; i != n_u_dofs; i++)
183  for (unsigned int j=0; j != n_u_dofs; ++j)
184  K(i,j) += JxW[qp] * -sigma * (dphi[i][qp] * dphi[j][qp]);
185  } // end of the quadrature point qp-loop
186 
187  return compute_jacobian;
188 }
void add_adjoint_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary, unsigned int q)
Adds a copy of the specified Dirichlet boundary to the system, corresponding to the adjoint problem d...
virtual bool element_time_derivative(bool request_jacobian, DiffContext &context)
Adds the time derivative contribution on elem to elem_residual.
Definition: heatsystem.C:78
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
void compute_jacobian(const NumericVector< Number > &, SparseMatrix< Number > &J, NonlinearImplicitSystem &system)
Definition: assembly.C:315
bool _analytic_jacobians
Definition: heatsystem.h:136
virtual void init_context(DiffContext &context)
Definition: heatsystem.C:48
NumberVectorValue Gradient
FEGenericBase< Real > FEBase
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
bool _averaged_model
Definition: heatsystem.h:86
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 >