Loading [MathJax]/extensions/tex2jax.js
libMesh
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
heatsystem.C
Go to the documentation of this file.
1 #include "heatsystem.h"
2 
3 #include "libmesh/dirichlet_boundaries.h"
4 #include "libmesh/dof_map.h"
5 #include "libmesh/fe_base.h"
6 #include "libmesh/fe_interface.h"
7 #include "libmesh/fem_context.h"
8 #include "libmesh/getpot.h"
9 #include "libmesh/mesh.h"
10 #include "libmesh/quadrature.h"
11 #include "libmesh/string_to_enum.h"
12 #include "libmesh/zero_function.h"
13 #include "libmesh/elem.h"
14 
15 using namespace libMesh;
16 
18 {
19  T_var = this->add_variable("T", static_cast<Order>(_fe_order),
20  Utility::string_to_enum<FEFamily>(_fe_family));
21 
22 #ifdef LIBMESH_ENABLE_DIRICHLET
23  const unsigned int dim = this->get_mesh().mesh_dimension();
24 
25  // Add dirichlet boundaries on all but the boundary element side
26  const boundary_id_type all_ids[6] = {0, 1, 2, 3, 4, 5};
27  std::set<boundary_id_type> nonyplus_bdys(all_ids, all_ids+(dim*2));
28  const boundary_id_type yplus_id = (dim == 3) ? 3 : 2;
29  nonyplus_bdys.erase(yplus_id);
30 
32 
33  // Most DirichletBoundary users will want to supply a "locally
34  // indexed" functor
35  this->get_dof_map().add_dirichlet_boundary
36  (DirichletBoundary (nonyplus_bdys, {T_var}, zero, LOCAL_VARIABLE_ORDER));
37 #endif // LIBMESH_ENABLE_DIRICHLET
38 
39  // Do the parent's initialization after variables are defined
41 
42  // The temperature is evolving, with a first-order time derivative
43  this->time_evolving(T_var, 1);
44 }
45 
46 
47 
49 {
50  FEMContext & c = cast_ref<FEMContext &>(context);
51 
52  const std::set<unsigned char> & elem_dims =
53  c.elem_dimensions();
54 
55  for (std::set<unsigned char>::const_iterator dim_it =
56  elem_dims.begin(); dim_it != elem_dims.end(); ++dim_it)
57  {
58  const unsigned char dim = *dim_it;
59 
60  FEBase * fe = nullptr;
61 
62  c.get_element_fe(T_var, fe, dim);
63 
64  fe->get_JxW(); // For integration
65  fe->get_dphi(); // For bilinear form
66  fe->get_xyz(); // For forcing
67  fe->get_phi(); // For forcing
68 
69  FEBase * side_fe = nullptr;
70  c.get_side_fe(T_var, side_fe, dim);
71  side_fe->get_nothing();
72  }
73 
74  FEMSystem::init_context(context);
75 }
76 
77 
78 bool HeatSystem::element_time_derivative (bool request_jacobian,
79  DiffContext & context)
80 {
81  FEMContext & c = cast_ref<FEMContext &>(context);
82 
83  const unsigned int mesh_dim =
85 
86  // First we get some references to cell-specific data that
87  // will be used to assemble the linear system.
88  const unsigned short dim = c.get_elem().dim();
89  FEBase * fe = nullptr;
90  c.get_element_fe(T_var, fe, dim);
91 
92  // Element Jacobian * quadrature weights for interior integration
93  const std::vector<Real> & JxW = fe->get_JxW();
94 
95  const std::vector<Point> & xyz = fe->get_xyz();
96 
97  const std::vector<std::vector<Real>> & phi = fe->get_phi();
98 
99  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
100 
101  // The number of local degrees of freedom in each variable
102  const unsigned int n_T_dofs = c.n_dof_indices(T_var);
103 
104  // The subvectors and submatrices we need to fill:
105  DenseSubMatrix<Number> & K = c.get_elem_jacobian(T_var, T_var);
107 
108  // Now we will build the element Jacobian and residual.
109  // Constructing the residual requires the solution and its
110  // gradient from the previous timestep. This must be
111  // calculated at each quadrature point by summing the
112  // solution degree-of-freedom values by the appropriate
113  // weight functions.
114  unsigned int n_qpoints = c.get_element_qrule().n_points();
115 
116  for (unsigned int qp=0; qp != n_qpoints; qp++)
117  {
118  // Compute the solution gradient at the Newton iterate
119  Gradient grad_T = c.interior_gradient(T_var, qp);
120 
121  const Number k = _k[dim];
122 
123  const Point & p = xyz[qp];
124 
125  // solution + laplacian depend on problem dimension
126  const Number u_exact = (mesh_dim == 2) ?
127  std::sin(libMesh::pi*p(0)) * std::sin(libMesh::pi*p(1)) :
128  std::sin(libMesh::pi*p(0)) * std::sin(libMesh::pi*p(1)) *
129  std::sin(libMesh::pi*p(2));
130 
131  // Only apply forcing to interior elements
132  const Number forcing = (dim == mesh_dim) ?
133  -k * u_exact * (dim * libMesh::pi * libMesh::pi) : 0;
134 
135  const Number JxWxNK = JxW[qp] * -k;
136 
137  for (unsigned int i=0; i != n_T_dofs; i++)
138  F(i) += JxWxNK * (grad_T * dphi[i][qp] + forcing * phi[i][qp]);
139  if (request_jacobian)
140  {
141  const Number JxWxNKxD = JxWxNK *
143 
144  for (unsigned int i=0; i != n_T_dofs; i++)
145  for (unsigned int j=0; j != n_T_dofs; ++j)
146  K(i,j) += JxWxNKxD * (dphi[i][qp] * dphi[j][qp]);
147  }
148  } // end of the quadrature point qp-loop
149 
150  return request_jacobian;
151 }
Real get_elem_solution_derivative() const
The derivative of the current elem_solution w.r.t.
Definition: diff_context.h:432
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
Definition: diff_context.h:274
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
void get_side_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for edge/face (2D/3D) finite element object for variable var for the largest dimension in th...
Definition: fem_context.h:317
ConstFunction that simply returns 0.
Definition: zero_function.h:38
const Elem & get_elem() const
Accessor for current Elem object.
Definition: fem_context.h:908
unsigned int dim
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 init_context(DiffContext &) override
Definition: fem_system.C:1370
unsigned int n_dof_indices() const
Total number of dof indices on the element.
Definition: diff_context.h:395
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
The libMesh namespace provides an interface to certain functionality in the library.
const Number zero
.
Definition: libmesh.h:287
const MeshBase & get_mesh() const
Definition: system.h:2354
Defines a dense subvector for use in finite element computations.
const std::vector< std::vector< OutputGradient > > & get_dphi() const
Definition: fe_base.h:230
Gradient interior_gradient(unsigned int var, unsigned int qp) const
Definition: fem_context.C:467
int8_t boundary_id_type
Definition: id_types.h:51
const System & get_system() const
Accessor for associated system.
Definition: diff_context.h:106
virtual void init_context(DiffContext &context)
Definition: heatsystem.C:48
Defines a dense submatrix for use in Finite Element-type computations.
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used...
Definition: fem_system.C:873
virtual_for_inffe const std::vector< Real > & get_JxW() const
Definition: fe_abstract.h:303
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62
unsigned int n_points() const
Definition: quadrature.h:131
virtual_for_inffe const std::vector< Point > & get_xyz() const
Definition: fe_abstract.h:280
void get_nothing() const
Definition: fe_abstract.h:269
const DenseVector< Number > & get_elem_residual() const
Const accessor for element residual.
Definition: diff_context.h:242
virtual unsigned short dim() const =0
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
unsigned int mesh_dimension() const
Definition: mesh_base.C:354
void get_element_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for interior finite element object for variable var for the largest dimension in the mesh...
Definition: fem_context.h:277
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
This class forms the foundation from which generic finite elements may be derived.
const QBase & get_element_qrule() const
Accessor for element interior quadrature rule for the dimension of the current _elem.
Definition: fem_context.h:802
const std::set< unsigned char > & elem_dimensions() const
Definition: fem_context.h:951
const Real pi
.
Definition: libmesh.h:281
const std::vector< std::vector< OutputShape > > & get_phi() const
Definition: fe_base.h:207