NEMoSys  0.63.0
A modular, extensible resource with robust automated mesh generation, mesh quality analysis, adaptive mesh refinement, and data transfer between arbitrary meshes.
meshSrch.H
Go to the documentation of this file.
1 /*******************************************************************************
2 * Promesh *
3 * Copyright (C) 2022, IllinoisRocstar LLC. All rights reserved. *
4 * *
5 * Promesh is the property of IllinoisRocstar LLC. *
6 * *
7 * IllinoisRocstar LLC *
8 * Champaign, IL *
9 * www.illinoisrocstar.com *
10 * promesh@illinoisrocstar.com *
11 *******************************************************************************/
12 /*******************************************************************************
13 * This file is part of Promesh *
14 * *
15 * This version of Promesh is free software: you can redistribute it and/or *
16 * modify it under the terms of the GNU Lesser General Public License as *
17 * published by the Free Software Foundation, either version 3 of the License, *
18 * or (at your option) any later version. *
19 * *
20 * Promesh is distributed in the hope that it will be useful, but WITHOUT ANY *
21 * WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS *
22 * FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more *
23 * details. *
24 * *
25 * You should have received a copy of the GNU Lesser General Public License *
26 * along with this program. If not, see <https://www.gnu.org/licenses/>. *
27 * *
28 *******************************************************************************/
29 #ifndef NEMOSYS_MESHSRCH_H_
30 #define NEMOSYS_MESHSRCH_H_
31 
32 #include <set>
33 
34 #include <vtkCellLocator.h>
35 
36 #include "nemosys_export.h"
37 #include "Mesh/meshBase.H"
38 
39 class NEMOSYS_EXPORT meshSrch : public meshBase {
40  // constructors and destructors
41  public:
42  meshSrch() = delete;
43  explicit meshSrch(meshBase *mb) : meshBase(*mb), upd_vcl(true) {}
44 
45  static meshSrch *Create(meshBase *mb) {
46  auto *ms = new meshSrch(mb);
47  return ms;
48  }
49 
50  static std::shared_ptr<meshSrch> CreateShared(meshBase *mb) {
51  std::shared_ptr<meshSrch> ms;
52  ms.reset(meshSrch::Create(mb));
53  return ms;
54  }
55 
56  static std::unique_ptr<meshSrch> CreateUnique(meshBase *mb) {
57  return std::unique_ptr<meshSrch>(meshSrch::Create(mb));
58  }
59 
60  ~meshSrch() override = default;
61 
62  // general access
63  public:
64  std::vector<std::vector<double>> getCellVec(nemId_t id) const override;
65  std::vector<double> getPoint(nemId_t id) const override;
66  std::vector<double> getCellCenter(nemId_t cellID) const override;
67 
68  int getCellType() const override { return 0; }
69  std::vector<std::vector<double>> getVertCrds() const override {
70  return std::vector<std::vector<double>>();
71  }
72  std::map<nemId_t, std::vector<double>> getCell(nemId_t id) const override {
73  return std::map<nemId_t, std::vector<double>>();
74  }
75  int getCellDataIdx(const std::string &name) override { return 0; }
76  std::vector<double> getCellLengths() const override {
77  return std::vector<double>();
78  }
79  std::vector<nemId_t> getConnectivities() const override {
80  return std::vector<nemId_t>();
81  }
82 
83  // check for special conditions
84  public:
85  bool chkDuplElm() const; // finds duplicate elements
86 
87  // point search methods
88  public:
89  // get coordinates and connectivities of the surface triangulation and returns
90  // ids for the nodes that reside on the triangulation within given tolerance.
91  void FindPntsOnTriSrf(const std::vector<double> &crds,
92  const std::vector<nemId_t> &conn,
93  std::set<nemId_t> &ids, double tol = 0.1e-15) const;
94  // get coordinates of the start and end points of an edge and returns ids for
95  // the nodes that reside on the edge within given tolerance.
96  void FindPntsOnEdge(std::vector<double> &crds, std::set<nemId_t> &ids,
97  double tol = 0.1e-15) const;
98 
99  // cell search methods
100  public:
101  // get vtkPolyData and return ids for the cells whose center reside inside the
102  // polyData.
103  void FindCellsInPolyData(vtkPolyData *polyData, std::vector<nemId_t> &ids,
104  bool query3Donly = true, double tol = 0.1e-15) const;
105  // get coordinates of a bounding box [xmin,xmax, ymin,ymax, zmin,zmax] and
106  // return ids for the cells whose center reside inside the bounding box.
107  void FindCellsWithinBounds(std::vector<double> &bb, std::vector<nemId_t> &ids,
108  bool fulImrsd = true);
109  // get coordinates and connectivities of the surface triangulation of a closed
110  // manifold and return ids for the cells whose center reside inside the
111  // manifold.
112  void FindCellsInTriSrf(const std::vector<std::vector<double>> &crds,
113  const std::vector<std::vector<vtkIdType>> &conns,
114  std::vector<nemId_t> &ids, bool query3Donly = true,
115  double tol = 0.1e-15) const;
116  // get center and radius of sphere and return ids for the cells whose center
117  // reside inside the sphere.
118  void FindCellsInSphere(const std::vector<double> &center, double radius,
119  std::vector<nemId_t> &ids, bool query3Donly = true,
120  double tol = 0.1e-15) const;
121 
122  // misc
123  public:
124  void inspectEdges(const std::string &ofname) const override {}
125  vtkSmartPointer<vtkDataSet> extractSurface() override {
126  return vtkSmartPointer<vtkDataSet>();
127  }
128 
129  void read(const std::string &fname = std::string()) override {}
130  using meshBase::write;
131  void write(const std::string &fname) const override {}
132 
133  // internal management
134  private:
135  void buildCellLocator();
136 
137  private:
138  bool upd_vcl;
139  vtkSmartPointer<vtkCellLocator> vcl;
140 };
141 
142 #endif // NEMOSYS_MESHSRCH_H_
int getCellType() const override
get cell type as an integer assumes all elements are the same type
Definition: meshSrch.H:68
std::vector< nemId_t > getConnectivities() const override
get connectivities.
Definition: meshSrch.H:79
void read(const std::string &fname=std::string()) override
abstract read method reserved for derived classes
Definition: meshSrch.H:129
A brief description of meshBase.
Definition: meshBase.H:64
std::size_t nemId_t
Definition: meshBase.H:51
bool upd_vcl
Definition: meshSrch.H:138
void inspectEdges(const std::string &ofname) const override
get edge lengths of dataSet
Definition: meshSrch.H:124
void write(const std::string &fname) const override
write the mesh to file named fname
Definition: meshSrch.H:131
vtkSmartPointer< vtkDataSet > extractSurface() override
extract the surface mesh
Definition: meshSrch.H:125
virtual std::vector< std::vector< double > > getCellVec(nemId_t id) const =0
get vector of coords of cell with id
std::vector< double > getCellLengths() const override
get diameter of circumsphere of each cell
Definition: meshSrch.H:76
virtual void write() const
write the mesh to file named after the private var &#39;filename&#39;.
Definition: meshBase.H:598
vtkSmartPointer< vtkCellLocator > vcl
Definition: meshSrch.H:139
static meshSrch * Create(meshBase *mb)
Definition: meshSrch.H:45
virtual std::vector< double > getCellCenter(nemId_t cellID) const =0
get center of a cell
int getCellDataIdx(const std::string &name) override
<>
Definition: meshSrch.H:75
static std::shared_ptr< meshSrch > CreateShared(meshBase *mb)
Definition: meshSrch.H:50
std::vector< std::vector< double > > getVertCrds() const override
get 3 vecs with x,y and z coords
Definition: meshSrch.H:69
std::map< nemId_t, std::vector< double > > getCell(nemId_t id) const override
get cell with id
Definition: meshSrch.H:72
static std::unique_ptr< meshSrch > CreateUnique(meshBase *mb)
Definition: meshSrch.H:56
virtual std::vector< double > getPoint(nemId_t id) const =0
get point with id
meshSrch(meshBase *mb)
Definition: meshSrch.H:43