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