Go to the documentation of this file.
    6 #include <libmesh/dof_map.h> 
    7 #include <libmesh/elem.h> 
    8 #include <libmesh/equation_systems.h> 
    9 #include <libmesh/fe_base.h> 
   10 #include <libmesh/fe_interface.h> 
   11 #include <libmesh/mesh.h> 
   12 #include <libmesh/mesh_generation.h> 
   13 #include <libmesh/numeric_vector.h> 
   14 #include <libmesh/system.h> 
   21   CPPUNIT_TEST( testU );                        \ 
   22   CPPUNIT_TEST( testGradU );                    \ 
   23   CPPUNIT_TEST( testGradUComp ); 
   33   const Real & x = p(0);
 
   34   const Real & y = (LIBMESH_DIM > 1) ? p(1) : 0;
 
   35   const Real & z = (LIBMESH_DIM > 2) ? p(2) : 0;
 
   37   return x + 0.25*y + 0.0625*z;
 
   67   const Real & x = p(0);
 
   68   const Real & y = (LIBMESH_DIM > 1) ? p(1) : 0;
 
   69   const Real & z = (LIBMESH_DIM > 2) ? p(2) : 0;
 
   75   return (x + 0.25*y + 0.0625*z)/denom;
 
   84   const Real & x = p(0);
 
   85   const Real & y = (LIBMESH_DIM > 1) ? p(1) : 0;
 
   86   const Real & z = (LIBMESH_DIM > 2) ? p(2) : 0;
 
   95   const Real denom = xpoly * ypoly * zpoly;
 
   97   const Real numer = (x + 0.25*y + 0.0625*z);
 
   99   Gradient grad_n = 1, grad_d = xderiv * ypoly * zpoly;
 
  103       grad_d(1) = xpoly * yderiv * zpoly;
 
  108       grad_d(2) = xpoly * ypoly * zderiv;
 
  111   Gradient grad = (grad_n - numer * grad_d / denom) / denom;
 
  119 template <Order order, FEFamily family, ElemType elem_type>
 
  120 class FETest : 
public CppUnit::TestCase {
 
  123   unsigned int _dim, _nx, _ny, 
_nz;
 
  135     const std::unique_ptr<Elem> test_elem = 
Elem::build(elem_type);
 
  136     _dim = test_elem->dim();
 
  137     const unsigned int ny = _dim > 1;
 
  138     const unsigned int nz = _dim > 2;
 
  140     unsigned char weight_index = 0;
 
  146         weight_index = cast_int<unsigned char>
 
  148         libmesh_assert_not_equal_to(weight_index, 0);
 
  158                                        0., 1., 0., ny, 0., nz,
 
  166             const unsigned int nv = elem->n_vertices();
 
  167             const unsigned int nn = elem->n_nodes();
 
  170             const unsigned int n_edges =
 
  171               (elem->type() == 
EDGE3) ? 1 : elem->n_edges();
 
  172             const unsigned int n_faces =
 
  173               (elem->type() == 
QUAD9) ? 1 : elem->n_faces();
 
  174             const unsigned int nve = std::min(nv + n_edges, nn);
 
  175             const unsigned int nvef = std::min(nve + n_faces, nn);
 
  177             for (
unsigned int i = 0; i != nv; ++i)
 
  178               elem->node_ref(i).set_extra_datum<
Real>(weight_index, 1.);
 
  179             for (
unsigned int i = nv; i != nve; ++i)
 
  180               elem->node_ref(i).set_extra_datum<
Real>(weight_index, 
rational_w);
 
  182             for (
unsigned int i = nve; i != nvef; ++i)
 
  183               elem->node_ref(i).set_extra_datum<
Real>(weight_index, w2);
 
  185             for (
unsigned int i = nvef; i != nn; ++i)
 
  186               elem->node_ref(i).set_extra_datum<
Real>(weight_index, w3);
 
  221     _elem = rng.begin() == rng.end() ? nullptr : *(rng.begin());
 
  226     _ny = (_dim > 1) ? _nx : 1;
 
  227     _nz = (_dim > 2) ? _nx : 1;
 
  254 #ifdef LIBMESH_ENABLE_EXCEPTIONS 
  255     for (
unsigned int i=0; i != _nx; ++i)
 
  256       for (
unsigned int j=0; j != _ny; ++j)
 
  257         for (
unsigned int k=0; k != _nz; ++k)
 
  259             Real x = (
Real(i)/_nx), y = 0, z = 0;
 
  262               p(1) = y = (
Real(j)/_ny);
 
  264               p(2) = z = (
Real(k)/_nz);
 
  268             std::vector<Point> master_points
 
  271             _fe->
reinit(_elem, &master_points);
 
  274             for (std::size_t d = 0; d != _dof_indices.size(); ++d)
 
  278               LIBMESH_ASSERT_FP_EQUAL
 
  283               LIBMESH_ASSERT_FP_EQUAL
 
  308 #ifdef LIBMESH_ENABLE_EXCEPTIONS 
  309     for (
unsigned int i=0; i != _nx; ++i)
 
  310       for (
unsigned int j=0; j != _ny; ++j)
 
  311         for (
unsigned int k=0; k != _nz; ++k)
 
  321             std::vector<Point> master_points
 
  324             _fe->
reinit(_elem, &master_points);
 
  327             for (std::size_t d = 0; d != _dof_indices.size(); ++d)
 
  355                   LIBMESH_ASSERT_FP_EQUAL(
libmesh_real(grad_u(2)), 0.0625,
 
  379 #ifdef LIBMESH_ENABLE_EXCEPTIONS 
  380     for (
unsigned int i=0; i != _nx; ++i)
 
  381       for (
unsigned int j=0; j != _ny; ++j)
 
  382         for (
unsigned int k=0; k != _nz; ++k)
 
  392             std::vector<Point> master_points
 
  395             _fe->
reinit(_elem, &master_points);
 
  397             Number grad_u_x = 0, grad_u_y = 0, grad_u_z = 0;
 
  398             for (std::size_t d = 0; d != _dof_indices.size(); ++d)
 
  444 #define INSTANTIATE_FETEST(order, family, elemtype)                     \ 
  445   class FETest_##order##_##family##_##elemtype : public FETest<order, family, elemtype> { \ 
  447   CPPUNIT_TEST_SUITE( FETest_##order##_##family##_##elemtype );         \ 
  449   CPPUNIT_TEST_SUITE_END();                                             \ 
  452   CPPUNIT_TEST_SUITE_REGISTRATION( FETest_##order##_##family##_##elemtype ); 
  454 #endif // #ifdef __fe_test_h__ 
  
Manages consistently variables, degrees of freedom, and coefficient vectors.
 
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
 
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
 
virtual SimpleRange< element_iterator > active_local_element_ptr_range() override
 
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Fills the vector di with the global degree of freedom indices for the element.
 
virtual bool contains_point(const Point &p, Real tol=TOLERANCE) const
 
const std::vector< std::vector< OutputShape > > & get_dphidx() const
 
The libMesh namespace provides an interface to certain functionality in the library.
 
virtual SimpleRange< element_iterator > active_element_ptr_range() override
 
This class forms the foundation from which generic finite elements may be derived.
 
void set_default_mapping_data(const unsigned char data)
Set the default master space to physical space mapping basis functions to be used on newly added elem...
 
MetaPhysicL::DualNumber< T, D > sqrt(const MetaPhysicL::DualNumber< T, D > &in)
 
static const Real TOLERANCE
 
const std::vector< std::vector< OutputGradient > > & get_dphi() const
 
Gradient rational_test_grad(const Point &p, const Parameters &, const std::string &, const std::string &)
 
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
 
static const Real rational_w
 
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
 
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.
 
std::unique_ptr< NumericVector< Number > > current_local_solution
All the values I need to compute my contribution to the simulation at hand.
 
libMesh::Parallel::Communicator * TestCommWorld
 
virtual void init()
Initialize all the systems.
 
A Point defines a location in LIBMESH_DIM dimensional Real space.
 
Number rational_test(const Point &p, const Parameters &, const std::string &, const std::string &)
 
unsigned int add_node_integer(const std::string &name, bool allocate_data=true)
Register an integer datum (of type dof_id_type) to be added to each node in the mesh.
 
This is the EquationSystems class.
 
const FEType & variable_type(const unsigned int i) const
 
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
 
const std::vector< std::vector< OutputShape > > & get_dphidy() const
 
unsigned int add_node_datum(const std::string &name, bool allocate_data=true)
Register a datum (of type T) to be added to each node in the mesh.
 
This is the base class from which all geometric element types are derived.
 
Gradient linear_test_grad(const Point &, const Parameters &, const std::string &, const std::string &)
 
const std::vector< std::vector< OutputShape > > & get_dphidz() const
 
const DofMap & get_dof_map() const
 
const std::vector< std::vector< OutputShape > > & get_phi() const
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr)=0
This is at the core of this class.
 
void project_solution(FunctionBase< Number > *f, FunctionBase< Gradient > *g=nullptr) const
Projects arbitrary functions onto the current solution.
 
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
 
std::vector< dof_id_type > _dof_indices
 
Number linear_test(const Point &p, const Parameters &, const std::string &, const std::string &)
 
This class provides the ability to map between arbitrary, user-defined strings and several data types...
 
Parameters parameters
Data structure holding arbitrary parameters.