Go to the documentation of this file.
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();
140 const std::vector<Real> *
JxW;
168 #endif // #ifdef LIBMESH_ENABLE_AMR
170 #endif // LIBMESH_HP_COARSENTEST_H
Manages consistently variables, degrees of freedom, and coefficient vectors.
std::unique_ptr< FEBase > fe_coarse
RealVectorValue RealGradient
The libMesh namespace provides an interface to certain functionality in the library.
RealTensorValue RealTensor
const std::vector< Real > * JxW
Mapping jacobians.
const std::vector< Point > * xyz_values
Quadrature locations.
VectorValue< Real > RealVectorValue
Useful typedefs to allow transparent switching between Real and Complex data types.
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
std::unique_ptr< QBase > qrule
The quadrature rule for the fine element.
DenseMatrix< Number > Ke
Linear system for projections.
Elem * coarse
The coarse element on which a solution projection is cached.
TensorValue< Real > RealTensorValue
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.
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
virtual ~HPCoarsenTest()=default
HPCoarsenTest()
Constructor.
DenseVector< Number > Uc
Coefficients for projected coarse and projected p-derefined solutions.
void add_projection(const System &, const Elem *, unsigned int var)
The helper function which adds individual fine element data to the coarse element projection.
std::unique_ptr< FEBase > fe
The finite element objects for fine and coarse elements.
std::vector< dof_id_type > dof_indices
Global DOF indices for fine elements.
const std::vector< std::vector< RealTensor > > * d2phi
Real p_weight
Because the coarsening test seems to always choose p refinement, we're providing an option to make h ...
const std::vector< std::vector< Real > > * phi_coarse
This is the base class from which all geometric element types are derived.
const std::vector< std::vector< RealGradient > > * dphi
const std::vector< std::vector< Real > > * phi
The shape functions and their derivatives.
const std::vector< std::vector< RealGradient > > * dphi_coarse
This class uses the error estimate given by different types of derefinement (h coarsening or p reduct...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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...
HPCoarsenTest & operator=(const HPCoarsenTest &)=delete
std::vector< Point > coarse_qpoints
const std::vector< std::vector< RealTensor > > * d2phi_coarse