Nektar++
Public Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
Nektar::NekMeshUtils::FaceMesh Class Reference

class for surface meshes on individual surfaces (paramter plane meshes) More...

#include <FaceMesh.h>

Public Member Functions

 FaceMesh (const int id, MeshSharedPtr m, const std::map< int, CurveMeshSharedPtr > &cmeshes, const int comp)
 Default constructor. More...
 
void Mesh ()
 mesh exectuation command More...
 
bool ValidateCurves ()
 validate the curve meshes More...
 
void ValidateLoops ()
 validate the curve meshes considering the loops More...
 

Private Member Functions

void OrientateCurves ()
 Get the boundries of the surface and extracts the nodes from the curve meshes in the correct order. More...
 
void Stretching ()
 Calculate the paramter plane streching factor. More...
 
void Smoothing ()
 performs node smoothing on face More...
 
void DiagonalSwap ()
 performs diagonal swapping of edges More...
 
void BuildLocalMesh ()
 build a local version of mesh elements More...
 
void OptimiseLocalMesh ()
 function which calls the optimisation routines More...
 
bool Validate ()
 Validate the surface mesh base on the octree and real dimensions of the edges. More...
 
void AddNewPoint (Array< OneD, NekDouble > uv)
 adds a new stiener point to the triangulation for meshing More...
 
void MakeBL ()
 adds a quad layer around any interior loops More...
 

Private Attributes

MeshSharedPtr m_mesh
 mesh pointer More...
 
CADSurfSharedPtr m_cadsurf
 CAD surface. More...
 
std::map< int, CurveMeshSharedPtrm_curvemeshes
 Map of the curve meshes which bound the surfaces. More...
 
std::vector< EdgeLoopSharedPtrm_edgeloops
 data structure containing the edges, their order and oreientation for the surface More...
 
int m_id
 id of the surface mesh More...
 
std::vector< std::vector< NodeSharedPtr > > orderedLoops
 list of boundary nodes in their order loops More...
 
std::vector< NodeSharedPtrm_stienerpoints
 list of stiener points in the triangulation More...
 
NekDouble m_str
 pplane stretching More...
 
std::vector< std::vector< NodeSharedPtr > > m_connec
 triangle connectiviities More...
 
NodeSet m_localNodes
 local set of nodes More...
 
EdgeSet m_localEdges
 local set of edges More...
 
std::vector< ElementSharedPtrm_localElements
 local list of elements More...
 
NodeSet m_inBoundary
 set of nodes which are in the boundary (easier to identify conflicts with) More...
 
int m_compId
 identity to put into element tags More...
 

Friends

class MemoryManager< FaceMesh >
 

Detailed Description

class for surface meshes on individual surfaces (paramter plane meshes)

Definition at line 51 of file FaceMesh.h.

Constructor & Destructor Documentation

◆ FaceMesh()

Nektar::NekMeshUtils::FaceMesh::FaceMesh ( const int  id,
MeshSharedPtr  m,
const std::map< int, CurveMeshSharedPtr > &  cmeshes,
const int  comp 
)
inline

Default constructor.

Definition at line 59 of file FaceMesh.h.

References AddNewPoint(), BuildLocalMesh(), DiagonalSwap(), m_cadsurf, m_edgeloops, m_id, m_mesh, MakeBL(), Mesh(), OptimiseLocalMesh(), OrientateCurves(), Smoothing(), Stretching(), Validate(), ValidateCurves(), and ValidateLoops().

63  : m_mesh(m), m_curvemeshes(cmeshes), m_id(id), m_compId(comp)
64 
65  {
66  m_cadsurf = m_mesh->m_cad->GetSurf(m_id);
67  m_edgeloops = m_cadsurf->GetEdges();
68  };
std::map< int, CurveMeshSharedPtr > m_curvemeshes
Map of the curve meshes which bound the surfaces.
Definition: FaceMesh.h:139
std::vector< EdgeLoopSharedPtr > m_edgeloops
data structure containing the edges, their order and oreientation for the surface ...
Definition: FaceMesh.h:142
int m_compId
identity to put into element tags
Definition: FaceMesh.h:162
int m_id
id of the surface mesh
Definition: FaceMesh.h:144
CADSurfSharedPtr m_cadsurf
CAD surface.
Definition: FaceMesh.h:137
MeshSharedPtr m_mesh
mesh pointer
Definition: FaceMesh.h:135

Member Function Documentation

◆ AddNewPoint()

void Nektar::NekMeshUtils::FaceMesh::AddNewPoint ( Array< OneD, NekDouble uv)
private

adds a new stiener point to the triangulation for meshing

Definition at line 1129 of file FaceMesh.cpp.

References class_topology::Node.

Referenced by FaceMesh().

1130 {
1131  // adds a new point but checks that there are no other points nearby first
1132  Array<OneD, NekDouble> np = m_cadsurf->P(uv);
1133  NekDouble npDelta = m_mesh->m_octree->Query(np);
1134 
1135  NodeSharedPtr n = std::shared_ptr<Node>(
1136  new Node(m_mesh->m_numNodes++, np[0], np[1], np[2]));
1137 
1138  bool add = true;
1139 
1140  for (int i = 0; i < orderedLoops.size(); i++)
1141  {
1142  for (int j = 0; j < orderedLoops[i].size(); j++)
1143  {
1144  NekDouble r = orderedLoops[i][j]->Distance(n);
1145 
1146  if (r < npDelta / 2.0)
1147  {
1148  add = false;
1149  break;
1150  }
1151  }
1152  }
1153 
1154  if (add)
1155  {
1156  for (int i = 0; i < m_stienerpoints.size(); i++)
1157  {
1158  NekDouble r = m_stienerpoints[i]->Distance(n);
1159 
1160  if (r < npDelta / 2.0)
1161  {
1162  add = false;
1163  break;
1164  }
1165  }
1166  }
1167 
1168  if (add)
1169  {
1170  n->SetCADSurf(m_cadsurf, uv);
1171  m_stienerpoints.push_back(n);
1172  }
1173 }
std::vector< std::vector< NodeSharedPtr > > orderedLoops
list of boundary nodes in their order loops
Definition: FaceMesh.h:146
std::shared_ptr< Node > NodeSharedPtr
Definition: CADVert.h:49
double NekDouble
std::vector< NodeSharedPtr > m_stienerpoints
list of stiener points in the triangulation
Definition: FaceMesh.h:148
CADSurfSharedPtr m_cadsurf
CAD surface.
Definition: FaceMesh.h:137
MeshSharedPtr m_mesh
mesh pointer
Definition: FaceMesh.h:135

◆ BuildLocalMesh()

void Nektar::NekMeshUtils::FaceMesh::BuildLocalMesh ( )
private

build a local version of mesh elements

Definition at line 902 of file FaceMesh.cpp.

References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::eTriangle, and Nektar::NekMeshUtils::GetElementFactory().

Referenced by FaceMesh().

903 {
904  /*************************
905  // build a local set of nodes edges and elemenets for optimstaion prior to
906  putting them into m_mesh
907  */
908 
909  for (int i = 0; i < m_connec.size(); i++)
910  {
911  ElmtConfig conf(LibUtilities::eTriangle, 1, false, false,
912  m_mesh->m_spaceDim != 3);
913 
914  vector<int> tags;
915  tags.push_back(m_compId);
917  LibUtilities::eTriangle, conf, m_connec[i], tags);
918  E->m_parentCAD = m_cadsurf;
919 
920  vector<NodeSharedPtr> nods = E->GetVertexList();
921  for (int j = 0; j < nods.size(); j++)
922  {
923  // nodes are already unique some will insert some wont
924  m_localNodes.insert(nods[j]);
925  }
926 
927  E->SetId(m_localElements.size());
928  m_localElements.push_back(E);
929  }
930 
931  for (int i = 0; i < m_localElements.size(); ++i)
932  {
933  for (int j = 0; j < m_localElements[i]->GetEdgeCount(); ++j)
934  {
935  pair<EdgeSet::iterator, bool> testIns;
936  EdgeSharedPtr ed = m_localElements[i]->GetEdge(j);
937  // look for edge in m_mesh edgeset from curves
938  EdgeSet::iterator s = m_mesh->m_edgeSet.find(ed);
939  if (!(s == m_mesh->m_edgeSet.end()))
940  {
941  ed = *s;
942  m_localElements[i]->SetEdge(j, *s);
943  }
944 
945  testIns = m_localEdges.insert(ed);
946 
947  if (testIns.second)
948  {
949  EdgeSharedPtr ed2 = *testIns.first;
950  ed2->m_elLink.push_back(
951  pair<ElementSharedPtr, int>(m_localElements[i], j));
952  }
953  else
954  {
955  EdgeSharedPtr e2 = *(testIns.first);
956  m_localElements[i]->SetEdge(j, e2);
957  e2->m_elLink.push_back(
958  pair<ElementSharedPtr, int>(m_localElements[i], j));
959  }
960  }
961  }
962 }
std::vector< ElementSharedPtr > m_localElements
local list of elements
Definition: FaceMesh.h:158
std::shared_ptr< Edge > EdgeSharedPtr
Shared pointer to an edge.
Definition: Edge.h:136
ElementFactory & GetElementFactory()
Definition: Element.cpp:44
int m_compId
identity to put into element tags
Definition: FaceMesh.h:162
std::vector< std::vector< NodeSharedPtr > > m_connec
triangle connectiviities
Definition: FaceMesh.h:152
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
EdgeSet m_localEdges
local set of edges
Definition: FaceMesh.h:156
std::shared_ptr< Element > ElementSharedPtr
Definition: Edge.h:49
CADSurfSharedPtr m_cadsurf
CAD surface.
Definition: FaceMesh.h:137
MeshSharedPtr m_mesh
mesh pointer
Definition: FaceMesh.h:135
NodeSet m_localNodes
local set of nodes
Definition: FaceMesh.h:154

◆ DiagonalSwap()

void Nektar::NekMeshUtils::FaceMesh::DiagonalSwap ( )
private

performs diagonal swapping of edges

Definition at line 506 of file FaceMesh.cpp.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::eTriangle, and Nektar::NekMeshUtils::GetElementFactory().

Referenced by FaceMesh().

507 {
508  map<int, int> idealConnec;
509  map<int, int> actualConnec;
510  map<int, vector<EdgeSharedPtr> > nodetoedge;
511  // figure out ideal node count and actual node count
512  EdgeSet::iterator eit;
513  for (eit = m_localEdges.begin(); eit != m_localEdges.end(); eit++)
514  {
515  nodetoedge[(*eit)->m_n1->m_id].push_back(*eit);
516  nodetoedge[(*eit)->m_n2->m_id].push_back(*eit);
517  }
518  NodeSet::iterator nit;
519  for (nit = m_localNodes.begin(); nit != m_localNodes.end(); nit++)
520  {
521  // this routine is broken and needs looking at
522  NodeSet::iterator f = m_inBoundary.find((*nit));
523  if (f == m_inBoundary.end()) // node isnt on curve so skip
524  {
525  // node is interior
526  idealConnec[(*nit)->m_id] = 6;
527  }
528  else
529  {
530  // need to identify the two other nodes on the boundary to find
531  // interior angle
532  /*vector<NodeSharedPtr> ns;
533  vector<EdgeSharedPtr> e = nodetoedge[(*nit)->m_id];
534  for (int i = 0; i < e.size(); i++)
535  {
536  if (!e[i]->onCurve)
537  continue; // the linking nodes are not going to exist on
538  // interior edges
539 
540  if (e[i]->m_n1 == (*nit))
541  ns.push_back(e[i]->m_n2);
542  else
543  ns.push_back(e[i]->m_n1);
544  }
545  ASSERTL0(ns.size() == 2,
546  "failed to find 2 nodes in the angle system");
547 
548  idealConnec[(*nit)->m_id] =
549  ceil((*nit)->Angle(ns[0], ns[1]) / 3.142 * 3) + 1;*/
550  idealConnec[(*nit)->m_id] = 4;
551  }
552  }
553  for (nit = m_localNodes.begin(); nit != m_localNodes.end(); nit++)
554  {
555  actualConnec[(*nit)->m_id] = nodetoedge[(*nit)->m_id].size();
556  }
557 
558  // edgeswapping fun times
559  // perfrom edge swap based on node defect and then angle
560  for (int q = 0; q < 4; q++)
561  {
562  int edgesStart = m_localEdges.size();
563  EdgeSet edges = m_localEdges;
564  m_localEdges.clear();
565 
566  int swappedEdges = 0;
567 
568  EdgeSet::iterator it;
569 
570  for (it = edges.begin(); it != edges.end(); it++)
571  {
572  EdgeSharedPtr e = *it;
573 
574  NodeSet::iterator f1 = m_inBoundary.find((*it)->m_n1);
575  NodeSet::iterator f2 = m_inBoundary.find((*it)->m_n2);
576  if (f1 != m_inBoundary.end() && f2 != m_inBoundary.end())
577  {
578  m_localEdges.insert(e);
579  continue;
580  }
581 
582  ElementSharedPtr tri1 = e->m_elLink[0].first;
583  ElementSharedPtr tri2 = e->m_elLink[1].first;
584 
585  NodeSharedPtr n1 = e->m_n1;
586  NodeSharedPtr n2 = e->m_n2;
587 
588  vector<NodeSharedPtr> nt = tri1->GetVertexList();
589 
590  // identify node a,b,c,d of the swapping
591  NodeSharedPtr A, B, C, D;
592  if (nt[0] != n1 && nt[0] != n2)
593  {
594  C = nt[0];
595  B = nt[1];
596  A = nt[2];
597  }
598  else if (nt[1] != n1 && nt[1] != n2)
599  {
600  C = nt[1];
601  B = nt[2];
602  A = nt[0];
603  }
604  else if (nt[2] != n1 && nt[2] != n2)
605  {
606  C = nt[2];
607  B = nt[0];
608  A = nt[1];
609  }
610  else
611  {
612  ASSERTL0(false, "failed to identify verticies in tri1");
613  }
614 
615  nt = tri2->GetVertexList();
616 
617  if (nt[0] != n1 && nt[0] != n2)
618  {
619  D = nt[0];
620  }
621  else if (nt[1] != n1 && nt[1] != n2)
622  {
623  D = nt[1];
624  }
625  else if (nt[2] != n1 && nt[2] != n2)
626  {
627  D = nt[2];
628  }
629  else
630  {
631  ASSERTL0(false, "failed to identify verticies in tri2");
632  }
633 
634  // determine signed area of alternate config
635  // //cout<< A->GetNumCADSurf() << " " << B->GetNumCADSurf() << " "
636  // // << C->GetNumCADSurf() << " " << D->GetNumCADSurf() << endl;
637  // ofstream file;
638  // file.open("pts.3D");
639  // file << "x y z value" << endl;
640  // file << A->m_x << " " << A->m_y << " " << A->m_z << endl;
641  // file << B->m_x << " " << B->m_y << " " << B->m_z << endl;
642  // file << C->m_x << " " << C->m_y << " " << C->m_z << endl;
643  // file << D->m_x << " " << D->m_y << " " << D->m_z << endl;
644  // file.close();
645  Array<OneD, NekDouble> ai, bi, ci, di;
646  ai = A->GetCADSurfInfo(m_id);
647  bi = B->GetCADSurfInfo(m_id);
648  ci = C->GetCADSurfInfo(m_id);
649  di = D->GetCADSurfInfo(m_id);
650 
651  NekDouble CDA, CBD;
652 
653  CDA = 0.5 * (-di[0] * ci[1] + ai[0] * ci[1] + ci[0] * di[1] -
654  ai[0] * di[1] - ci[0] * ai[1] + di[0] * ai[1]);
655 
656  CBD = 0.5 * (-bi[0] * ci[1] + di[0] * ci[1] + ci[0] * bi[1] -
657  di[0] * bi[1] - ci[0] * di[1] + bi[0] * di[1]);
658 
659  // if signed area of the swapping triangles is less than zero
660  // that configuration is invalid and swap cannot be performed
661  if (!(CDA > 0.001 && CBD > 0.001))
662  {
663  m_localEdges.insert(e);
664  continue;
665  }
666 
667  bool swap = false; // assume do not swap
668 
669  if (q < 2)
670  {
671  int nodedefectbefore = 0;
672  nodedefectbefore +=
673  abs(actualConnec[A->m_id] - idealConnec[A->m_id]);
674  nodedefectbefore +=
675  abs(actualConnec[B->m_id] - idealConnec[B->m_id]);
676  nodedefectbefore +=
677  abs(actualConnec[C->m_id] - idealConnec[C->m_id]);
678  nodedefectbefore +=
679  abs(actualConnec[D->m_id] - idealConnec[D->m_id]);
680 
681  int nodedefectafter = 0;
682  nodedefectafter +=
683  abs(actualConnec[A->m_id] - 1 - idealConnec[A->m_id]);
684  nodedefectafter +=
685  abs(actualConnec[B->m_id] - 1 - idealConnec[B->m_id]);
686  nodedefectafter +=
687  abs(actualConnec[C->m_id] + 1 - idealConnec[C->m_id]);
688  nodedefectafter +=
689  abs(actualConnec[D->m_id] + 1 - idealConnec[D->m_id]);
690 
691  if (nodedefectafter < nodedefectbefore)
692  {
693  swap = true;
694  }
695  }
696  else
697  {
698  NekDouble minanglebefore = C->Angle(A, B);
699  minanglebefore = min(minanglebefore, A->Angle(B, C));
700  minanglebefore = min(minanglebefore, B->Angle(A, C));
701  minanglebefore = min(minanglebefore, B->Angle(A, D));
702  minanglebefore = min(minanglebefore, A->Angle(B, D));
703  minanglebefore = min(minanglebefore, D->Angle(A, B));
704 
705  NekDouble minangleafter = C->Angle(B, D);
706  minangleafter = min(minangleafter, D->Angle(B, C));
707  minangleafter = min(minangleafter, B->Angle(C, D));
708  minangleafter = min(minangleafter, C->Angle(A, D));
709  minangleafter = min(minangleafter, A->Angle(C, D));
710  minangleafter = min(minangleafter, D->Angle(A, C));
711 
712  if (minangleafter > minanglebefore)
713  {
714  swap = true;
715  }
716  }
717 
718  if (swap)
719  {
720  actualConnec[A->m_id]--;
721  actualConnec[B->m_id]--;
722  actualConnec[C->m_id]++;
723  actualConnec[D->m_id]++;
724 
725  // make the 4 other edges
726  EdgeSharedPtr CA, AD, DB, BC, CAt, ADt, DBt, BCt;
727  CAt = std::shared_ptr<Edge>(new Edge(C, A));
728  ADt = std::shared_ptr<Edge>(new Edge(A, D));
729  DBt = std::shared_ptr<Edge>(new Edge(D, B));
730  BCt = std::shared_ptr<Edge>(new Edge(B, C));
731 
732  vector<EdgeSharedPtr> es = tri1->GetEdgeList();
733  for (int i = 0; i < 3; i++)
734  {
735  if (es[i] == CAt)
736  {
737  CA = es[i];
738  }
739  if (es[i] == BCt)
740  {
741  BC = es[i];
742  }
743  }
744  es = tri2->GetEdgeList();
745  for (int i = 0; i < 3; i++)
746  {
747  if (es[i] == DBt)
748  {
749  DB = es[i];
750  }
751  if (es[i] == ADt)
752  {
753  AD = es[i];
754  }
755  }
756 
757  // now sort out links for the 4 edges surrounding the patch
758  vector<pair<ElementSharedPtr, int> > links;
759 
760  links = CA->m_elLink;
761  CA->m_elLink.clear();
762  for (int i = 0; i < links.size(); i++)
763  {
764  if (links[i].first->GetId() == tri1->GetId())
765  {
766  continue;
767  }
768  CA->m_elLink.push_back(links[i]);
769  }
770 
771  links = BC->m_elLink;
772  BC->m_elLink.clear();
773  for (int i = 0; i < links.size(); i++)
774  {
775  if (links[i].first->GetId() == tri1->GetId())
776  {
777  continue;
778  }
779  BC->m_elLink.push_back(links[i]);
780  }
781 
782  links = AD->m_elLink;
783  AD->m_elLink.clear();
784  for (int i = 0; i < links.size(); i++)
785  {
786  if (links[i].first->GetId() == tri2->GetId())
787  {
788  continue;
789  }
790  AD->m_elLink.push_back(links[i]);
791  }
792 
793  links = DB->m_elLink;
794  DB->m_elLink.clear();
795  for (int i = 0; i < links.size(); i++)
796  {
797  if (links[i].first->GetId() == tri2->GetId())
798  {
799  continue;
800  }
801  DB->m_elLink.push_back(links[i]);
802  }
803 
804  EdgeSharedPtr newe = std::shared_ptr<Edge>(new Edge(C, D));
805 
806  vector<NodeSharedPtr> t1, t2;
807  t1.push_back(B);
808  t1.push_back(D);
809  t1.push_back(C);
810  t2.push_back(A);
811  t2.push_back(C);
812  t2.push_back(D);
813 
814  ElmtConfig conf(LibUtilities::eTriangle, 1, false, false,
815  m_mesh->m_spaceDim != 3);
816  vector<int> tags = tri1->GetTagList();
817 
818  int id1 = tri1->GetId();
819  int id2 = tri2->GetId();
820 
822  LibUtilities::eTriangle, conf, t1, tags);
823  tags = tri2->GetTagList();
825  LibUtilities::eTriangle, conf, t2, tags);
826 
827  ntri1->SetId(id1);
828  ntri2->SetId(id2);
829  ntri1->m_parentCAD = m_cadsurf;
830  ntri2->m_parentCAD = m_cadsurf;
831 
832  vector<EdgeSharedPtr> t1es = ntri1->GetEdgeList();
833  for (int i = 0; i < 3; i++)
834  {
835  if (t1es[i] == DB)
836  {
837  ntri1->SetEdge(i, DB);
838  DB->m_elLink.push_back(
839  pair<ElementSharedPtr, int>(ntri1, i));
840  }
841  else if (t1es[i] == BC)
842  {
843  ntri1->SetEdge(i, BC);
844  BC->m_elLink.push_back(
845  pair<ElementSharedPtr, int>(ntri1, i));
846  }
847  else if (t1es[i] == newe)
848  {
849  ntri1->SetEdge(i, newe);
850  newe->m_elLink.push_back(
851  pair<ElementSharedPtr, int>(ntri1, i));
852  }
853  else
854  {
855  ASSERTL0(false, "weird edge in new tri 1");
856  }
857  }
858  vector<EdgeSharedPtr> t2es = ntri2->GetEdgeList();
859  for (int i = 0; i < 3; i++)
860  {
861  if (t2es[i] == CA)
862  {
863  ntri2->SetEdge(i, CA);
864  CA->m_elLink.push_back(
865  pair<ElementSharedPtr, int>(ntri2, i));
866  }
867  else if (t2es[i] == AD)
868  {
869  ntri2->SetEdge(i, AD);
870  AD->m_elLink.push_back(
871  pair<ElementSharedPtr, int>(ntri2, i));
872  }
873  else if (t2es[i] == newe)
874  {
875  ntri2->SetEdge(i, newe);
876  newe->m_elLink.push_back(
877  pair<ElementSharedPtr, int>(ntri2, i));
878  }
879  else
880  {
881  ASSERTL0(false, "weird edge in new tri 2");
882  }
883  }
884 
885  m_localEdges.insert(newe);
886 
887  m_localElements[id1] = ntri1;
888  m_localElements[id2] = ntri2;
889 
890  swappedEdges++;
891  }
892  else
893  {
894  m_localEdges.insert(e);
895  }
896  }
897 
898  ASSERTL0(m_localEdges.size() == edgesStart, "mismatch edge count");
899  }
900 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::unordered_set< EdgeSharedPtr, EdgeHash > EdgeSet
Definition: Edge.h:159
std::vector< ElementSharedPtr > m_localElements
local list of elements
Definition: FaceMesh.h:158
std::shared_ptr< Edge > EdgeSharedPtr
Shared pointer to an edge.
Definition: Edge.h:136
ElementFactory & GetElementFactory()
Definition: Element.cpp:44
std::shared_ptr< Node > NodeSharedPtr
Definition: CADVert.h:49
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
EdgeSet m_localEdges
local set of edges
Definition: FaceMesh.h:156
std::shared_ptr< Element > ElementSharedPtr
Definition: Edge.h:49
int m_id
id of the surface mesh
Definition: FaceMesh.h:144
double NekDouble
CADSurfSharedPtr m_cadsurf
CAD surface.
Definition: FaceMesh.h:137
MeshSharedPtr m_mesh
mesh pointer
Definition: FaceMesh.h:135
NodeSet m_localNodes
local set of nodes
Definition: FaceMesh.h:154
NodeSet m_inBoundary
set of nodes which are in the boundary (easier to identify conflicts with)
Definition: FaceMesh.h:160

◆ MakeBL()

void Nektar::NekMeshUtils::FaceMesh::MakeBL ( )
private

adds a quad layer around any interior loops

Referenced by FaceMesh().

◆ Mesh()

void Nektar::NekMeshUtils::FaceMesh::Mesh ( )

mesh exectuation command

Definition at line 145 of file FaceMesh.cpp.

References ASSERTL0.

Referenced by FaceMesh().

146 {
147  Stretching();
148  OrientateCurves();
149 
150  int numPoints = 0;
151  for (int i = 0; i < orderedLoops.size(); i++)
152  {
153  numPoints += orderedLoops[i].size();
154  for (int j = 0; j < orderedLoops[i].size(); j++)
155  {
156  m_inBoundary.insert(orderedLoops[i][j]);
157  }
158  }
159 
160  stringstream ss;
161  ss << "3 points required for triangulation, " << numPoints << " in loop"
162  << endl;
163  ss << "curves: ";
164  for (int i = 0; i < m_edgeloops.size(); i++)
165  {
166  for (int j = 0; j < m_edgeloops[i]->edges.size(); j++)
167  {
168  ss << m_edgeloops[i]->edges[j]->GetId() << " ";
169  }
170  }
171 
172  ASSERTL0(numPoints > 2, "number of verts in face is less than 3");
173 
174  // create interface to triangle thirdparty library
175  TriangleInterfaceSharedPtr pplanemesh =
177 
178  vector<Array<OneD, NekDouble> > centers;
179  for (int i = 0; i < m_edgeloops.size(); i++)
180  {
181  centers.push_back(m_edgeloops[i]->center);
182  }
183 
184  pplanemesh->Assign(orderedLoops, centers, m_id, m_str);
185 
186  pplanemesh->Mesh();
187 
188  pplanemesh->Extract(m_connec);
189 
190  bool repeat = true;
191  int meshcounter = 1;
192 
193  // continuously remesh until all triangles conform to the spacing in the
194  // octree
195  while (repeat)
196  {
197  repeat = Validate();
198  if (!repeat)
199  {
200  break;
201  }
202  m_connec.clear();
203  pplanemesh->AssignStiener(m_stienerpoints);
204  pplanemesh->Mesh();
205  pplanemesh->Extract(m_connec);
206  meshcounter++;
207  // break;
208  }
209 
210  // build a local version of the mesh (one set of triangles). this is done
211  // so edge connectivity infomration can be used for optimisation
212  BuildLocalMesh();
213 
215 
216  // make new elements and add to list from list of nodes and connectivity
217  // from triangle removing unnesercary infomration from the elements
218  for (int i = 0; i < m_localElements.size(); i++)
219  {
220  vector<EdgeSharedPtr> e = m_localElements[i]->GetEdgeList();
221  for (int j = 0; j < e.size(); j++)
222  {
223  e[j]->m_elLink.clear();
224  }
225  m_mesh->m_element[2].push_back(m_localElements[i]);
226  }
227 
228  if (m_mesh->m_verbose)
229  {
230  cout << "\r "
231  " "
232  " ";
233  cout << scientific << "\r\t\tFace " << m_id << endl
234  << "\t\t\tNodes: " << m_localNodes.size() << endl
235  << "\t\t\tEdges: " << m_localEdges.size() << endl
236  << "\t\t\tTriangles: " << m_localElements.size() << endl
237  << "\t\t\tLoops: " << m_edgeloops.size() << endl
238  << "\t\t\tSTR: " << m_str << endl
239  << endl;
240  }
241 }
std::vector< std::vector< NodeSharedPtr > > orderedLoops
list of boundary nodes in their order loops
Definition: FaceMesh.h:146
bool Validate()
Validate the surface mesh base on the octree and real dimensions of the edges.
Definition: FaceMesh.cpp:1018
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
void Stretching()
Calculate the paramter plane streching factor.
Definition: FaceMesh.cpp:964
std::vector< ElementSharedPtr > m_localElements
local list of elements
Definition: FaceMesh.h:158
std::vector< EdgeLoopSharedPtr > m_edgeloops
data structure containing the edges, their order and oreientation for the surface ...
Definition: FaceMesh.h:142
std::vector< std::vector< NodeSharedPtr > > m_connec
triangle connectiviities
Definition: FaceMesh.h:152
EdgeSet m_localEdges
local set of edges
Definition: FaceMesh.h:156
int m_id
id of the surface mesh
Definition: FaceMesh.h:144
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void OrientateCurves()
Get the boundries of the surface and extracts the nodes from the curve meshes in the correct order...
Definition: FaceMesh.cpp:1175
std::vector< NodeSharedPtr > m_stienerpoints
list of stiener points in the triangulation
Definition: FaceMesh.h:148
std::shared_ptr< TriangleInterface > TriangleInterfaceSharedPtr
NekDouble m_str
pplane stretching
Definition: FaceMesh.h:150
void BuildLocalMesh()
build a local version of mesh elements
Definition: FaceMesh.cpp:902
void OptimiseLocalMesh()
function which calls the optimisation routines
Definition: FaceMesh.cpp:243
MeshSharedPtr m_mesh
mesh pointer
Definition: FaceMesh.h:135
NodeSet m_localNodes
local set of nodes
Definition: FaceMesh.h:154
NodeSet m_inBoundary
set of nodes which are in the boundary (easier to identify conflicts with)
Definition: FaceMesh.h:160

◆ OptimiseLocalMesh()

void Nektar::NekMeshUtils::FaceMesh::OptimiseLocalMesh ( )
private

function which calls the optimisation routines

Definition at line 243 of file FaceMesh.cpp.

Referenced by FaceMesh().

244 {
245  // each optimisation algorithm is based on the work in chapter 19
246  DiagonalSwap();
247 
248  Smoothing();
249 
250  DiagonalSwap();
251 
252  Smoothing();
253 }
void DiagonalSwap()
performs diagonal swapping of edges
Definition: FaceMesh.cpp:506
void Smoothing()
performs node smoothing on face
Definition: FaceMesh.cpp:255

◆ OrientateCurves()

void Nektar::NekMeshUtils::FaceMesh::OrientateCurves ( )
private

Get the boundries of the surface and extracts the nodes from the curve meshes in the correct order.

Definition at line 1175 of file FaceMesh.cpp.

References Nektar::StdRegions::eForwards.

Referenced by FaceMesh().

1176 {
1177  // this could be a second run on orentate so clear some info
1178  orderedLoops.clear();
1179 
1180  // create list of bounding loop nodes
1181  for (int i = 0; i < m_edgeloops.size(); i++)
1182  {
1183  vector<NodeSharedPtr> cE;
1184  for (int j = 0; j < m_edgeloops[i]->edges.size(); j++)
1185  {
1186  int cid = m_edgeloops[i]->edges[j]->GetId();
1187  vector<NodeSharedPtr> edgePoints =
1188  m_curvemeshes[cid]->GetMeshPoints();
1189 
1190  int numPoints = m_curvemeshes[cid]->GetNumPoints();
1191 
1192  if (m_edgeloops[i]->edgeo[j] == CADOrientation::eForwards)
1193  {
1194  for (int k = 0; k < numPoints - 1; k++)
1195  {
1196  cE.push_back(edgePoints[k]);
1197  }
1198  }
1199  else
1200  {
1201  for (int k = numPoints - 1; k > 0; k--)
1202  {
1203  cE.push_back(edgePoints[k]);
1204  }
1205  }
1206  }
1207  orderedLoops.push_back(cE);
1208  }
1209 }
std::vector< std::vector< NodeSharedPtr > > orderedLoops
list of boundary nodes in their order loops
Definition: FaceMesh.h:146
std::map< int, CurveMeshSharedPtr > m_curvemeshes
Map of the curve meshes which bound the surfaces.
Definition: FaceMesh.h:139
std::vector< EdgeLoopSharedPtr > m_edgeloops
data structure containing the edges, their order and oreientation for the surface ...
Definition: FaceMesh.h:142

◆ Smoothing()

void Nektar::NekMeshUtils::FaceMesh::Smoothing ( )
private

performs node smoothing on face

Definition at line 255 of file FaceMesh.cpp.

References ASSERTL0, and class_topology::Node.

Referenced by FaceMesh().

256 {
257  EdgeSet::iterator eit;
258  NodeSet::iterator nit;
259 
260  Array<OneD, NekDouble> bounds = m_cadsurf->GetBounds();
261 
262  map<int, vector<EdgeSharedPtr> > connectingedges;
263 
264  map<int, vector<ElementSharedPtr> > connectingelements;
265 
266  for (eit = m_localEdges.begin(); eit != m_localEdges.end(); eit++)
267  {
268  connectingedges[(*eit)->m_n1->m_id].push_back(*eit);
269  connectingedges[(*eit)->m_n2->m_id].push_back(*eit);
270  }
271 
272  for (int i = 0; i < m_localElements.size(); i++)
273  {
274  vector<NodeSharedPtr> v = m_localElements[i]->GetVertexList();
275  for (int j = 0; j < 3; j++)
276  {
277  connectingelements[v[j]->m_id].push_back(m_localElements[i]);
278  }
279  }
280 
281  // perform 4 runs of elastic relaxation based on the octree
282  for (int q = 0; q < 4; q++)
283  {
284  for (nit = m_localNodes.begin(); nit != m_localNodes.end(); nit++)
285  {
286  NodeSet::iterator f = m_inBoundary.find((*nit));
287 
288  // node is on curve so skip
289  if (f != m_inBoundary.end())
290  {
291  continue;
292  }
293 
294  // this can be real nodes or dummy nodes depending on the system
295  vector<NodeSharedPtr> connodes;
296 
297  vector<EdgeSharedPtr> edges = connectingedges[(*nit)->m_id];
298  vector<ElementSharedPtr> els = connectingelements[(*nit)->m_id];
299 
300  vector<NodeSharedPtr> nodesystem;
301  vector<NekDouble> lamp;
302 
303  for (int i = 0; i < edges.size(); i++)
304  {
305  vector<NekDouble> lambda;
306 
307  NodeSharedPtr J;
308  if (*nit == edges[i]->m_n1)
309  {
310  J = edges[i]->m_n2;
311  }
312  else if (*nit == edges[i]->m_n2)
313  {
314  J = edges[i]->m_n1;
315  }
316  else
317  {
318  ASSERTL0(false, "could not find node");
319  }
320 
321  Array<OneD, NekDouble> ui = (*nit)->GetCADSurfInfo(m_id);
322  Array<OneD, NekDouble> uj = J->GetCADSurfInfo(m_id);
323 
324  for (int j = 0; j < els.size(); j++)
325  {
326  vector<NodeSharedPtr> v = els[j]->GetVertexList();
327 
328  // elememt is adjacent to J therefore no intersection on IJ
329  if (v[0] == J || v[1] == J || v[2] == J)
330  {
331  continue;
332  }
333 
334  // need to find other edge
335  EdgeSharedPtr AtoB;
336  bool found = false;
337  vector<EdgeSharedPtr> es = els[j]->GetEdgeList();
338  for (int k = 0; k < es.size(); k++)
339  {
340  if (!(es[k]->m_n1 == *nit || es[k]->m_n2 == *nit))
341  {
342  found = true;
343  AtoB = es[k];
344  break;
345  }
346  }
347  ASSERTL0(found, "failed to find edge to test");
348 
349  Array<OneD, NekDouble> A = AtoB->m_n1->GetCADSurfInfo(m_id);
350  Array<OneD, NekDouble> B = AtoB->m_n2->GetCADSurfInfo(m_id);
351 
352  NekDouble lam = ((A[0] - uj[0]) * (B[1] - A[1]) -
353  (A[1] - uj[1]) * (B[0] - A[0])) /
354  ((ui[0] - uj[0]) * (B[1] - A[1]) -
355  (ui[1] - uj[1]) * (B[0] - A[0]));
356 
357  if (!(lam < 0) && !(lam > 1))
358  lambda.push_back(lam);
359  }
360 
361  if (lambda.size() > 0)
362  {
363  sort(lambda.begin(), lambda.end());
364  // make a new dummy node based on the system
365  Array<OneD, NekDouble> ud(2);
366  ud[0] = uj[0] + lambda[0] * (ui[0] - uj[0]);
367  ud[1] = uj[1] + lambda[0] * (ui[1] - uj[1]);
368  Array<OneD, NekDouble> locd = m_cadsurf->P(ud);
369  NodeSharedPtr dn = std::shared_ptr<Node>(
370  new Node(0, locd[0], locd[1], locd[2]));
371  dn->SetCADSurf(m_cadsurf, ud);
372 
373  nodesystem.push_back(dn);
374  lamp.push_back(lambda[0]);
375  }
376  else
377  {
378  nodesystem.push_back(J);
379  lamp.push_back(1.0);
380  }
381  }
382 
383  Array<OneD, NekDouble> u0(2);
384  u0[0] = 0.0;
385  u0[1] = 0.0;
386 
387  for (int i = 0; i < nodesystem.size(); i++)
388  {
389  Array<OneD, NekDouble> uj = nodesystem[i]->GetCADSurfInfo(m_id);
390  u0[0] += uj[0] / nodesystem.size();
391  u0[1] += uj[1] / nodesystem.size();
392  }
393 
394  /*Array<OneD, NekDouble> pu0 = m_cadsurf->P(u0);
395  NekDouble di = m_mesh->m_octree->Query(pu0);
396  Array<OneD, NekDouble> F(2, 0.0), dF(4, 0.0);
397  for (int i = 0; i < nodesystem.size(); i++)
398  {
399  Array<OneD, NekDouble> rj = nodesystem[i]->GetLoc();
400  Array<OneD, NekDouble> uj = nodesystem[i]->GetCADSurfInfo(m_id);
401  NekDouble dj = m_mesh->m_octree->Query(rj);
402  NekDouble d = (di + dj) / 2.0;
403  NekDouble wij = sqrt((rj[0] - pu0[0]) * (rj[0] - pu0[0]) +
404  (rj[1] - pu0[1]) * (rj[1] - pu0[1]) +
405  (rj[2] - pu0[2]) * (rj[2] - pu0[2])) -
406  d;
407 
408  NekDouble umag = sqrt((uj[0] - u0[0]) * (uj[0] - u0[0]) +
409  (uj[1] - u0[1]) * (uj[1] - u0[1]));
410 
411  F[0] += wij * (uj[0] - u0[0]) / umag;
412  F[1] += wij * (uj[1] - u0[1]) / umag;
413 
414  Array<OneD, NekDouble> d1 = m_cadsurf->D1(u0);
415 
416  Array<OneD, NekDouble> dw(2, 0.0);
417  dw[0] = -2.0 *
418  ((rj[0] - pu0[0]) * d1[3] + (rj[1] - pu0[1]) * d1[4] +
419  (rj[2] - pu0[2]) * d1[5]);
420  dw[1] = -2.0 *
421  ((rj[0] - pu0[0]) * d1[6] + (rj[1] - pu0[1]) * d1[7] +
422  (rj[2] - pu0[2]) * d1[8]);
423 
424  Array<OneD, NekDouble> drhs(4, 0.0);
425  drhs[0] = 2.0 * ((uj[0] - u0[0]) * (uj[0] - u0[0]) - umag) /
426  umag / umag;
427  drhs[1] =
428  2.0 * ((uj[0] - u0[0]) * (uj[1] - u0[1])) / umag / umag;
429  drhs[2] = drhs[1];
430  drhs[3] = 2.0 * ((uj[1] - u0[1]) * (uj[1] - u0[1]) - umag) /
431  umag / umag;
432 
433  dF[0] += dw[0] * (uj[0] - u0[0]) / umag + wij * drhs[0];
434 
435  dF[1] += dw[0] * (uj[1] - u0[1]) / umag + wij * drhs[1];
436 
437  dF[2] += dw[1] * (uj[0] - u0[0]) / umag + wij * drhs[2];
438 
439  dF[3] += dw[1] * (uj[1] - u0[1]) / umag + wij * drhs[3];
440  }
441 
442  NekDouble det = dF[0] * dF[3] - dF[1] * dF[2];
443  NekDouble tmp = dF[0];
444  dF[0] = dF[3] / det;
445  dF[3] = tmp / det;
446  dF[1] *= -1.0 / det;
447  dF[2] *= -1.0 / det;
448 
449  u0[0] -= (dF[0] * F[0] + dF[2] * F[1]);
450  u0[1] -= (dF[1] * F[0] + dF[3] * F[1]);*/
451 
452  bool inbounds = true;
453  if (u0[0] < bounds[0])
454  {
455  inbounds = false;
456  }
457  else if (u0[0] > bounds[1])
458  {
459  inbounds = false;
460  }
461  else if (u0[1] < bounds[2])
462  {
463  inbounds = false;
464  }
465  else if (u0[1] > bounds[3])
466  {
467  inbounds = false;
468  }
469 
470  if (!inbounds)
471  {
472  continue;
473  }
474 
475  /*Array<OneD, NekDouble> FN(2, 0.0);
476  pu0 = m_cadsurf->P(u0);
477  di = m_mesh->m_octree->Query(pu0);
478  for (int i = 0; i < nodesystem.size(); i++)
479  {
480  Array<OneD, NekDouble> rj = nodesystem[i]->GetLoc();
481  Array<OneD, NekDouble> uj = nodesystem[i]->GetCADSurfInfo(m_id);
482  NekDouble dj = m_mesh->m_octree->Query(rj);
483  NekDouble d = (di + dj) / 2.0;
484  NekDouble wij = sqrt((rj[0] - pu0[0]) * (rj[0] - pu0[0]) +
485  (rj[1] - pu0[1]) * (rj[1] - pu0[1]) +
486  (rj[2] - pu0[2]) * (rj[2] - pu0[2])) - d;
487 
488  NekDouble umag = sqrt((uj[0] - u0[0]) * (uj[0] - u0[0]) +
489  (uj[1] - u0[1]) * (uj[1] - u0[1]));
490 
491  FN[0] += wij * (uj[0] - u0[0]) / umag;
492  FN[1] += wij * (uj[1] - u0[1]) / umag;
493  }
494 
495  if (F[0] * F[0] + F[1] * F[1] < FN[0] * FN[0] + FN[1] * FN[1])
496  {
497  continue;
498  }*/
499 
500  Array<OneD, NekDouble> l2 = m_cadsurf->P(u0);
501  (*nit)->Move(l2, m_id, u0);
502  }
503  }
504 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::vector< ElementSharedPtr > m_localElements
local list of elements
Definition: FaceMesh.h:158
std::shared_ptr< Edge > EdgeSharedPtr
Shared pointer to an edge.
Definition: Edge.h:136
std::shared_ptr< Node > NodeSharedPtr
Definition: CADVert.h:49
EdgeSet m_localEdges
local set of edges
Definition: FaceMesh.h:156
int m_id
id of the surface mesh
Definition: FaceMesh.h:144
double NekDouble
CADSurfSharedPtr m_cadsurf
CAD surface.
Definition: FaceMesh.h:137
NodeSet m_localNodes
local set of nodes
Definition: FaceMesh.h:154
NodeSet m_inBoundary
set of nodes which are in the boundary (easier to identify conflicts with)
Definition: FaceMesh.h:160

◆ Stretching()

void Nektar::NekMeshUtils::FaceMesh::Stretching ( )
private

Calculate the paramter plane streching factor.

Definition at line 964 of file FaceMesh.cpp.

Referenced by FaceMesh().

965 {
966  // define a sampling and calculate the aspect ratio of the paramter plane
967  m_str = 0.0;
968  Array<OneD, NekDouble> bnds = m_cadsurf->GetBounds();
969 
970  NekDouble dxu = int(bnds[1] - bnds[0] < bnds[3] - bnds[2]
971  ? 40
972  : (bnds[1] - bnds[0]) / (bnds[3] - bnds[2]) * 40);
973  NekDouble dxv = int(bnds[3] - bnds[2] < bnds[1] - bnds[0]
974  ? 40
975  : (bnds[3] - bnds[2]) / (bnds[1] - bnds[0]) * 40);
976 
977  NekDouble du = (bnds[1] - bnds[0]) / dxu;
978  NekDouble dv = (bnds[3] - bnds[2]) / dxv;
979 
980  int ct = 0;
981 
982  for (int i = 0; i < dxu; i++)
983  {
984  for (int j = 0; j < dxv; j++)
985  {
986  Array<OneD, NekDouble> uv(2);
987  uv[0] = bnds[0] + i * du;
988  uv[1] = bnds[2] + j * dv;
989  if (i == dxu - 1)
990  {
991  uv[0] = bnds[1];
992  }
993  if (j == dxv - 1)
994  {
995  uv[1] = bnds[3];
996  }
997  Array<OneD, NekDouble> r = m_cadsurf->D1(uv);
998 
999  NekDouble ru = sqrt(r[3] * r[3] + r[4] * r[4] + r[5] * r[5]);
1000  NekDouble rv = sqrt(r[6] * r[6] + r[7] * r[7] + r[8] * r[8]);
1001 
1002  ru *= du;
1003  rv *= dv;
1004 
1005  if (rv < 1E-8)
1006  {
1007  continue;
1008  }
1009 
1010  m_str += ru / rv;
1011  ct++;
1012  }
1013  }
1014 
1015  m_str /= ct;
1016 }
double NekDouble
NekDouble m_str
pplane stretching
Definition: FaceMesh.h:150
CADSurfSharedPtr m_cadsurf
CAD surface.
Definition: FaceMesh.h:137

◆ Validate()

bool Nektar::NekMeshUtils::FaceMesh::Validate ( )
private

Validate the surface mesh base on the octree and real dimensions of the edges.

Definition at line 1018 of file FaceMesh.cpp.

Referenced by FaceMesh().

1019 {
1020  // check all edges in the current mesh for length against the octree
1021  // if the octree is not conformed to add a new point inside the triangle
1022  // if no new points are added meshing can stop
1023  int pointBefore = m_stienerpoints.size();
1024  for (int i = 0; i < m_connec.size(); i++)
1025  {
1026  Array<OneD, NekDouble> r(3), a(3);
1027 
1028  vector<Array<OneD, NekDouble> > info;
1029 
1030  for (int j = 0; j < 3; j++)
1031  {
1032  info.push_back(m_connec[i][j]->GetCADSurfInfo(m_id));
1033  }
1034 
1035  r[0] = m_connec[i][0]->Distance(m_connec[i][1]);
1036  r[1] = m_connec[i][1]->Distance(m_connec[i][2]);
1037  r[2] = m_connec[i][2]->Distance(m_connec[i][0]);
1038 
1039  a[0] = m_connec[i][0]->Angle(m_connec[i][1]->GetLoc(),
1040  m_connec[i][2]->GetLoc(),
1041  m_cadsurf->N(info[0]));
1042  a[1] = m_connec[i][1]->Angle(m_connec[i][2]->GetLoc(),
1043  m_connec[i][0]->GetLoc(),
1044  m_cadsurf->N(info[1]));
1045  a[2] = m_connec[i][2]->Angle(m_connec[i][0]->GetLoc(),
1046  m_connec[i][1]->GetLoc(),
1047  m_cadsurf->N(info[2]));
1048 
1049  NekDouble d1 = m_mesh->m_octree->Query(m_connec[i][0]->GetLoc());
1050  NekDouble d2 = m_mesh->m_octree->Query(m_connec[i][1]->GetLoc());
1051  NekDouble d3 = m_mesh->m_octree->Query(m_connec[i][2]->GetLoc());
1052 
1053  Array<OneD, NekDouble> uvc(2);
1054  uvc[0] = (info[0][0] + info[1][0] + info[2][0]) / 3.0;
1055  uvc[1] = (info[0][1] + info[1][1] + info[2][1]) / 3.0;
1056 
1057  Array<OneD, NekDouble> locc = m_cadsurf->P(uvc);
1058  NekDouble d4 = m_mesh->m_octree->Query(locc);
1059 
1060  NekDouble d = (d1 + d2 + d3 + d4) / 4.0;
1061 
1062  vector<bool> valid(3);
1063  valid[0] = r[0] < d * 1.5;
1064  valid[1] = r[1] < d * 1.5;
1065  valid[2] = r[2] < d * 1.5;
1066 
1067  vector<bool> angValid(3);
1068  angValid[0] = a[0] / M_PI * 180.0 > 20.0 && a[0] / M_PI * 180.0 < 120.0;
1069  angValid[1] = a[1] / M_PI * 180.0 > 20.0 && a[1] / M_PI * 180.0 < 120.0;
1070  angValid[2] = a[2] / M_PI * 180.0 > 20.0 && a[2] / M_PI * 180.0 < 120.0;
1071 
1072  int numValid = 0;
1073  int numAngValid = 0;
1074  for (int j = 0; j < 3; j++)
1075  {
1076  if (valid[j])
1077  {
1078  numValid++;
1079  }
1080  if (angValid[j])
1081  {
1082  numAngValid++;
1083  }
1084  }
1085 
1086  // if numvalid is zero no work to be done
1087  /*if (numValid != 3)
1088  {
1089  AddNewPoint(uvc);
1090  }*/
1091 
1092  if (numValid != 3 || numAngValid != 3)
1093  {
1094  // break the bad edge
1095  /*int a=0, b=0;
1096  if(!valid[0])
1097  {
1098  a = 0;
1099  b = 1;
1100  }
1101  else if(!valid[1])
1102  {
1103  a = 1;
1104  b = 2;
1105  }
1106  else if(!valid[2])
1107  {
1108  a = 2;
1109  b = 0;
1110  }
1111 
1112  Array<OneD, NekDouble> uvn(2);
1113  uvn[0] = (info[a][0] + info[b][0]) / 2.0;
1114  uvn[1] = (info[a][1] + info[b][1]) / 2.0;*/
1115  AddNewPoint(uvc);
1116  }
1117  }
1118 
1119  if (m_stienerpoints.size() == pointBefore)
1120  {
1121  return false;
1122  }
1123  else
1124  {
1125  return true;
1126  }
1127 }
std::vector< std::vector< NodeSharedPtr > > m_connec
triangle connectiviities
Definition: FaceMesh.h:152
int m_id
id of the surface mesh
Definition: FaceMesh.h:144
double NekDouble
std::vector< NodeSharedPtr > m_stienerpoints
list of stiener points in the triangulation
Definition: FaceMesh.h:148
void AddNewPoint(Array< OneD, NekDouble > uv)
adds a new stiener point to the triangulation for meshing
Definition: FaceMesh.cpp:1129
CADSurfSharedPtr m_cadsurf
CAD surface.
Definition: FaceMesh.h:137
MeshSharedPtr m_mesh
mesh pointer
Definition: FaceMesh.h:135

◆ ValidateCurves()

bool Nektar::NekMeshUtils::FaceMesh::ValidateCurves ( )

validate the curve meshes

Definition at line 46 of file FaceMesh.cpp.

References CG_Iterations::loc.

Referenced by FaceMesh().

47 {
48  vector<int> curvesInSurface;
49  for (int i = 0; i < m_edgeloops.size(); i++)
50  {
51  for (int j = 0; j < m_edgeloops[i]->edges.size(); j++)
52  {
53  curvesInSurface.push_back(m_edgeloops[i]->edges[j]->GetId());
54  }
55  }
56 
57  bool error = false;
58 
59  for (int i = 0; i < curvesInSurface.size(); i++)
60  {
61  vector<EdgeSharedPtr> es =
62  m_curvemeshes[curvesInSurface[i]]->GetMeshEdges();
63 
64  for (int j = i; j < curvesInSurface.size(); j++)
65  {
66  if (i == j)
67  {
68  continue;
69  }
70 
71  vector<EdgeSharedPtr> es2 =
72  m_curvemeshes[curvesInSurface[j]]->GetMeshEdges();
73 
74  for (int l = 0; l < es.size(); l++)
75  {
76  Array<OneD, NekDouble> P1 = es[l]->m_n1->GetCADSurfInfo(m_id);
77  Array<OneD, NekDouble> P2 = es[l]->m_n2->GetCADSurfInfo(m_id);
78  for (int k = 0; k < es2.size(); k++)
79  {
80  if (es[l]->m_n1 == es2[k]->m_n1 ||
81  es[l]->m_n1 == es2[k]->m_n2 ||
82  es[l]->m_n2 == es2[k]->m_n1 ||
83  es[l]->m_n2 == es2[k]->m_n2)
84  {
85  continue;
86  }
87 
88  Array<OneD, NekDouble> P3 =
89  es2[k]->m_n1->GetCADSurfInfo(m_id);
90  Array<OneD, NekDouble> P4 =
91  es2[k]->m_n2->GetCADSurfInfo(m_id);
92 
93  NekDouble den = (P4[0] - P3[0]) * (P2[1] - P1[1]) -
94  (P2[0] - P1[0]) * (P4[1] - P3[1]);
95 
96  if (fabs(den) < 1e-8)
97  {
98  continue;
99  }
100 
101  NekDouble t = ((P1[0] - P3[0]) * (P4[1] - P3[1]) -
102  (P4[0] - P3[0]) * (P1[1] - P3[1])) /
103  den;
104  NekDouble u =
105  (P1[0] - P3[0] + t * (P2[0] - P1[0])) / (P4[0] - P3[0]);
106 
107  if (t < 1.0 && t > 0.0 && u < 1.0 && u > 0.0)
108  {
109  Array<OneD, NekDouble> uv(2);
110  uv[0] = P1[0] + t * (P2[0] - P1[0]);
111  uv[1] = P1[1] + t * (P2[1] - P1[1]);
112  Array<OneD, NekDouble> loc = m_cadsurf->P(uv);
113  cout << endl
114  << "Curve mesh error at " << loc[0] << " "
115  << loc[1] << " " << loc[2] << " on face " << m_id
116  << endl;
117  error = true;
118  }
119  }
120  }
121  }
122  }
123  return error;
124 }
std::map< int, CurveMeshSharedPtr > m_curvemeshes
Map of the curve meshes which bound the surfaces.
Definition: FaceMesh.h:139
std::vector< EdgeLoopSharedPtr > m_edgeloops
data structure containing the edges, their order and oreientation for the surface ...
Definition: FaceMesh.h:142
int m_id
id of the surface mesh
Definition: FaceMesh.h:144
double NekDouble
CADSurfSharedPtr m_cadsurf
CAD surface.
Definition: FaceMesh.h:137

◆ ValidateLoops()

void Nektar::NekMeshUtils::FaceMesh::ValidateLoops ( )

validate the curve meshes considering the loops

Definition at line 126 of file FaceMesh.cpp.

Referenced by FaceMesh().

127 {
128  OrientateCurves();
129 
130  for (int i = 0; i < orderedLoops.size(); i++)
131  {
132  int numPoints = orderedLoops[i].size();
133  if(numPoints == 2)
134  {
135  //force a remesh of the curves
136  for(int j = 0; j < m_edgeloops[i]->edges.size(); j++)
137  {
138  int cid = m_edgeloops[i]->edges[j]->GetId();
139  m_curvemeshes[cid]->ReMesh();
140  }
141  }
142  }
143 }
std::vector< std::vector< NodeSharedPtr > > orderedLoops
list of boundary nodes in their order loops
Definition: FaceMesh.h:146
std::map< int, CurveMeshSharedPtr > m_curvemeshes
Map of the curve meshes which bound the surfaces.
Definition: FaceMesh.h:139
std::vector< EdgeLoopSharedPtr > m_edgeloops
data structure containing the edges, their order and oreientation for the surface ...
Definition: FaceMesh.h:142
void OrientateCurves()
Get the boundries of the surface and extracts the nodes from the curve meshes in the correct order...
Definition: FaceMesh.cpp:1175

Friends And Related Function Documentation

◆ MemoryManager< FaceMesh >

friend class MemoryManager< FaceMesh >
friend

Definition at line 54 of file FaceMesh.h.

Member Data Documentation

◆ m_cadsurf

CADSurfSharedPtr Nektar::NekMeshUtils::FaceMesh::m_cadsurf
private

CAD surface.

Definition at line 137 of file FaceMesh.h.

Referenced by FaceMesh().

◆ m_compId

int Nektar::NekMeshUtils::FaceMesh::m_compId
private

identity to put into element tags

Definition at line 162 of file FaceMesh.h.

◆ m_connec

std::vector<std::vector<NodeSharedPtr> > Nektar::NekMeshUtils::FaceMesh::m_connec
private

triangle connectiviities

Definition at line 152 of file FaceMesh.h.

◆ m_curvemeshes

std::map<int, CurveMeshSharedPtr> Nektar::NekMeshUtils::FaceMesh::m_curvemeshes
private

Map of the curve meshes which bound the surfaces.

Definition at line 139 of file FaceMesh.h.

◆ m_edgeloops

std::vector<EdgeLoopSharedPtr> Nektar::NekMeshUtils::FaceMesh::m_edgeloops
private

data structure containing the edges, their order and oreientation for the surface

Definition at line 142 of file FaceMesh.h.

Referenced by FaceMesh().

◆ m_id

int Nektar::NekMeshUtils::FaceMesh::m_id
private

id of the surface mesh

Definition at line 144 of file FaceMesh.h.

Referenced by FaceMesh().

◆ m_inBoundary

NodeSet Nektar::NekMeshUtils::FaceMesh::m_inBoundary
private

set of nodes which are in the boundary (easier to identify conflicts with)

Definition at line 160 of file FaceMesh.h.

◆ m_localEdges

EdgeSet Nektar::NekMeshUtils::FaceMesh::m_localEdges
private

local set of edges

Definition at line 156 of file FaceMesh.h.

◆ m_localElements

std::vector<ElementSharedPtr> Nektar::NekMeshUtils::FaceMesh::m_localElements
private

local list of elements

Definition at line 158 of file FaceMesh.h.

◆ m_localNodes

NodeSet Nektar::NekMeshUtils::FaceMesh::m_localNodes
private

local set of nodes

Definition at line 154 of file FaceMesh.h.

◆ m_mesh

MeshSharedPtr Nektar::NekMeshUtils::FaceMesh::m_mesh
private

mesh pointer

Definition at line 135 of file FaceMesh.h.

Referenced by FaceMesh().

◆ m_stienerpoints

std::vector<NodeSharedPtr> Nektar::NekMeshUtils::FaceMesh::m_stienerpoints
private

list of stiener points in the triangulation

Definition at line 148 of file FaceMesh.h.

◆ m_str

NekDouble Nektar::NekMeshUtils::FaceMesh::m_str
private

pplane stretching

Definition at line 150 of file FaceMesh.h.

◆ orderedLoops

std::vector<std::vector<NodeSharedPtr> > Nektar::NekMeshUtils::FaceMesh::orderedLoops
private

list of boundary nodes in their order loops

Definition at line 146 of file FaceMesh.h.