Go to the documentation of this file.
   18 #include "libmesh/getpot.h" 
   20 #include "curl_curl_system.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" 
   29 #include "libmesh/zero_function.h" 
   35                                const std::string & name_in,
 
   36                                const unsigned int number_in) :
 
   44   GetPot infile(
"vector_fe_ex3.in");
 
   63   FEMSystem::init_data();
 
   69   std::set<boundary_id_type> boundary_ids(all_ids, all_ids+4);
 
   71   std::vector<unsigned int> vars;
 
   72   vars.push_back(
u_var);
 
   81   FEMContext & c = cast_ref<FEMContext &>(context);
 
  107   FEMContext & c = cast_ref<FEMContext &>(context);
 
  117   const std::vector<Real> & JxW = fe->
get_JxW();
 
  120   const std::vector<std::vector<RealGradient>> & phi = fe->
get_phi();
 
  124   const std::vector<std::vector<RealGradient>> & curl_phi = fe->
get_curl_phi();
 
  126   const std::vector<Point> & qpoint = fe->
get_xyz();
 
  144   for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  157       for (
unsigned int i=0; i != n_u_dofs; i++)
 
  159           Fu(i) += (curl_u*curl_phi[i][qp] + u*phi[i][qp] - f*phi[i][qp])*JxW[qp];
 
  162           if (request_jacobian)
 
  163             for (
unsigned int j=0; j != n_u_dofs; j++)
 
  164               Kuu(i,j) += (curl_phi[j][qp]*curl_phi[i][qp] +
 
  165                            phi[j][qp]*phi[i][qp])*JxW[qp];
 
  169   return request_jacobian;
 
  176   FEMContext & c = cast_ref<FEMContext &>(context);
 
  186   const std::vector<Real> & JxW = side_fe->
get_JxW();
 
  189   const std::vector<std::vector<RealGradient>> & phi = side_fe->
get_phi();
 
  194   const std::vector<Point> & normals = side_fe->
get_normals();
 
  196   const std::vector<Point> & qpoint = side_fe->
get_xyz();
 
  200   const Real penalty = 1.e10;
 
  207   for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  216       Gradient Ncu = (u - u_exact).cross(N);
 
  218       for (
unsigned int i=0; i != n_u_dofs; i++)
 
  220           Fu(i) += penalty*Ncu*(phi[i][qp].cross(N))*JxW[qp];
 
  222           if (request_jacobian)
 
  223             for (
unsigned int j=0; j != n_u_dofs; j++)
 
  224               Kuu(i,j) += penalty*(phi[j][qp].cross(N))*(phi[i][qp].cross(N))*JxW[qp];
 
  228   return request_jacobian;
 
  
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
 
virtual void init_data()
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
 
RealGradient forcing(const Point &p)
 
const std::vector< Point > & get_normals() const
 
Defines a dense submatrix for use in Finite Element-type computations.
 
unsigned int n_points() const
 
The libMesh namespace provides an interface to certain functionality in the library.
 
RealGradient forcing(Real x, Real y)
 
This class forms the foundation from which generic finite elements may be derived.
 
int extra_quadrature_order
A member int that can be employed to indicate increased or reduced quadrature order.
 
const std::vector< Real > & get_JxW() const
 
CurlCurlSystem(EquationSystems &es, const std::string &name, const unsigned int number)
 
CurlCurlExactSolution soln
 
void interior_curl(unsigned int var, unsigned int qp, OutputType &curl_u) const
 
This class provides a specific system class.
 
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.
 
bool print_jacobians
Set print_jacobians to true to print J whenever it is assembled.
 
A Point defines a location in LIBMESH_DIM dimensional Real space.
 
unsigned int n_dof_indices() const
Total number of dof indices on the element.
 
ConstFunction that simply returns 0.
 
This class provides all data required for a physics package (e.g.
 
virtual void init_context(DiffContext &context)
 
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.
 
Real verify_analytic_jacobians
If verify_analytic_jacobian is equal to zero (as it is by default), no numeric jacobians will be calc...
 
const std::vector< std::vector< OutputShape > > & get_curl_phi() const
 
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
 
RealGradient exact_solution(const Point &p)
 
Number side_value(unsigned int var, unsigned int qp) const
 
virtual bool element_time_derivative(bool request_jacobian, DiffContext &context)
Adds the time derivative contribution on elem to elem_residual.
 
const QBase & get_side_qrule() const
Accessor for element side quadrature rule for the dimension of the current _elem.
 
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 std::vector< std::vector< OutputShape > > & get_phi() const
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
void init_dirichlet_bcs()
 
virtual bool side_time_derivative(bool request_jacobian, DiffContext &context)
Adds the time derivative contribution on side of elem to elem_residual.
 
Number interior_value(unsigned int var, unsigned int qp) const
 
This class provides all data required for a physics package (e.g.