libMesh
Functions
meshbcid.C File Reference

Go to the source code of this file.

Functions

void usage_error (const char *progname)
 
int main (int argc, char **argv)
 

Function Documentation

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 48 of file meshbcid.C.

49 {
50  LibMeshInit init(argc, argv);
51 
52  GetPot cl(argc, argv);
53 
54  unsigned char dim = -1;
55  if (!cl.search("--dim"))
56  {
57  libMesh::err << "No --dim argument found!" << std::endl;
58  usage_error(argv[0]);
59  }
60  dim = cl.next(dim);
61 
62  Mesh mesh(init.comm(), dim);
63 
64  if (!cl.search("--input"))
65  {
66  libMesh::err << "No --input argument found!" << std::endl;
67  usage_error(argv[0]);
68  }
69  const char * meshname = cl.next("mesh.xda");
70 
71  mesh.read(meshname);
72  libMesh::out << "Loaded mesh " << meshname << std::endl;
73 
74  if (!cl.search("--newbcid"))
75  {
76  libMesh::err << "No --bcid argument found!" << std::endl;
77  usage_error(argv[0]);
78  }
79  boundary_id_type bcid = 0;
80  bcid = cl.next(bcid);
81 
82  Point minpt(-std::numeric_limits<Real>::max());
83 #if LIBMESH_DIM > 1
84  minpt(1) = -std::numeric_limits<Real>::max();
85 #endif
86 #if LIBMESH_DIM > 2
87  minpt(2) = -std::numeric_limits<Real>::max();
88 #endif
89  Point maxpt = -minpt;
90 
91  BoundingBox normals(minpt, maxpt),
92  points(minpt, maxpt);
93 
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));
98 
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));
103 
104 #if LIBMESH_DIM > 1
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));
109 
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));
114 #endif
115 
116 #if LIBMESH_DIM > 2
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));
121 
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));
126 #endif
127 
128  libMesh::out << "min point = " << points.min() << std::endl;
129  libMesh::out << "max point = " << points.max() << std::endl;
130  libMesh::out << "min normal = " << normals.min() << std::endl;
131  libMesh::out << "max normal = " << normals.max() << std::endl;
132 
133  bool matcholdbcid = false;
134  boundary_id_type oldbcid = 0;
135  if (cl.search("--oldbcid"))
136  {
137  matcholdbcid = true;
138  oldbcid = cl.next(oldbcid);
139  if (oldbcid < 0)
140  oldbcid = BoundaryInfo::invalid_id;
141  }
142 
143  std::unique_ptr<FEBase> fe = FEBase::build(dim, FEType(FIRST,LAGRANGE));
144  QGauss qface(dim-1, CONSTANT);
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();
148 
149  for (auto & elem : mesh.element_ptr_range())
150  {
151  unsigned int n_sides = elem->n_sides();
152 
153  // Container to catch ids handed back from BoundaryInfo
154  std::vector<boundary_id_type> ids;
155 
156  for (unsigned short s=0; s != n_sides; ++s)
157  {
158  if (elem->neighbor_ptr(s))
159  continue;
160 
161  fe->reinit(elem,s);
162  const Point & p = face_points[0];
163  const Point & n = face_normals[0];
164 
165  //libMesh::out << "elem = " << elem->id() << std::endl;
166  //libMesh::out << "centroid = " << elem->centroid() << std::endl;
167  //libMesh::out << "p = " << p << std::endl;
168  //libMesh::out << "n = " << n << std::endl;
169 
170  if (points.contains_point(p) &&
171  normals.contains_point(n))
172  {
173  // Get the list of boundary ids for this side
174  mesh.get_boundary_info().boundary_ids(elem, s, ids);
175 
176  // There should be at most one value present, otherwise the
177  // logic here won't work.
178  libmesh_assert(ids.size() <= 1);
179 
180  // A convenient name for the side's ID.
181  boundary_id_type b_id = ids.empty() ? BoundaryInfo::invalid_id : ids[0];
182 
183  if (matcholdbcid && b_id != oldbcid)
184  continue;
185 
187  mesh.get_boundary_info().add_side(elem, s, bcid);
188  //libMesh::out << "Set element " << elem->id() << " side " << s <<
189  // " to boundary " << bcid << std::endl;
190  }
191  }
192  }
193 
194  // We might have removed *every* instance of a given id, and if that
195  // happened then we should make sure that file formats which write
196  // out id sets do not write out the removed id.
198 
199  std::string outputname;
200  if (cl.search("--output"))
201  {
202  outputname = cl.next("mesh.xda");
203  }
204  else
205  {
206  outputname = "new.";
207  outputname += meshname;
208  }
209 
210 
211  mesh.write(outputname.c_str());
212  libMesh::out << "Wrote mesh " << outputname << std::endl;
213 
214  return 0;
215 }

References libMesh::BoundaryInfo::add_side(), libMesh::BoundaryInfo::boundary_ids(), libMesh::FEGenericBase< OutputType >::build(), libMesh::CONSTANT, libMesh::BoundingBox::contains_point(), dim, libMesh::MeshBase::element_ptr_range(), libMesh::err, libMesh::FIRST, libMesh::MeshBase::get_boundary_info(), libMesh::TriangleWrapper::init(), libMesh::BoundaryInfo::invalid_id, libMesh::LAGRANGE, libMesh::libmesh_assert(), libMesh::BoundingBox::max(), mesh, libMesh::BoundingBox::min(), libMesh::out, libMesh::MeshBase::read(), libMesh::BoundaryInfo::regenerate_id_sets(), libMesh::BoundaryInfo::remove_side(), usage_error(), and libMesh::MeshBase::write().

◆ usage_error()

void usage_error ( const char *  progname)

Definition at line 39 of file meshbcid.C.

40 {
41  libMesh::out << "Usage: " << progname
42  << " --dim d --input inputmesh --output outputmesh --newbcid idnum --tests --moretests"
43  << std::endl;
44 
45  exit(1);
46 }

References libMesh::out.

Referenced by main().

libMesh::BoundaryInfo::boundary_ids
std::vector< boundary_id_type > boundary_ids(const Node *node) const
Definition: boundary_info.C:985
libMesh::Mesh
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
libMesh::MeshBase::read
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.
libMesh::MeshBase::get_boundary_info
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:132
libMesh::QGauss
This class implements specific orders of Gauss quadrature.
Definition: quadrature_gauss.h:39
libMesh::MeshBase::write
virtual void write(const std::string &name)=0
libMesh::BoundingBox
Defines a Cartesian bounding box by the two corner extremum.
Definition: bounding_box.h:40
usage_error
void usage_error(const char *progname)
Definition: meshbcid.C:39
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::boundary_id_type
int8_t boundary_id_type
Definition: id_types.h:51
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
libMesh::MeshBase::element_ptr_range
virtual SimpleRange< element_iterator > element_ptr_range()=0
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::CONSTANT
Definition: enum_order.h:41
libMesh::BoundaryInfo::remove_side
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.
Definition: boundary_info.C:1462
libMesh::LibMeshInit
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:83
libMesh::BoundaryInfo::regenerate_id_sets
void regenerate_id_sets()
Clears and regenerates the cached sets of ids.
Definition: boundary_info.C:159
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::err
OStreamProxy err
libMesh::LAGRANGE
Definition: enum_fe_family.h:36
libMesh::out
OStreamProxy out
libMesh::BoundaryInfo::add_side
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.
Definition: boundary_info.C:886
libMesh::FIRST
Definition: enum_order.h:42