31 #include "libmesh/libmesh.h"    32 #include "libmesh/mesh.h"    33 #include "libmesh/mesh_generation.h"    34 #include "libmesh/mesh_refinement.h"    35 #include "libmesh/quadrature_gauss.h"    36 #include "libmesh/quadrature_composite.h"    37 #include "libmesh/fe.h"    38 #include "libmesh/elem.h"    39 #include "libmesh/parallel.h"    40 #include "libmesh/getpot.h"    65 int main (
int argc, 
char ** argv)
    75 #if !defined(LIBMESH_HAVE_TRIANGLE) || !defined(LIBMESH_HAVE_TETGEN) || !defined(LIBMESH_ENABLE_AMR)    76   libmesh_example_requires(
false, 
"--disable-strict-lgpl --enable-amr");
    80   libmesh_example_requires(2 <= LIBMESH_DIM, 
"2D support");
    89     const unsigned int dim =
    93       std::cout << 
"Running in 3D" << std::endl;
    95     mesh.
read ((
dim==2) ? 
"mesh.xda" : 
"hybrid_3d.xda");
   103   const int n_refinements =
   123 #if defined(LIBMESH_HAVE_TRIANGLE) && defined(LIBMESH_HAVE_TETGEN)   125   std::vector<Real> vertex_distance;
   134   const std::vector<Point> & q_points = fe->get_xyz();
   135   const std::vector<Real>  & JxW      = fe->get_JxW();
   137   for (
const auto & elem : 
mesh.active_local_element_ptr_range())
   139       vertex_distance.clear();
   141       for (
unsigned int v=0; v<elem->n_vertices(); v++)
   142         vertex_distance.push_back (
distance(elem->point(v)));
   144       qrule.
init (*elem, vertex_distance);
   147                   &(qrule.get_points()),
   148                   &(qrule.get_weights()));
   155       for (std::size_t qp=0; qp<q_points.size(); qp++)
   156         int_val += JxW[qp] * 
integrand(q_points[qp]);
   161   libMesh::out << 
"\n***********************************\n"   162                << 
" int_val   = " << int_val << std::endl
   164                << 
"\n***********************************\n" class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
T command_line_next(std::string name, T default_value)
Use GetPot's search()/next() functions to get following arguments from the command line...
virtual void read(const std::string &name, void *mesh_data=nullptr, bool skip_renumber_nodes_and_elements=false, bool skip_find_neighbors=false)=0
Interfaces for reading/writing a mesh to/from a file. 
const Parallel::Communicator & comm() const
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g. 
virtual void init(const Elem &elem, const std::vector< Real > &vertex_distance_func, unsigned int p_level=0) override
Overrides the base class init() function, and uses the ElemCutter to subdivide the element into "insi...
The libMesh namespace provides an interface to certain functionality in the library. 
Real distance(const Point &p)
Real integrand(const Point &p)
This is the MeshBase class. 
void libmesh_ignore(const Args &...)
void integrate_function(const MeshBase &mesh)
Implements (adaptive) mesh refinement algorithms for a MeshBase. 
void print_info(std::ostream &os=libMesh::out, const unsigned int verbosity=0, const bool global=true) const
Prints relevant information about the mesh. 
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary. 
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type. 
This class implements generic composite quadrature rules. 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
unsigned int mesh_dimension() const
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default. 
int main(int argc, char **argv)
A Point defines a location in LIBMESH_DIM dimensional Real space. 
void uniformly_refine(unsigned int n=1)
Uniformly refines the mesh n times.