Go to the documentation of this file.
    4 #include "libmesh/sparse_matrix.h" 
    5 #include "libmesh/numeric_vector.h" 
    6 #include "libmesh/dense_matrix.h" 
    7 #include "libmesh/dense_vector.h" 
    8 #include "libmesh/fe.h" 
    9 #include "libmesh/fe_interface.h" 
   10 #include "libmesh/fe_base.h" 
   11 #include "libmesh/elem_assembly.h" 
   12 #include "libmesh/quadrature_gauss.h" 
   15 #include "libmesh/transient_rb_theta_expansion.h" 
   16 #include "libmesh/transient_rb_assembly_expansion.h" 
   45     const unsigned int u_var = 0;
 
   47     FEBase * elem_fe = 
nullptr;
 
   50     const std::vector<Real> & JxW = elem_fe->get_JxW();
 
   52     const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
 
   60     for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
   61       for (
unsigned int i=0; i != n_u_dofs; i++)
 
   62         for (
unsigned int j=0; j != n_u_dofs; j++)
 
   72     const unsigned int u_var = 0;
 
   74     FEBase * elem_fe = 
nullptr;
 
   77     const std::vector<Real> & JxW = elem_fe->get_JxW();
 
   81     const std::vector<std::vector<RealGradient>> & dphi = elem_fe->get_dphi();
 
   89     for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
   90       for (
unsigned int i=0; i != n_u_dofs; i++)
 
   91         for (
unsigned int j=0; j != n_u_dofs; j++)
 
  102     const unsigned int u_var = 0;
 
  104     FEBase * elem_fe = 
nullptr;
 
  107     const std::vector<Real> & JxW = elem_fe->get_JxW();
 
  109     const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
 
  111     const std::vector<std::vector<RealGradient>> & dphi = elem_fe->get_dphi();
 
  119     for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  120       for (
unsigned int i=0; i != n_u_dofs; i++)
 
  121         for (
unsigned int j=0; j != n_u_dofs; j++)
 
  131     const unsigned int u_var = 0;
 
  133     FEBase * elem_fe = 
nullptr;
 
  136     const std::vector<Real> & JxW = elem_fe->get_JxW();
 
  138     const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
 
  140     const std::vector<std::vector<RealGradient>> & dphi = elem_fe->get_dphi();
 
  148     for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  149       for (
unsigned int i=0; i != n_u_dofs; i++)
 
  150         for (
unsigned int j=0; j != n_u_dofs; j++)
 
  160     const unsigned int u_var = 0;
 
  162     FEBase * elem_fe = 
nullptr;
 
  165     const std::vector<Real> & JxW = elem_fe->get_JxW();
 
  167     const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
 
  175     for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  176       for (
unsigned int i=0; i != n_u_dofs; i++)
 
  195     const unsigned int u_var = 0;
 
  197     FEBase * elem_fe = 
nullptr;
 
  200     const std::vector<Real> & JxW = elem_fe->get_JxW();
 
  202     const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
 
  213     if ((
min_x <= centroid(0)) && (centroid(0) <= 
max_x) &&
 
  214         (
min_y <= centroid(1)) && (centroid(1) <= 
max_y))
 
  215       for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  216         for (
unsigned int i=0; i != n_u_dofs; i++)
 
  264     L0(0.72, 0.88, 0.72, 0.88), 
 
  265     L1(0.12, 0.28, 0.72, 0.88),
 
  266     L2(0.12, 0.28, 0.12, 0.28),
 
  267     L3(0.72, 0.88, 0.12, 0.28)
 
  
void attach_F_assembly(ElemAssembly *Fq_assembly)
Attach ElemAssembly object for the right-hand side (both interior and boundary assembly).
 
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.
 
virtual void interior_assembly(FEMContext &c)
Perform the element interior assembly.
 
This class stores the set of RBTheta functor objects that define the "parameter-dependent expansion" ...
 
This class provides an encapsulated access to all static public member functions of finite element cl...
 
const std::vector< dof_id_type > & get_dof_indices() const
Accessor for element dof indices.
 
virtual void attach_output_theta(std::vector< std::unique_ptr< RBTheta >> &theta_q_l)
Attach a vector of pointers to functor objects that define one of the outputs.
 
unsigned int n_points() const
 
virtual Number evaluate(const RBParameters &mu)
Evaluate the functor object for the given parameter.
 
RealVectorValue RealGradient
 
const Elem & get_elem() const
Accessor for current Elem object.
 
This class is part of the rbOOmit framework.
 
CDRBAssemblyExpansion()
Constructor.
 
OutputAssembly(Real min_x_in, Real max_x_in, Real min_y_in, Real max_y_in)
 
virtual void attach_M_theta(RBTheta *theta_q_m)
Attach a pointer to a functor object that defines one of the theta_q_m terms.
 
virtual void interior_assembly(FEMContext &c)
Perform the element interior assembly.
 
virtual Point centroid() const
 
void attach_M_assembly(ElemAssembly *A_q_assembly)
Attach ElemAssembly object for the time-derivative (both interior and boundary assembly).
 
virtual void interior_assembly(FEMContext &c)
Perform the element interior assembly.
 
FEGenericBase< Real > FEBase
 
virtual void interior_assembly(FEMContext &c)
Perform the element interior assembly.
 
A Point defines a location in LIBMESH_DIM dimensional Real space.
 
virtual void attach_F_theta(RBTheta *theta_q_f)
Attach a pointer to a functor object that defines one of the theta_q_a terms.
 
void attach_A_assembly(ElemAssembly *Aq_assembly)
Attach ElemAssembly object for the left-hand side (both interior and boundary assembly).
 
virtual Number evaluate(const RBParameters &)
Evaluate the functor object for the given parameter.
 
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.
 
virtual void interior_assembly(FEMContext &c)
Perform the element interior assembly.
 
virtual void interior_assembly(FEMContext &c)
Perform the element interior assembly.
 
ElemAssembly provides a per-element (interior and boundary) assembly functionality.
 
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
 
This class is part of the rbOOmit framework.
 
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
 
CDRBThetaExpansion()
Constructor.
 
virtual void attach_output_assembly(std::vector< std::unique_ptr< ElemAssembly >> &output_assembly)
Attach ElemAssembly object for an output (both interior and boundary assembly).
 
virtual void attach_A_theta(RBTheta *theta_q_a)
Attach a pointer to a functor object that defines one of the theta_q_a terms.
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
virtual Number evaluate(const RBParameters &mu)
Evaluate the functor object for the given parameter.
 
This extends RBAssemblyExpansion to provide an assembly expansion for the case of time-dependent PDEs...
 
Real get_value(const std::string ¶m_name) const
Get the value of the specific parameter.
 
This class provides all data required for a physics package (e.g.