#include "heatsystem.hpp" #include "libmesh/dirichlet_boundaries.h" #include "libmesh/dof_map.h" #include "libmesh/elem.h" #include "libmesh/fe_base.h" #include "libmesh/fe_interface.h" #include "libmesh/fem_context.h" #include "libmesh/getpot.h" #include "libmesh/mesh.h" #include "libmesh/quadrature.h" #include "libmesh/string_to_enum.h" #include "libmesh/zero_function.h" using namespace libMesh; // user defined literal for libmesh's Real constexpr auto operator"" _r(long double n) { return libMesh::Real(n); } void HeatSystem::init_data() { T_var = this->add_variable("T", static_cast(_fe_order), Utility::string_to_enum(_fe_family)); ZeroFunction zero; this->get_dof_map().add_dirichlet_boundary( DirichletBoundary({0, 1, 2}, {T_var}, &zero)); // Do the parent's initialization after variables are defined FEMSystem::init_data(); this->time_evolving(0); } void HeatSystem::init_context(DiffContext& context) { auto& c = libmesh_cast_ref(context); const auto& elem_dims = c.elem_dimensions(); for (const auto& dim : elem_dims) { FEBase* fe = nullptr; c.get_element_fe(T_var, fe, dim); fe->get_JxW(); // For integration fe->get_dphi(); // For bilinear form fe->get_xyz(); // For forcing fe->get_phi(); // For forcing } FEMSystem::init_context(context); } bool HeatSystem::element_time_derivative(bool request_jacobian, DiffContext& context) { auto& c = libmesh_cast_ref(context); const auto& mesh_dim = c.get_system().get_mesh().mesh_dimension(); // First we get some references to cell-specific data that // will be used to assemble the linear system. const auto& dim = c.get_elem().dim(); FEBase* fe = nullptr; c.get_element_fe(T_var, fe, dim); // Element Jacobian * quadrature weights for interior integration const auto& JxW = fe->get_JxW(); const auto& xyz = fe->get_xyz(); const auto& phi = fe->get_phi(); const auto& dphi = fe->get_dphi(); // The number of local degrees of freedom in each variable const auto& n_T_dofs = c.get_dof_indices(T_var).size(); // The subvectors and submatrices we need to fill: auto& K = c.get_elem_jacobian(T_var, T_var); auto& F = c.get_elem_residual(T_var); // Now we will build the element Jacobian and residual. // Constructing the residual requires the solution and its // gradient from the previous timestep. This must be // calculated at each quadrature point by summing the // solution degree-of-freedom values by the appropriate // weight functions. auto n_qpoints = c.get_element_qrule().n_points(); for (std::size_t qp{0}; qp != n_qpoints; qp++) { // Compute the solution gradient at the Newton iterate Gradient grad_T = c.interior_gradient(T_var, qp); const auto& k = _k[dim]; const auto& u_exact = +1.0_r; // Only apply forcing to interior elements const auto& forcing = (dim == mesh_dim) ? -k * u_exact * (dim * libMesh::pi * libMesh::pi) : 0; const auto& JxWxNK = JxW[qp] * -k; for (std::size_t i{0}; i != n_T_dofs; i++) F(i) += JxWxNK * (grad_T * dphi[i][qp] + forcing * phi[i][qp]); if (request_jacobian) { const auto& JxWxNKxD = JxWxNK * context.get_elem_solution_derivative(); for (std::size_t i{0}; i != n_T_dofs; i++) for (std::size_t j{0}; j != n_T_dofs; ++j) K(i, j) += JxWxNKxD * (dphi[i][qp] * dphi[j][qp]); } } // end of the quadrature point qp-loop return request_jacobian; }