20 #include "libmesh/getpot.h"    22 #include "heatsystem.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"    42     std::ifstream i(
"heat.in");
    43     libmesh_error_msg_if(!i, 
'[' << this->
processor_id() << 
"] Can't find heat.in; exiting early.");
    45   GetPot infile(
"heat.in");
    46   _k = infile(
"k", 1.0);
    53   ZeroFunction<Number> 
zero;
    55   ConstFunction<Number> one(1.0);
    57 #ifdef LIBMESH_ENABLE_DIRICHLET    61     (DirichletBoundary ({0,1,2,3}, {  0 }, one,
    68   FEMSystem::init_data();
    73   FEMContext & c = cast_ref<FEMContext &>(context);
    75   FEBase * elem_fe = 
nullptr;
    76   c.get_element_fe(0, elem_fe);
    86   FEBase * side_fe = 
nullptr;
    87   c.get_side_fe(0, side_fe);
    88   side_fe->get_nothing();
    95       const System & sys = c.get_system();
   102       c.add_localized_vector(adjoint_solution0, sys);
   109       c.add_localized_vector(adjoint_solution1, sys);
   112   FEMSystem::init_context(context);
   116                                           DiffContext & context)
   120   FEMContext & c = cast_ref<FEMContext &>(context);
   124   FEBase * elem_fe = 
nullptr;
   125   c.get_element_fe(0, elem_fe);
   128   const std::vector<Real> & JxW = elem_fe->get_JxW();
   131   const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
   132   const std::vector<std::vector<RealGradient>> & dphi = elem_fe->get_dphi();
   135   const unsigned int n_u_dofs = c.n_dof_indices(0);
   138   DenseSubMatrix<Number> & K = c.get_elem_jacobian(0, 0);
   139   DenseSubVector<Number> & F = c.get_elem_residual(0);
   142   const std::vector<Point > & q_point = elem_fe->get_xyz();
   150   unsigned int n_qpoints = c.get_element_qrule().n_points();
   158   for (
unsigned int qp=0; qp != n_qpoints; qp++)
   161       Gradient grad_T = c.interior_gradient(0, qp);
   164       const Real x = q_point[qp](0);
   176       for (
unsigned int i=0; i != n_u_dofs; i++)
   178 F(i) += JxW[qp] * ( ( -sigma * (grad_T * dphi[i][qp]) ) + (f * phi[i][qp]) );
   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]);
 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. 
virtual void time_evolving(unsigned int var, unsigned int order)
Tells the DiffSystem that variable var is evolving with respect to time. 
void compute_jacobian(const NumericVector< Number > &, SparseMatrix< Number > &J, NonlinearImplicitSystem &system)
virtual void init_context(DiffContext &context)
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...
NumericVector< Number > & get_adjoint_solution(unsigned int i=0)
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system. 
processor_id_type processor_id() const
System(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor. 
const DofMap & get_dof_map() const
template class LIBMESH_EXPORT NumericVector< Number >