Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Member Functions | List of all members
Nektar::Utilities::ProcessProjectCAD Class Reference

#include <ProcessProjectCAD.h>

Inheritance diagram for Nektar::Utilities::ProcessProjectCAD:
[legend]

Public Member Functions

 ProcessProjectCAD (NekMeshUtils::MeshSharedPtr m)
 
virtual ~ProcessProjectCAD ()
 
virtual void Process ()
 Write mesh to output file. More...
 
- Public Member Functions inherited from Nektar::NekMeshUtils::ProcessModule
NEKMESHUTILS_EXPORT ProcessModule (MeshSharedPtr p_m)
 
- Public Member Functions inherited from Nektar::NekMeshUtils::Module
NEKMESHUTILS_EXPORT Module (MeshSharedPtr p_m)
 
NEKMESHUTILS_EXPORT void RegisterConfig (std::string key, std::string value=std::string())
 Register a configuration option with a module. More...
 
NEKMESHUTILS_EXPORT void PrintConfig ()
 Print out all configuration options for a module. More...
 
NEKMESHUTILS_EXPORT void SetDefaults ()
 Sets default configuration options for those which have not been set. More...
 
NEKMESHUTILS_EXPORT MeshSharedPtr GetMesh ()
 
virtual NEKMESHUTILS_EXPORT void ProcessVertices ()
 Extract element vertices. More...
 
virtual NEKMESHUTILS_EXPORT void ProcessEdges (bool ReprocessEdges=true)
 Extract element edges. More...
 
virtual NEKMESHUTILS_EXPORT void ProcessFaces (bool ReprocessFaces=true)
 Extract element faces. More...
 
virtual NEKMESHUTILS_EXPORT void ProcessElements ()
 Generate element IDs. More...
 
virtual NEKMESHUTILS_EXPORT void ProcessComposites ()
 Generate composites. More...
 
virtual NEKMESHUTILS_EXPORT void ClearElementLinks ()
 

Static Public Member Functions

static std::shared_ptr< Modulecreate (NekMeshUtils::MeshSharedPtr m)
 Creates an instance of this class. More...
 

Static Public Attributes

static NekMeshUtils::ModuleKey className
 

Private Member Functions

bool findAndProject (bgi::rtree< boxI, bgi::quadratic< 16 > > &rtree, Array< OneD, NekDouble > in, int &surf)
 
bool IsNotValid (std::vector< NekMeshUtils::ElementSharedPtr > &els)
 

Additional Inherited Members

- Protected Member Functions inherited from Nektar::NekMeshUtils::Module
NEKMESHUTILS_EXPORT void ReorderPrisms (PerMap &perFaces)
 Reorder node IDs so that prisms and tetrahedra are aligned correctly. More...
 
NEKMESHUTILS_EXPORT void PrismLines (int prism, PerMap &perFaces, std::set< int > &prismsDone, std::vector< ElementSharedPtr > &line)
 
- Protected Attributes inherited from Nektar::NekMeshUtils::Module
MeshSharedPtr m_mesh
 Mesh object. More...
 
std::map< std::string, ConfigOptionm_config
 List of configuration values. More...
 

Detailed Description

Definition at line 59 of file ProcessProjectCAD.h.

Constructor & Destructor Documentation

◆ ProcessProjectCAD()

Nektar::Utilities::ProcessProjectCAD::ProcessProjectCAD ( NekMeshUtils::MeshSharedPtr  m)

Definition at line 63 of file ProcessProjectCAD.cpp.

References Nektar::NekMeshUtils::Module::m_config.

63  : ProcessModule(m)
64 {
65  m_config["file"] =
66  ConfigOption(false, "", "CAD file");
67  m_config["order"] =
68  ConfigOption(false, "4", "CAD file");
69 }
Represents a command-line configuration option.
NEKMESHUTILS_EXPORT ProcessModule(MeshSharedPtr p_m)
std::map< std::string, ConfigOption > m_config
List of configuration values.

◆ ~ProcessProjectCAD()

Nektar::Utilities::ProcessProjectCAD::~ProcessProjectCAD ( )
virtual

Definition at line 71 of file ProcessProjectCAD.cpp.

72 {
73 }

Member Function Documentation

◆ create()

static std::shared_ptr<Module> Nektar::Utilities::ProcessProjectCAD::create ( NekMeshUtils::MeshSharedPtr  m)
inlinestatic

Creates an instance of this class.

Definition at line 63 of file ProcessProjectCAD.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

64  {
66  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

◆ findAndProject()

bool Nektar::Utilities::ProcessProjectCAD::findAndProject ( bgi::rtree< boxI, bgi::quadratic< 16 > > &  rtree,
Array< OneD, NekDouble in,
int &  surf 
)
private

Definition at line 75 of file ProcessProjectCAD.cpp.

References Nektar::NekMeshUtils::Module::m_mesh.

Referenced by Process().

78 {
79  boost::ignore_unused(surf);
80 
81  point q(in[0], in[1], in[2]);
82  vector<boxI> result;
83  rtree.query(bgi::intersects(q), back_inserter(result));
84 
85  if(result.size() == 0)
86  {
87  //along a projecting edge the node is too far from any surface boxes
88  //this is hardly surprising but rare, return false and linearise
89  return false;
90  }
91 
92  int minsurf = 0;
93  NekDouble minDist = numeric_limits<double>::max();
94  NekDouble dist;
95 
96  for (int j = 0; j < result.size(); j++)
97  {
98  m_mesh->m_cad->GetSurf(result[j].second)->locuv(in, dist);
99 
100  if (dist < minDist)
101  {
102  minDist = dist;
103  minsurf = result[j].second;
104  }
105  }
106 
107  Array<OneD, NekDouble> uv = m_mesh->m_cad->GetSurf(minsurf)->locuv(in, dist);
108 
109  return true;
110 }
double NekDouble
bg::model::point< double, 3, bg::cs::cartesian > point
Definition: BLMesh.cpp:53

◆ IsNotValid()

bool Nektar::Utilities::ProcessProjectCAD::IsNotValid ( std::vector< NekMeshUtils::ElementSharedPtr > &  els)
private

Definition at line 112 of file ProcessProjectCAD.cpp.

References ASSERTL0, Nektar::LibUtilities::ePrism, Nektar::LibUtilities::eTetrahedron, and Nektar::NekConstants::kNekZeroTol.

Referenced by Process().

113 {
114  //short algebraic method to figure out the vailidy of elements
115  for(int i = 0; i < els.size(); i++)
116  {
117  if(els[i]->GetShapeType() == LibUtilities::ePrism)
118  {
119  vector<NodeSharedPtr> ns = els[i]->GetVertexList();
120  for(int j = 0; j < 6; j++)
121  {
122  NekDouble a2 = 0.5 * (1 + prismU1[j]);
123  NekDouble b1 = 0.5 * (1 - prismV1[j]);
124  NekDouble b2 = 0.5 * (1 + prismV1[j]);
125  NekDouble c2 = 0.5 * (1 + prismW1[j]);
126  NekDouble d = 0.5 * (prismU1[j] + prismW1[j]);
127 
128  Array<OneD, NekDouble> jac(9,0.0);
129 
130  jac[0] = -0.5 * b1 * ns[0]->m_x + 0.5 * b1 * ns[1]->m_x +
131  0.5 * b2 * ns[2]->m_x - 0.5 * b2 * ns[3]->m_x;
132  jac[1] = -0.5 * b1 * ns[0]->m_y + 0.5 * b1 * ns[1]->m_y +
133  0.5 * b2 * ns[2]->m_y - 0.5 * b2 * ns[3]->m_y;
134  jac[2] = -0.5 * b1 * ns[0]->m_z + 0.5 * b1 * ns[1]->m_z +
135  0.5 * b2 * ns[2]->m_z - 0.5 * b2 * ns[3]->m_z;
136 
137  jac[3] = 0.5 * d * ns[0]->m_x - 0.5 * a2 * ns[1]->m_x +
138  0.5 * a2 * ns[2]->m_x - 0.5 * d * ns[3]->m_x -
139  0.5 * c2 * ns[4]->m_x + 0.5 * c2 * ns[5]->m_x;
140  jac[4] = 0.5 * d * ns[0]->m_y - 0.5 * a2 * ns[1]->m_y +
141  0.5 * a2 * ns[2]->m_y - 0.5 * d * ns[3]->m_y -
142  0.5 * c2 * ns[4]->m_y + 0.5 * c2 * ns[5]->m_y;
143  jac[5] = 0.5 * d * ns[0]->m_z - 0.5 * a2 * ns[1]->m_z +
144  0.5 * a2 * ns[2]->m_z - 0.5 * d * ns[3]->m_z -
145  0.5 * c2 * ns[4]->m_z + 0.5 * c2 * ns[5]->m_z;
146 
147  jac[6] = -0.5 * b1 * ns[0]->m_x - 0.5 * b2 * ns[3]->m_x +
148  0.5 * b1 * ns[4]->m_x + 0.5 * b2 * ns[5]->m_x;
149  jac[7] = -0.5 * b1 * ns[0]->m_y - 0.5 * b2 * ns[3]->m_y +
150  0.5 * b1 * ns[4]->m_y + 0.5 * b2 * ns[5]->m_y;
151  jac[8] = -0.5 * b1 * ns[0]->m_z - 0.5 * b2 * ns[3]->m_z +
152  0.5 * b1 * ns[4]->m_z + 0.5 * b2 * ns[5]->m_z;
153 
154  NekDouble jc = jac[0] * (jac[4] * jac[8] - jac[5] * jac[7]) -
155  jac[3] * (jac[1] * jac[8] - jac[2] * jac[7]) +
156  jac[6] * (jac[1] * jac[5] - jac[2] * jac[4]);
157 
158  if (jc < NekConstants::kNekZeroTol)
159  {
160  return true;
161  }
162  }
163  }
164  else if(els[i]->GetShapeType() == LibUtilities::eTetrahedron)
165  {
166  vector<NodeSharedPtr> ns = els[i]->GetVertexList();
167  for(int j = 0; j < 4; j++)
168  {
169  Array<OneD, NekDouble> jac(9,0.0);
170 
171  jac[0] = 0.5*(ns[1]->m_x - ns[0]->m_x);
172  jac[1] = 0.5*(ns[1]->m_y - ns[0]->m_y);
173  jac[2] = 0.5*(ns[1]->m_z - ns[0]->m_z);
174  jac[3] = 0.5*(ns[2]->m_x - ns[0]->m_x);
175  jac[4] = 0.5*(ns[2]->m_y - ns[0]->m_y);
176  jac[5] = 0.5*(ns[2]->m_z - ns[0]->m_z);
177  jac[6] = 0.5*(ns[3]->m_x - ns[0]->m_x);
178  jac[7] = 0.5*(ns[3]->m_y - ns[0]->m_y);
179  jac[8] = 0.5*(ns[3]->m_z - ns[0]->m_z);
180 
181  NekDouble jc = jac[0] * (jac[4] * jac[8] - jac[5] * jac[7]) -
182  jac[3] * (jac[1] * jac[8] - jac[2] * jac[7]) +
183  jac[6] * (jac[1] * jac[5] - jac[2] * jac[4]);
184 
185  if (jc < NekConstants::kNekZeroTol)
186  {
187  return true;
188  }
189  }
190  }
191  else
192  {
193  ASSERTL0(false, "not coded");
194  }
195  }
196 
197  return false;
198 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
const NekDouble prismV1[6]
const NekDouble prismW1[6]
static const NekDouble kNekZeroTol
double NekDouble
const NekDouble prismU1[6]

◆ Process()

void Nektar::Utilities::ProcessProjectCAD::Process ( )
virtual

Write mesh to output file.

Implements Nektar::NekMeshUtils::Module.

Definition at line 200 of file ProcessProjectCAD.cpp.

References ASSERTL0, Nektar::NekMeshUtils::Module::ClearElementLinks(), Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::FieldUtils::eProcessModule, findAndProject(), Nektar::FieldUtils::GetModuleFactory(), IsNotValid(), CG_Iterations::loc, Nektar::NekMeshUtils::Module::m_config, Nektar::NekMeshUtils::Module::m_mesh, class_topology::Node, Nektar::LibUtilities::PointsManager(), and Nektar::LibUtilities::PrintProgressbar().

201 {
202  if (m_mesh->m_verbose)
203  {
204  cout << "ProcessAssignCAD: Processing... " << endl;
205  }
206 
207  cout << "ProcessAssignCAD: Warning: This module is designed for use with "
208  "starccm+ meshes only, it also requires that the star mesh was "
209  "created in a certain way" << endl;
210 
211  if(!m_config["order"].beenSet)
212  {
213  cout << "Mesh order not set will assume 4" << endl;
214  }
215 
216  m_mesh->m_nummode = m_config["order"].as<int>() + 1;
217 
218  //load instance of the CAD model
220  ModuleKey(eProcessModule, "loadcad"), m_mesh);
221  module->RegisterConfig("filename", m_config["file"].as<string>());
222  module->SetDefaults();
223  module->Process();
224 
225  //turn the bounding box of the CAD surfaces into a k-d tree
226  vector<boxI> boxes;
227  for (int i = 1; i <= m_mesh->m_cad->GetNumSurf(); i++)
228  {
229  LibUtilities::PrintProgressbar(i, m_mesh->m_cad->GetNumSurf(), "building surface bboxes", i-1);
230  Array<OneD, NekDouble> bx = m_mesh->m_cad->GetSurf(i)->BoundingBox();
231  boxes.push_back(make_pair(box(point(bx[0], bx[1], bx[2]), point(bx[3], bx[4], bx[5])), i));
232  }
233  cout << endl;
234 
235  cout << "building admin data structures" << endl;
236 
237  bgi::rtree<boxI, bgi::quadratic<16> > rtree(boxes);
238 
239  NodeSet surfNodes;
240 
241  //find unique nodes on the surface
242  for(int i = 0; i < m_mesh->m_element[2].size(); i++)
243  {
244  ElementSharedPtr el = m_mesh->m_element[2][i];
245  vector<NodeSharedPtr> ns = el->GetVertexList();
246  for(int j = 0; j < ns.size(); j++)
247  {
248  surfNodes.insert(ns[j]);
249  }
250  }
251 
252  map<NodeSharedPtr, vector<ElementSharedPtr> > surfNodeToEl;
253 
254  //link surface nodes to their 3D element
255  for(int i = 0; i < m_mesh->m_element[3].size(); i++)
256  {
257  if(m_mesh->m_element[3][i]->HasBoundaryLinks())
258  {
259  vector<NodeSharedPtr> ns = m_mesh->m_element[3][i]->GetVertexList();
260  for(int j = 0; j < ns.size(); j++)
261  {
262  if(surfNodes.count(ns[j]) > 0)
263  {
264  surfNodeToEl[ns[j]].push_back(m_mesh->m_element[3][i]);
265  }
266  }
267  }
268  }
269 
270  //link the surface node to a value for the shortest connecting edge to it
271  map<NodeSharedPtr, NekDouble> minConEdge;
272  for (int i = 0; i < m_mesh->m_element[2].size(); i++)
273  {
274  ElementSharedPtr el = m_mesh->m_element[2][i];
275  vector<NodeSharedPtr> ns = el->GetVertexList();
276  NekDouble l1 = ns[0]->Distance(ns[1]);
277  NekDouble l2 = ns[1]->Distance(ns[2]);
278  NekDouble l3 = ns[2]->Distance(ns[0]);
279 
280  if (minConEdge.count(ns[0]))
281  {
282  NekDouble l = minConEdge[ns[0]];
283  minConEdge[ns[0]] = min(l, min(l1, l3));
284  }
285  else
286  {
287  minConEdge.insert(make_pair(ns[0], min(l1, l3)));
288  }
289  if (minConEdge.count(ns[1]))
290  {
291  NekDouble l = minConEdge[ns[1]];
292  minConEdge[ns[1]] = min(l, min(l1, l1));
293  }
294  else
295  {
296  minConEdge.insert(make_pair(ns[1], min(l1, l1)));
297  }
298  if (minConEdge.count(ns[2]))
299  {
300  NekDouble l = minConEdge[ns[2]];
301  minConEdge[ns[2]] = min(l, min(l2, l3));
302  }
303  else
304  {
305  minConEdge.insert(make_pair(ns[2], min(l2, l3)));
306  }
307  }
308 
309  map<int, vector<int> > finds;
310 
311  cout << "searching tree" << endl;
312 
313  NekDouble maxNodeCor = 0;
314 
315  // this is a set of nodes which have a CAD failure
316  // if touched in the HO stage they should be ignored and linearised
317  NodeSet lockedNodes;
318 
319  //find nodes surface and parametric location
320  int ct = 0;
321  for (auto i = surfNodes.begin(); i != surfNodes.end(); i++, ct++)
322  {
323  LibUtilities::PrintProgressbar(ct, surfNodes.size(), "projecting verts", ct-1);
324 
325  point q((*i)->m_x, (*i)->m_y, (*i)->m_z);
326  vector<boxI> result;
327  rtree.query(bgi::intersects(q), back_inserter(result));
328 
329  ASSERTL0(result.size() > 0, "node thinks its not in any boxes");
330 
331  NekDouble tol = minConEdge[*i] * 0.5;
332  tol = min(tol, 5e-4);
333  tol = max(tol, 1e-6);
334 
335  vector<int> distId;
336  vector<NekDouble> distList;
337 
338  for (int j = 0; j < result.size(); j++)
339  {
340  NekDouble dist;
341  m_mesh->m_cad->GetSurf(result[j].second)->locuv((*i)->GetLoc(), dist);
342  distList.push_back(dist);
343  distId.push_back(result[j].second);
344  }
345 
346  bool repeat = true;
347  while (repeat)
348  {
349  repeat = false;
350  for (int j = 0; j < distId.size() - 1; j++)
351  {
352  if (distList[j+1] < distList[j])
353  {
354  repeat = true;
355  swap(distList[j+1], distList[j]);
356  swap(distId[j+1], distId[j]);
357  }
358  }
359  }
360 
361  int pos = 0;
362  for (int j = 0; j < distId.size(); j++)
363  {
364  if (distList[j] < tol)
365  {
366  pos++;
367  }
368  }
369 
370  distId.resize(pos);
371 
372  finds[pos].push_back(0);
373 
374  if (pos == 0)
375  {
376  lockedNodes.insert(*i);
377  cout << endl << "WARNING: surface unknown " << distList[0] << " " << tol << endl;
378  }
379  else
380  {
381  NekDouble shift;
382  bool st = false;
383  for (int j = 0; j < distId.size(); j++)
384  {
385  if (distList[j] > tol)
386  {
387  continue;
388  }
389  if (m_mesh->m_cad->GetSurf(distId[j])->IsPlanar())
390  {
391  continue;
392  }
393 
394  shift = distList[j];
395  Array<OneD, NekDouble> uvt(2);
396  NekDouble dist = 0;
397  CADSurfSharedPtr s = m_mesh->m_cad->GetSurf(distId[j]);
398  Array<OneD, NekDouble> l = (*i)->GetLoc();
399  uvt = s->locuv(l, dist);
400 
401  NekDouble tmpX = (*i)->m_x;
402  NekDouble tmpY = (*i)->m_y;
403  NekDouble tmpZ = (*i)->m_z;
404 
405  (*i)->m_x = l[0];
406  (*i)->m_y = l[1];
407  (*i)->m_z = l[2];
408 
409  if(IsNotValid(surfNodeToEl[*i]))
410  {
411  (*i)->m_x = tmpX;
412  (*i)->m_y = tmpY;
413  (*i)->m_z = tmpZ;
414  cout << "reset element projection" << endl;
415  break;
416  }
417 
418  st = true;
419  break;
420  }
421 
422  if (!st)
423  {
424  lockedNodes.insert(*i);
425  continue;
426  }
427 
428  for (int j = 0; j < distId.size(); j++)
429  {
430  if (distList[j] > tol)
431  {
432  continue;
433  }
434  if (m_mesh->m_cad->GetSurf(distId[j])->IsPlanar())
435  {
436  continue;
437  }
438 
439  CADSurfSharedPtr s = m_mesh->m_cad->GetSurf(distId[j]);
440  Array<OneD, NekDouble> uv;
441  NekDouble dist=0;
442  Array<OneD, NekDouble> loc = (*i)->GetLoc();
443  uv = s->locuv(loc, dist);
444  (*i)->SetCADSurf(s, uv);
445  }
446  maxNodeCor = max(maxNodeCor, shift);
447  }
448  }
449 
450  cout << endl;
451 
452  cout << "max surface node correction " << maxNodeCor << endl;
453 
454  //make edges of surface mesh unique
456  EdgeSet surfEdges;
457  vector<ElementSharedPtr> &elmt = m_mesh->m_element[2];
458  map<int, int> surfIdToLoc;
459  for (int i = 0; i < elmt.size(); i++)
460  {
461  surfIdToLoc.insert(make_pair(elmt[i]->GetId(), i));
462  for (int j = 0; j < elmt[i]->GetEdgeCount(); ++j)
463  {
464  pair<EdgeSet::iterator,bool> testIns;
465  EdgeSharedPtr ed = elmt[i]->GetEdge(j);
466  testIns = surfEdges.insert(ed);
467 
468  if (testIns.second)
469  {
470  EdgeSharedPtr ed2 = *testIns.first;
471  ed2->m_elLink.push_back(
472  pair<ElementSharedPtr,int>(elmt[i],j));
473  }
474  else
475  {
476  EdgeSharedPtr e2 = *(testIns.first);
477  elmt[i]->SetEdge(j, e2);
478 
479  // Update edge to element map.
480  e2->m_elLink.push_back( pair<ElementSharedPtr,int>(elmt[i],j));
481  }
482  }
483  }
484 
485  int order = m_config["order"].as<int>();
486 
487  map<int, vector<int> > eds;
488 
489  LibUtilities::PointsKey ekey(order+1, LibUtilities::eGaussLobattoLegendre);
490  Array<OneD, NekDouble> gll;
491  LibUtilities::PointsManager()[ekey]->GetPoints(gll);
492 
493  //make surface edges high-order
494  for (auto i = surfEdges.begin(); i != surfEdges.end(); i++)
495  {
496  if (lockedNodes.count((*i)->m_n1) || lockedNodes.count((*i)->m_n2))
497  {
498  continue;
499  }
500  vector<CADSurfSharedPtr> v1 = (*i)->m_n1->GetCADSurfs();
501  vector<CADSurfSharedPtr> v2 = (*i)->m_n2->GetCADSurfs();
502 
503  vector<int> vi1, vi2;
504  for (size_t j = 0; j < v1.size(); ++j)
505  {
506  vi1.push_back(v1[j]->GetId());
507  }
508  for (size_t j = 0; j < v2.size(); ++j)
509  {
510  vi2.push_back(v2[j]->GetId());
511  }
512 
513  sort(vi1.begin(), vi1.end());
514  sort(vi2.begin(), vi2.end());
515 
516  vector<int> cmn;
517  set_intersection(vi1.begin(), vi1.end(), vi2.begin(), vi2.end(), back_inserter(cmn));
518  eds[cmn.size()].push_back(0);
519 
520  (*i)->m_curveType = LibUtilities::eGaussLobattoLegendre;
521 
522  if (cmn.size() == 1 || cmn.size() == 2)
523  {
524  for (int j = 0; j < cmn.size(); j++)
525  {
526  if (m_mesh->m_cad->GetSurf(cmn[j])->IsPlanar())
527  {
528  // if its planar dont care
529  continue;
530  }
531 
532  Array<OneD, NekDouble> uvb = (*i)->m_n1->GetCADSurfInfo(cmn[j]);
533  Array<OneD, NekDouble> uve = (*i)->m_n2->GetCADSurfInfo(cmn[j]);
534 
535  // can compare the loction of the projection to the
536  // corresponding position of the straight sided edge
537  // if the two differ by more than the length of the edge
538  // something has gone wrong
539  NekDouble len = (*i)->m_n1->Distance((*i)->m_n2);
540 
541  for (int k = 1; k < order+1 - 1; k++)
542  {
543  Array<OneD, NekDouble> uv(2);
544  uv[0] = uvb[0] * (1.0 - gll[k]) / 2.0 + uve[0] * (1.0 + gll[k]) / 2.0;
545  uv[1] = uvb[1] * (1.0 - gll[k]) / 2.0 + uve[1] * (1.0 + gll[k]) / 2.0;
546  Array<OneD, NekDouble> loc;
547  loc = m_mesh->m_cad->GetSurf(cmn[j])->P(uv);
548  Array<OneD, NekDouble> locT(3);
549  locT[0] = (*i)->m_n1->m_x * (1.0 - gll[k]) / 2.0 +
550  (*i)->m_n2->m_x * (1.0 + gll[k]) / 2.0;
551  locT[1] = (*i)->m_n1->m_y * (1.0 - gll[k]) / 2.0 +
552  (*i)->m_n2->m_y * (1.0 + gll[k]) / 2.0;
553  locT[2] = (*i)->m_n1->m_z * (1.0 - gll[k]) / 2.0 +
554  (*i)->m_n2->m_z * (1.0 + gll[k]) / 2.0;
555 
556  NekDouble d = sqrt((locT[0] - loc[0]) * (locT[0] - loc[0]) +
557  (locT[1] - loc[1]) * (locT[1] - loc[1]) +
558  (locT[2] - loc[2]) * (locT[2] - loc[2]));
559 
560  if (d > len)
561  {
562  (*i)->m_edgeNodes.clear();
563  break;
564  }
565 
566  NodeSharedPtr nn = std::shared_ptr<Node>(
567  new Node(0, loc[0], loc[1], loc[2]));
568 
569  (*i)->m_edgeNodes.push_back(nn);
570  }
571 
572  if ((*i)->m_edgeNodes.size() != 0)
573  {
574  // it suceeded on this surface so skip the other possibility
575  break;
576  }
577  }
578  }
579  else if (cmn.size() == 0)
580  {
581  // projection, if the projection requires more than two surfaces
582  // including the edge nodes, then, in theory projection shouldnt be
583  // used
584 
585  set<int> sused;
586  for (int k = 1; k < order+1 - 1; k++)
587  {
588  Array<OneD, NekDouble> locT(3);
589  locT[0] = (*i)->m_n1->m_x * (1.0 - gll[k]) / 2.0 +
590  (*i)->m_n2->m_x * (1.0 + gll[k]) / 2.0;
591  locT[1] = (*i)->m_n1->m_y * (1.0 - gll[k]) / 2.0 +
592  (*i)->m_n2->m_y * (1.0 + gll[k]) / 2.0;
593  locT[2] = (*i)->m_n1->m_z * (1.0 - gll[k]) / 2.0 +
594  (*i)->m_n2->m_z * (1.0 + gll[k]) / 2.0;
595 
596  int s;
597  if(!findAndProject(rtree, locT, s))
598  {
599  (*i)->m_edgeNodes.clear();
600  break;
601  }
602  sused.insert(s);
603 
604  if (sused.size() > 2)
605  {
606  (*i)->m_edgeNodes.clear();
607  break;
608  }
609 
610  NodeSharedPtr nn = std::shared_ptr<Node>(
611  new Node(0, locT[0], locT[1], locT[2]));
612 
613  (*i)->m_edgeNodes.push_back(nn);
614  }
615  }
616  }
617 
618 }
std::shared_ptr< CADSurf > CADSurfSharedPtr
Definition: CADCurve.h:51
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
int PrintProgressbar(const int position, const int goal, const std::string message, int lastprogress=-1)
Prints a progressbar.
Definition: Progressbar.hpp:67
std::unordered_set< EdgeSharedPtr, EdgeHash > EdgeSet
Definition: Edge.h:159
std::shared_ptr< Module > ModuleSharedPtr
std::shared_ptr< Edge > EdgeSharedPtr
Shared pointer to an edge.
Definition: Edge.h:136
std::shared_ptr< Node > NodeSharedPtr
Definition: CADVert.h:49
bool IsNotValid(std::vector< NekMeshUtils::ElementSharedPtr > &els)
std::unordered_set< NodeSharedPtr, NodeHash > NodeSet
Definition: Node.h:447
bg::model::box< point > box
Definition: BLMesh.cpp:54
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
std::pair< ModuleType, std::string > ModuleKey
std::shared_ptr< Element > ElementSharedPtr
Definition: Edge.h:49
PointsManagerT & PointsManager(void)
double NekDouble
std::map< std::string, ConfigOption > m_config
List of configuration values.
virtual NEKMESHUTILS_EXPORT void ClearElementLinks()
bg::model::point< double, 3, bg::cs::cartesian > point
Definition: BLMesh.cpp:53
bool findAndProject(bgi::rtree< boxI, bgi::quadratic< 16 > > &rtree, Array< OneD, NekDouble > in, int &surf)
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
ModuleFactory & GetModuleFactory()

Member Data Documentation

◆ className

ModuleKey Nektar::Utilities::ProcessProjectCAD::className
static
Initial value:

Definition at line 67 of file ProcessProjectCAD.h.