Go to the documentation of this file.
26 #include "libmesh/libmesh.h"
28 #include "libmesh/boundary_info.h"
29 #include "libmesh/bounding_box.h"
30 #include "libmesh/elem.h"
31 #include "libmesh/fe.h"
32 #include "libmesh/getpot.h"
33 #include "libmesh/mesh.h"
34 #include "libmesh/point.h"
35 #include "libmesh/quadrature_gauss.h"
42 <<
" --dim d --input inputmesh --output outputmesh --newbcid idnum --tests --moretests"
48 int main(
int argc,
char ** argv)
52 GetPot cl(argc, argv);
54 unsigned char dim = -1;
55 if (!cl.search(
"--dim"))
64 if (!cl.search(
"--input"))
66 libMesh::err <<
"No --input argument found!" << std::endl;
69 const char * meshname = cl.next(
"mesh.xda");
74 if (!cl.search(
"--newbcid"))
76 libMesh::err <<
"No --bcid argument found!" << std::endl;
82 Point minpt(-std::numeric_limits<Real>::max());
84 minpt(1) = -std::numeric_limits<Real>::max();
87 minpt(2) = -std::numeric_limits<Real>::max();
94 if (cl.search(
"--minnormalx"))
95 normals.min()(0) = cl.next(normals.min()(0));
96 if (cl.search(
"--maxnormalx"))
97 normals.max()(0) = cl.next(normals.max()(0));
99 if (cl.search(
"--minpointx"))
100 points.
min()(0) = cl.next(points.
min()(0));
101 if (cl.search(
"--maxpointx"))
102 points.
max()(0) = cl.next(points.
max()(0));
105 if (cl.search(
"--minnormaly"))
106 normals.min()(1) = cl.next(normals.min()(1));
107 if (cl.search(
"--maxnormaly"))
108 normals.max()(1) = cl.next(normals.max()(1));
110 if (cl.search(
"--minpointy"))
111 points.
min()(1) = cl.next(points.
min()(1));
112 if (cl.search(
"--maxpointy"))
113 points.
max()(1) = cl.next(points.
max()(1));
117 if (cl.search(
"--minnormalz"))
118 normals.min()(2) = cl.next(normals.min()(2));
119 if (cl.search(
"--maxnormalz"))
120 normals.max()(2) = cl.next(normals.max()(2));
122 if (cl.search(
"--minpointz"))
123 points.
min()(2) = cl.next(points.
min()(2));
124 if (cl.search(
"--maxpointz"))
125 points.
max()(2) = cl.next(points.
max()(2));
130 libMesh::out <<
"min normal = " << normals.min() << std::endl;
131 libMesh::out <<
"max normal = " << normals.max() << std::endl;
133 bool matcholdbcid =
false;
135 if (cl.search(
"--oldbcid"))
138 oldbcid = cl.next(oldbcid);
145 fe->attach_quadrature_rule(&qface);
146 const std::vector<Point> & face_points = fe->get_xyz();
147 const std::vector<Point> & face_normals = fe->get_normals();
151 unsigned int n_sides = elem->n_sides();
154 std::vector<boundary_id_type> ids;
156 for (
unsigned short s=0; s != n_sides; ++s)
158 if (elem->neighbor_ptr(s))
162 const Point & p = face_points[0];
163 const Point & n = face_normals[0];
171 normals.contains_point(n))
183 if (matcholdbcid && b_id != oldbcid)
199 std::string outputname;
200 if (cl.search(
"--output"))
202 outputname = cl.next(
"mesh.xda");
207 outputname += meshname;
212 libMesh::out <<
"Wrote mesh " << outputname << std::endl;
std::vector< boundary_id_type > boundary_ids(const Node *node) const
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.
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
This class implements specific orders of Gauss quadrature.
virtual void write(const std::string &name)=0
Defines a Cartesian bounding box by the two corner extremum.
The libMesh namespace provides an interface to certain functionality in the library.
void usage_error(const char *progname)
int main(int argc, char **argv)
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
virtual SimpleRange< element_iterator > element_ptr_range()=0
const Point & max() const
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
A Point defines a location in LIBMESH_DIM dimensional Real space.
void remove_side(const Elem *elem, const unsigned short int side)
Removes all boundary conditions associated with side side of element elem, if any exist.
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
void regenerate_id_sets()
Clears and regenerates the cached sets of ids.
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
bool contains_point(const Point &) const
static const boundary_id_type invalid_id
Number used for internal use.
void add_side(const dof_id_type elem, const unsigned short int side, const boundary_id_type id)
Add side side of element number elem with boundary id id to the boundary information data structure.
const Point & min() const