Go to the documentation of this file.
   18 #include "libmesh/getpot.h" 
   22 #include "libmesh/dof_map.h" 
   23 #include "libmesh/fe_base.h" 
   24 #include "libmesh/fe_interface.h" 
   25 #include "libmesh/fem_context.h" 
   26 #include "libmesh/mesh.h" 
   27 #include "libmesh/quadrature.h" 
   28 #include "libmesh/string_to_enum.h" 
   35                              const std::string & name_in,
 
   36                              const unsigned int number_in) :
 
   44   GetPot infile(
"laplace.in");
 
   61   FEMSystem::init_data();
 
   66 #ifdef LIBMESH_ENABLE_DIRICHLET 
   68   std::set<boundary_id_type> boundary_ids(all_ids, all_ids+6);
 
   70   std::vector<unsigned int> vars;
 
   71   vars.push_back(
u_var);
 
   85 #endif // LIBMESH_ENABLE_DIRICHLET 
   90   FEMContext & c = cast_ref<FEMContext &>(context);
 
  114   FEMContext & c = cast_ref<FEMContext &>(context);
 
  124   const std::vector<Real> & JxW = fe->
get_JxW();
 
  127   const std::vector<std::vector<RealGradient>> & phi = fe->
get_phi();
 
  131   const std::vector<std::vector<RealTensor>> & grad_phi = fe->
get_dphi();
 
  133   const std::vector<Point> & qpoint = fe->
get_xyz();
 
  150   for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  162       for (
unsigned int i=0; i != n_u_dofs; i++)
 
  164           Fu(i) += (grad_u.
contract(grad_phi[i][qp]) - f*phi[i][qp])*JxW[qp];
 
  167           if (request_jacobian)
 
  168             for (
unsigned int j=0; j != n_u_dofs; j++)
 
  169               Kuu(i,j) += grad_phi[j][qp].contract(grad_phi[i][qp])*JxW[qp];
 
  173   return request_jacobian;
 
  239   const Real eps = 1.e-3;
 
  
const DenseVector< Number > & get_elem_residual() const
Const accessor for element residual.
 
const QBase & get_element_qrule() const
Accessor for element interior quadrature rule for the dimension of the current _elem.
 
bool print_element_jacobians
Set print_element_jacobians to true to print each J_elem contribution.
 
virtual void time_evolving(unsigned int var)
Tells the DiffSystem that variable var is evolving with respect to time.
 
const std::vector< Point > & get_xyz() const
 
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
 
Defines a dense submatrix for use in Finite Element-type computations.
 
unsigned int n_points() const
 
RealVectorValue RealGradient
 
The libMesh namespace provides an interface to certain functionality in the library.
 
This class forms the foundation from which generic finite elements may be derived.
 
const std::vector< Real > & get_JxW() const
 
const std::vector< std::vector< OutputGradient > > & get_dphi() const
 
LaplaceSystem(EquationSystems &es, const std::string &name_in, const unsigned int number_in)
 
LaplaceExactSolution exact_solution
 
virtual bool element_time_derivative(bool request_jacobian, DiffContext &context)
Adds the time derivative contribution on elem to elem_residual.
 
This class provides a specific system class.
 
RealGradient forcing(const Point &p)
 
unsigned int add_variable(const std::string &var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
 
void init_dirichlet_bcs()
 
bool print_jacobians
Set print_jacobians to true to print J whenever it is assembled.
 
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
 
A Point defines a location in LIBMESH_DIM dimensional Real space.
 
CompareTypes< T, T2 >::supertype contract(const TypeTensor< T2 > &) const
Multiply 2 tensors together to return a scalar, i.e.
 
unsigned int n_dof_indices() const
Total number of dof indices on the element.
 
This class provides all data required for a physics package (e.g.
 
This is the EquationSystems class.
 
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.
 
Gradient interior_gradient(unsigned int var, unsigned int qp) const
 
Real verify_analytic_jacobians
If verify_analytic_jacobian is equal to zero (as it is by default), no numeric jacobians will be calc...
 
virtual void init_context(DiffContext &context)
 
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
 
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
 
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...
 
const DofMap & get_dof_map() const
 
virtual void init_data()
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
 
const std::vector< std::vector< OutputShape > > & get_phi() const
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
This class provides all data required for a physics package (e.g.