20 #ifndef LIBMESH_HP_COARSENTEST_H 21 #define LIBMESH_HP_COARSENTEST_H 24 #include "libmesh/dense_matrix.h" 25 #include "libmesh/dense_vector.h" 26 #include "libmesh/hp_selector.h" 27 #include "libmesh/id_types.h" 28 #include "libmesh/libmesh_common.h" 29 #include "libmesh/fe.h" 30 #include "libmesh/quadrature.h" 36 #ifdef LIBMESH_ENABLE_AMR 45 template <
typename T>
class TensorValue;
76 libmesh_experimental();
150 const std::vector<Real> *
JxW;
183 #endif // #ifdef LIBMESH_ENABLE_AMR 185 #endif // LIBMESH_HP_COARSENTEST_H RealVectorValue RealGradient
const std::vector< std::vector< Real > > * phi
The shape functions and their derivatives.
void add_projection(const System &, const Elem *, unsigned int var)
The helper function which adds individual fine element data to the coarse element projection...
HPCoarsenTest & operator=(const HPCoarsenTest &)=delete
virtual void select_refinement(System &system) override
This pure virtual function must be redefined in derived classes to take a mesh flagged for h refineme...
VectorValue< Real > RealVectorValue
Useful typedefs to allow transparent switching between Real and Complex data types.
Subclasses of this abstract base class choose between h refining and p elevation. ...
std::unique_ptr< FEBase > fe_coarse
const std::vector< std::vector< Real > > * phi_coarse
int _extra_order
Extra order to use for quadrature rule.
This is the base class from which all geometric element types are derived.
const std::vector< std::vector< RealGradient > > * dphi_coarse
RealTensorValue RealTensor
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
The libMesh namespace provides an interface to certain functionality in the library.
const std::vector< Point > * xyz_values
Quadrature locations.
DenseVector< Number > Uc
Coefficients for projected coarse and projected p-derefined solutions.
const std::vector< Real > * JxW
Mapping jacobians.
TensorValue< Real > RealTensorValue
Useful typedefs to allow transparent switching between Real and Complex data types.
virtual ~HPCoarsenTest()=default
HPCoarsenTest()
Constructor.
Manages consistently variables, degrees of freedom, and coefficient vectors.
This class uses the error estimate given by different types of derefinement (h coarsening or p reduct...
std::unique_ptr< FEBase > fe
The finite element objects for fine and coarse elements.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const std::vector< std::vector< RealTensor > > * d2phi_coarse
const std::vector< std::vector< RealTensor > > * d2phi
const std::vector< std::vector< RealGradient > > * dphi
std::vector< dof_id_type > dof_indices
Global DOF indices for fine elements.
Real p_weight
Because the coarsening test seems to always choose p refinement, we're providing an option to make h ...
std::vector< Point > coarse_qpoints
Elem * coarse
The coarse element on which a solution projection is cached.
std::unique_ptr< QBase > qrule
The quadrature rule for the fine element.
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
void extra_quadrature_order(const int extraorder)
Increases or decreases the order of the quadrature rule used for numerical integration.
DenseMatrix< Number > Ke
Linear system for projections.