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"
13 #include "libmesh/boundary_info.h"
16 #include "libmesh/rb_theta.h"
17 #include "libmesh/rb_assembly_expansion.h"
18 #include "libmesh/rb_parametrized_function.h"
19 #include "libmesh/rb_eim_construction.h"
20 #include "libmesh/rb_eim_theta.h"
59 return 1. + curvature*p(0);
71 return 1. + curvature*p(0);
83 return 1./(1. + curvature*p(0));
99 std::vector<boundary_id_type> bc_ids;
101 for (std::vector<boundary_id_type>::const_iterator b =
102 bc_ids.begin(); b != bc_ids.end(); ++b)
103 if (*b == 1 || *b == 2 || *b == 3 || *b == 4)
105 const unsigned int u_var = 0;
107 FEBase * side_fe =
nullptr;
110 const std::vector<Real> & JxW_side = side_fe->get_JxW();
112 const std::vector<std::vector<Real>> & phi_side = side_fe->get_phi();
120 for (
unsigned int qp=0; qp != n_sidepoints; qp++)
121 for (
unsigned int i=0; i != n_u_dofs; i++)
122 for (
unsigned int j=0; j != n_u_dofs; j++)
142 std::vector<boundary_id_type> bc_ids;
144 for (std::vector<boundary_id_type>::const_iterator b =
145 bc_ids.begin(); b != bc_ids.end(); ++b)
146 if (*b == 1 || *b == 3)
148 const unsigned int u_var = 0;
150 FEBase * side_fe =
nullptr;
153 const std::vector<Real> & JxW_side = side_fe->get_JxW();
155 const std::vector<std::vector<Real>> & phi_side = side_fe->get_phi();
157 const std::vector<Point> & xyz = side_fe->get_xyz();
165 for (
unsigned int qp=0; qp != n_sidepoints; qp++)
167 Real x_hat = xyz[qp](0);
169 for (
unsigned int i=0; i != n_u_dofs; i++)
170 for (
unsigned int j=0; j != n_u_dofs; j++)
171 c.
get_elem_jacobian()(i,j) += JxW_side[qp] * x_hat * phi_side[j][qp]*phi_side[i][qp];
189 std::vector<boundary_id_type> bc_ids;
191 for (std::vector<boundary_id_type>::const_iterator b =
192 bc_ids.begin(); b != bc_ids.end(); ++b)
193 if (*b == 2 || *b == 4)
195 const unsigned int u_var = 0;
197 FEBase * side_fe =
nullptr;
200 const std::vector<Real> & JxW_side = side_fe->get_JxW();
202 const std::vector<std::vector<Real>> & phi_side = side_fe->get_phi();
212 for (
unsigned int qp=0; qp != n_sidepoints; qp++)
214 for (
unsigned int i=0; i != n_u_dofs; i++)
215 for (
unsigned int j=0; j != n_u_dofs; j++)
222 for (
unsigned int qp=0; qp != n_sidepoints; qp++)
224 for (
unsigned int i=0; i != n_u_dofs; i++)
225 for (
unsigned int j=0; j != n_u_dofs; j++)
241 return mu.
get_value(
"kappa") * RBEIMTheta::evaluate(mu);
248 unsigned int basis_function_index_in) :
250 basis_function_index_in)
256 const unsigned int u_var = 0;
259 const unsigned int Gx_var = 0;
260 const unsigned int Gy_var = 1;
261 const unsigned int Gz_var = 2;
263 FEBase * elem_fe =
nullptr;
266 const std::vector<Real> & JxW = elem_fe->get_JxW();
268 const std::vector<std::vector<RealGradient>> & dphi = elem_fe->get_dphi();
273 std::vector<Number> eim_values_Gx;
279 std::vector<Number> eim_values_Gy;
285 std::vector<Number> eim_values_Gz;
292 for (
unsigned int i=0; i != n_u_dofs; i++)
293 for (
unsigned int j=0; j != n_u_dofs; j++)
294 c.
get_elem_jacobian()(i,j) += JxW[qp] * (eim_values_Gx[qp]*dphi[i][qp](0)*dphi[j][qp](0) +
295 eim_values_Gy[qp]*dphi[i][qp](1)*dphi[j][qp](1) +
296 eim_values_Gz[qp]*dphi[i][qp](2)*dphi[j][qp](2));
311 const unsigned int u_var = 0;
313 FEBase * elem_fe =
nullptr;
316 const std::vector<Real> & JxW = elem_fe->get_JxW();
318 const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
326 for (
unsigned int qp=0; qp != n_qpoints; qp++)
327 for (
unsigned int i=0; i != n_u_dofs; i++)
345 const unsigned int u_var = 0;
347 FEBase * elem_fe =
nullptr;
350 const std::vector<Real> & JxW = elem_fe->get_JxW();
352 const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
354 const std::vector<Point> & xyz = elem_fe->get_xyz();
362 for (
unsigned int qp=0; qp != n_qpoints; qp++)
364 Real x_hat = xyz[qp](0);
366 for (
unsigned int i=0; i != n_u_dofs; i++)
377 const unsigned int u_var = 0;
379 FEBase * elem_fe =
nullptr;
382 const std::vector<Real> & JxW = elem_fe->get_JxW();
384 const std::vector<std::vector<RealGradient>> & dphi = elem_fe->get_dphi();
392 for (
unsigned int qp=0; qp != n_qpoints; qp++)
393 for (
unsigned int i=0; i != n_u_dofs; i++)
394 for (
unsigned int j=0; j != n_u_dofs; j++)
404 FEBase * elem_fe =
nullptr;
407 const std::vector<Real> & JxW = elem_fe->get_JxW();
409 const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
415 for (
unsigned int qp=0; qp != n_qpoints; qp++)
416 for (
unsigned int i=0; i != n_dofs; i++)
417 for (
unsigned int j=0; j != n_dofs; j++)
This class provides functionality required to define an RBTheta object that arises from an "Empirical...
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.
std::vector< boundary_id_type > boundary_ids(const Node *node) const
unsigned char side
Current side for side_* to examine.
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
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.
Defines a dense submatrix for use in Finite Element-type computations.
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.
virtual Number evaluate(const RBParameters &)
Evaluate the functor object for the given parameter.
virtual Number evaluate(const RBParameters &mu)
Evaluate the functor object for the given parameter.
Ex6AssemblyExpansion(RBConstruction &rb_con)
Constructor.
This class provides functionality required to define an assembly object that arises from an "Empirica...
virtual void interior_assembly(FEMContext &c)
Perform the element interior assembly.
virtual Number evaluate(const RBParameters &mu, const Point &p, const Elem &)
Evaluate this parametrized function for the parameter value mu at the point p.
This class stores the set of ElemAssembly functor objects that define the "parameter-independent expa...
virtual void boundary_assembly(FEMContext &c)
Perform the element boundary assembly.
FEGenericBase< Real > FEBase
virtual Number evaluate(const RBParameters &mu)
Evaluate the functor object for the given parameter.
virtual void interior_assembly(FEMContext &c)
Perform the element interior assembly.
const MeshBase & get_mesh() const
const std::vector< Point > & get_points() const
A Point defines a location in LIBMESH_DIM dimensional Real space.
A simple functor class that provides a RBParameter-dependent function.
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.
virtual Number evaluate(const RBParameters &mu, const Point &p, const Elem &)
Evaluate this parametrized function for the parameter value mu at the point p.
void attach_A_assembly(ElemAssembly *Aq_assembly)
Attach ElemAssembly object for the left-hand side (both interior and boundary assembly).
This class is part of the rbOOmit framework.
virtual void evaluate_basis_function(unsigned int var, const Elem &element, const std::vector< Point > &points, std::vector< Number > &values)
Evaluate variable var_number of this object's EIM basis function at the points points,...
This class is part of the rbOOmit framework.
virtual void interior_assembly(FEMContext &c)
Perform the element interior assembly.
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.
ElemAssembly provides a per-element (interior and boundary) assembly functionality.
virtual void interior_assembly(FEMContext &c)
Perform the element interior assembly.
virtual void boundary_assembly(FEMContext &c)
Perform the element boundary assembly.
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
This class is part of the rbOOmit framework.
This is the base class from which all geometric element types are derived.
const QBase & get_side_qrule() const
Accessor for element side quadrature rule for the dimension of the current _elem.
virtual Number evaluate(const RBParameters &mu)
Evaluate this RBEIMTheta object at the parameter mu.
virtual Number evaluate(const RBParameters &mu, const Point &p, const Elem &)
Evaluate this parametrized function for the parameter value mu at the point p.
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...
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
ThetaEIM(RBEIMEvaluation &rb_eim_eval_in, unsigned int index_in)
This class stores the set of RBTheta functor objects that define the "parameter-dependent expansion" ...
virtual Number evaluate(const RBParameters &mu)
Evaluate the functor object for the given parameter.
AssemblyEIM(RBEIMConstruction &rb_eim_con_in, unsigned int basis_function_index_in)
virtual void interior_assembly(FEMContext &c)
Perform the element interior assembly.
Ex6ThetaExpansion()
Constructor.
virtual void boundary_assembly(FEMContext &c)
Perform the element boundary assembly.
This class is part of the rbOOmit framework.
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.