Go to the documentation of this file.
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"
64 int main (
int argc,
char ** argv)
74 #if !defined(LIBMESH_HAVE_TRIANGLE) || !defined(LIBMESH_HAVE_TETGEN) || !defined(LIBMESH_ENABLE_AMR)
75 libmesh_example_requires(
false,
"--disable-strict-lgpl --enable-amr");
79 libmesh_example_requires(2 <= LIBMESH_DIM,
"2D support");
90 if (argc == 3 && std::atoi(argv[2]) == 3)
96 mesh.
read ((
dim==2) ?
"mesh.xda" :
"hybrid_3d.xda");
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();
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"
void integrate_function(const MeshBase &mesh)
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
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.
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...
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
The libMesh namespace provides an interface to certain functionality in the library.
const Parallel::Communicator & comm() const
unsigned int mesh_dimension() const
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
Implements (adaptive) mesh refinement algorithms for a MeshBase.
This is the MeshBase class.
void libmesh_ignore(const Args &...)
Real integrand(const Point &p)
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
int main(int argc, char **argv)
A Point defines a location in LIBMESH_DIM dimensional Real space.
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Real distance(const Point &p)
This class implements generic composite quadrature rules.
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
MetaPhysicL::DualNumber< T, D > norm(const MetaPhysicL::DualNumber< T, D > &in)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void uniformly_refine(unsigned int n=1)
Uniformly refines the mesh n times.