20 #include "libmesh/getpot.h" 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" 41 T_var =
dynamic_cast<FEMSystem&
>(sys).add_variable (
"T", static_cast<Order>(1),
LAGRANGE);
43 GetPot infile(
"heat.in");
44 _k = infile(
"k", 1.0);
50 ZeroFunction<Number>
zero;
52 ConstFunction<Number> one(1.0);
55 sys.get_dof_map().add_dirichlet_boundary(DirichletBoundary ({0,1,2,3}, {
T_var}, &one));
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);
65 FEMContext & c = cast_ref<FEMContext &>(context);
67 FEBase * elem_fe =
nullptr;
68 c.get_element_fe(
T_var, elem_fe);
78 FEBase * side_fe =
nullptr;
79 c.get_side_fe(
T_var, side_fe);
80 side_fe->get_nothing();
87 const System & sys = c.get_system();
91 const_cast<System &
>(sys).get_adjoint_solution(0);
94 c.add_localized_vector(adjoint_solution0, sys);
98 const_cast<System &
>(sys).get_adjoint_solution(1);
101 c.add_localized_vector(adjoint_solution1, sys);
107 DiffContext & context)
111 FEMContext & c = cast_ref<FEMContext &>(context);
115 FEBase * elem_fe =
nullptr;
116 c.get_element_fe(
T_var, elem_fe);
119 const std::vector<Real> & JxW = elem_fe->get_JxW();
122 const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
123 const std::vector<std::vector<RealGradient>> & dphi = elem_fe->get_dphi();
126 const unsigned int n_u_dofs = c.n_dof_indices(
T_var);
129 DenseSubMatrix<Number> & K = c.get_elem_jacobian(
T_var,
T_var);
130 DenseSubVector<Number> & F = c.get_elem_residual(
T_var);
133 const std::vector<Point > & q_point = elem_fe->get_xyz();
141 unsigned int n_qpoints = c.get_element_qrule().n_points();
149 for (
unsigned int qp=0; qp != n_qpoints; qp++)
155 const Real x = q_point[qp](0);
160 for (
unsigned int i=0; i != n_u_dofs; i++)
162 F(i) += JxW[qp] * ( ( -sigma * (grad_T * dphi[i][qp]) ) + (f * phi[i][qp]) );
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]);
176 return std::make_unique<SigmaPhysics>(*this);
virtual void init_context(DiffContext &context) override
virtual void init_data(System &sys)
virtual void time_evolving(unsigned int var, unsigned int order)
Tells the DiffSystem that variable var is evolving with respect to time.
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)
NumberVectorValue Gradient
FEGenericBase< Real > FEBase
virtual std::unique_ptr< DifferentiablePhysics > clone_physics() override
Copy of this object.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
template class LIBMESH_EXPORT NumericVector< Number >