Nektar++
PreProcessing.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: ProcessJac.h
4 //
5 // For more information, please see: http://www.nektar.info/
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Calculate jacobians of elements.
32 //
33 ////////////////////////////////////////////////////////////////////////////////
34 
35 #include "ProcessVarOpti.h"
36 
40 
41 using namespace std;
42 
43 namespace Nektar
44 {
45 namespace Utilities
46 {
47 
48 map<LibUtilities::ShapeType, DerivUtilSharedPtr> ProcessVarOpti::BuildDerivUtil(
49  int o)
50 {
51  // build Vandermonde information
52  map<LibUtilities::ShapeType, DerivUtilSharedPtr> ret;
53 
54  // Typedef for points types used in the variational optimser. First entry is
55  // the evaluation points; second entry is the distributions used for the
56  // full element (includes element boundary)
57  typedef std::pair<LibUtilities::PointsType, LibUtilities::PointsType>
58  PTypes;
59 
60  map<LibUtilities::ShapeType, PTypes> typeMap;
61 
62  if (m_mesh->m_nummode + o <= 11)
63  {
64  typeMap[LibUtilities::eTriangle] =
68  typeMap[LibUtilities::ePrism] =
70  }
71 
74  // typeMap[LibUtilities::eHexahedron] =
75  // PTypes(LibUtilities::eNodalHexElec, LibUtilities::eNodalHexElec);
76 
77  for (auto &it : typeMap)
78  {
79  PTypes pType = it.second;
80  DerivUtilSharedPtr der = std::shared_ptr<DerivUtil>(new DerivUtil());
81 
82  LibUtilities::PointsKey pkey1(m_mesh->m_nummode, pType.second);
83  LibUtilities::PointsKey pkey2(m_mesh->m_nummode + o, pType.first);
84 
85  const int pDim = pkey1.GetPointsDim();
86  const int order = m_mesh->m_nummode - 1;
87 
88  Array<OneD, Array<OneD, NekDouble> > u1(pDim), u2(pDim);
89 
90  switch (pDim)
91  {
92  case 2:
93  {
94  LibUtilities::PointsManager()[pkey1]->GetPoints(u1[0], u1[1]);
95  LibUtilities::PointsManager()[pkey2]->GetPoints(u2[0], u2[1]);
96  break;
97  }
98  case 3:
99  {
100  LibUtilities::PointsManager()[pkey1]->GetPoints(u1[0], u1[1],
101  u1[2]);
102  LibUtilities::PointsManager()[pkey2]->GetPoints(u2[0], u2[1],
103  u2[2]);
104  break;
105  }
106  }
107 
108  der->ptsStd = u1[0].num_elements();
109  der->pts = u2[0].num_elements();
110 
111  LibUtilities::NodalUtil *nodalUtil = NULL;
112 
113  if (it.first == LibUtilities::eTriangle)
114  {
115  nodalUtil =
116  new LibUtilities::NodalUtilTriangle(order, u1[0], u1[1]);
117  }
118  else if (it.first == LibUtilities::eQuadrilateral)
119  {
120  nodalUtil = new LibUtilities::NodalUtilQuad(order, u1[0], u1[1]);
121  }
122  else if (it.first == LibUtilities::eTetrahedron)
123  {
124  nodalUtil = new LibUtilities::NodalUtilTetrahedron(order, u1[0],
125  u1[1], u1[2]);
126  }
127  else if (it.first == LibUtilities::ePrism)
128  {
129  nodalUtil =
130  new LibUtilities::NodalUtilPrism(order, u1[0], u1[1], u1[2]);
131  }
132  else if (it.first == LibUtilities::eHexahedron)
133  {
134  nodalUtil =
135  new LibUtilities::NodalUtilHex(order, u1[0], u1[1], u1[2]);
136  }
137  else
138  {
139  ASSERTL0(false,
140  "Unknown element type for derivative utility setup");
141  }
142 
143  NekMatrix<NekDouble> interp = *nodalUtil->GetInterpolationMatrix(u2);
144  NekMatrix<NekDouble> Vandermonde = *nodalUtil->GetVandermonde();
145  NekMatrix<NekDouble> VandermondeI = Vandermonde;
146  VandermondeI.Invert();
147 
148  for (int i = 0; i < pDim; ++i)
149  {
150  der->VdmDStd[i] =
151  *nodalUtil->GetVandermondeForDeriv(i) * VandermondeI;
152  der->VdmD[i] = interp * der->VdmDStd[i];
153  }
154 
156  LibUtilities::PointsManager()[pkey2]->GetW();
157  NekVector<NekDouble> quadWi(qds);
158  der->quadW = quadWi;
159 
160  ret[it.first] = der;
161  //delete nodalUtil;
162  }
163 
164  return ret;
165 }
166 
167 
168 
170 {
171  const vector<int> & value_vector;
172 
173  NodeComparator(const vector<int> & val_vec):
174  value_vector(val_vec) {}
175 
176  bool operator()(int i1, int i2)
177  {
178  return value_vector[i1] > value_vector[i2];
179  }
180 };
181 
182 
183 
184 vector<vector<NodeSharedPtr> > ProcessVarOpti::GetColouredNodes(
185  vector<ElUtilSharedPtr> elLock)
186 {
187 
188  // create set of nodes to be ignored and hence not included in the coloursets
189  NodeSet ignoredNodes;
190  for (int i = 0; i < elLock.size(); i++)
191  {
192  vector<NodeSharedPtr> nodes;
193  elLock[i]->GetEl()->GetCurvedNodes(nodes);
194  for (int j = 0; j < nodes.size(); j++)
195  {
196  ignoredNodes.insert(nodes[j]);
197  }
198  }
199 
200  // create set of nodes which are at the boundary and hence not included in the colourset
201  NodeSet boundaryNodes;
202 
203  switch (m_mesh->m_spaceDim)
204  {
205  case 2:
206  {
207  for(auto &edge : m_mesh->m_edgeSet)
208  {
209  if(edge->m_elLink.size() == 2)
210  {
211  continue;
212  }
213 
214  boundaryNodes.insert(edge->m_n1);
215  boundaryNodes.insert(edge->m_n2);
216  for(int i = 0; i < edge->m_edgeNodes.size(); i++)
217  {
218  boundaryNodes.insert(edge->m_edgeNodes[i]);
219  }
220  }
221  break;
222  }
223  case 3:
224  {
225  if(!m_mesh->m_cad)
226  {
227  for (auto &face : m_mesh->m_faceSet)
228  {
229  if (face->m_elLink.size() == 2)
230  {
231  continue;
232  }
233 
234  vector<NodeSharedPtr> vs = face->m_vertexList;
235  for (int j = 0; j < vs.size(); j++)
236  {
237  boundaryNodes.insert(vs[j]);
238  }
239 
240  vector<EdgeSharedPtr> es = face->m_edgeList;
241  for (int j = 0; j < es.size(); j++)
242  {
243  for (int k = 0; k < es[j]->m_edgeNodes.size(); k++)
244  {
245  boundaryNodes.insert(es[j]->m_edgeNodes[k]);
246  }
247  }
248 
249  for (int i = 0; i < face->m_faceNodes.size(); i++)
250  {
251  boundaryNodes.insert(face->m_faceNodes[i]);
252  }
253  }
254  }
255  else
256  {
257  //if we have CAD therefore the only fixed nodes exist on vertices only
258  for (auto &node : m_mesh->m_vertexSet)
259  {
260  if(node->GetNumCadCurve() > 1)
261  {
262  boundaryNodes.insert(node);
263  }
264  }
265  }
266  break;
267  }
268  default:
269  ASSERTL0(false,"space dim issue");
270  }
271 
272 
273  //create vector of free nodes which "remain", hence will be included in the coloursets
274  vector<NodeSharedPtr> remainEdgeVertex;
275  vector<NodeSharedPtr> remainFace;
276  vector<NodeSharedPtr> remainVolume;
277  m_res->nDoF = 0;
278 
279  // check if vertex nodes are in boundary or ignored nodes, otherwise add to EDGE-VERTEX remain nodes
280  for (auto &node : m_mesh->m_vertexSet)
281  {
282  if (boundaryNodes.find(node) == boundaryNodes.end() &&
283  ignoredNodes.find(node) == ignoredNodes.end())
284  {
285  remainEdgeVertex.push_back(node);
286  if (node->GetNumCadCurve() == 1)
287  {
288  m_res->nDoF++;
289  }
290  else if (node->GetNumCADSurf() == 1)
291  {
292  m_res->nDoF += 2;
293  }
294  else
295  {
296  m_res->nDoF += m_mesh->m_spaceDim;
297  }
298  }
299  }
300 
301  // check if edge nodes are in boundary or ignored nodes, otherwise add to EDGE-VERTEX remain nodes
302  for (auto &edge : m_mesh->m_edgeSet)
303  {
304  vector<NodeSharedPtr> &n = edge->m_edgeNodes;
305  for (int j = 0; j < n.size(); j++)
306  {
307  if (boundaryNodes.find(n[j]) == boundaryNodes.end() &&
308  ignoredNodes.find(n[j]) == ignoredNodes.end())
309  {
310  remainEdgeVertex.push_back(n[j]);
311  if (n[j]->GetNumCadCurve() == 1)
312  {
313  m_res->nDoF++;
314  }
315  else if (n[j]->GetNumCADSurf() == 1)
316  {
317  m_res->nDoF += 2;
318  }
319  else
320  {
321  m_res->nDoF += m_mesh->m_spaceDim;
322  }
323  }
324  }
325  }
326 
327  // check if face nodes are in boundary or ignored nodes, otherwise add to FACE remain nodes
328  for (auto &face : m_mesh->m_faceSet)
329  {
330  vector<NodeSharedPtr> &n = face->m_faceNodes;
331  for (int j = 0; j < n.size(); j++)
332  {
333  if (boundaryNodes.find(n[j]) == boundaryNodes.end() &&
334  ignoredNodes.find(n[j]) == ignoredNodes.end())
335  {
336  remainFace.push_back(n[j]);
337  if (n[j]->GetNumCADSurf() == 1)
338  {
339  m_res->nDoF += 2;
340  }
341  else
342  {
343  m_res->nDoF += m_mesh->m_spaceDim;
344  }
345  }
346  }
347  }
348 
349  // check if volume nodes are in boundary or ignored nodes, otherwise add to VOLUME remain nodes
350  for (int i = 0; i < m_mesh->m_element[m_mesh->m_expDim].size(); i++)
351  {
352  vector<NodeSharedPtr> ns =
353  m_mesh->m_element[m_mesh->m_expDim][i]->GetVolumeNodes();
354  for (int j = 0; j < ns.size(); j++)
355  {
356  if (boundaryNodes.find(ns[j]) == boundaryNodes.end() &&
357  ignoredNodes.find(ns[j]) == ignoredNodes.end())
358  {
359  remainVolume.push_back(ns[j]);
360  m_res->nDoF += m_mesh->m_spaceDim;
361  }
362  }
363  }
364 
365  // size of all free nodes to be included in the coloursets
366  m_res->n = remainEdgeVertex.size()
367  + remainFace.size() + remainVolume.size();
368 
369  // data structure for coloursets, that will ultimately contain all free nodes
370  vector<vector<NodeSharedPtr> > ret;
371  vector<vector<NodeSharedPtr> > retPart;
372 
373  // edge and vertex nodes
374  // create vector num_el of number of associated elements of each node
375  vector<int> num_el(remainEdgeVertex.size());
376  for (int i = 0; i < remainEdgeVertex.size(); i++)
377  {
378  //try to find node within all elements
379  auto it = m_nodeElMap.find(remainEdgeVertex[i]->m_id);
380  vector<ElUtilSharedPtr> &elUtils = it->second;
381  num_el[i] = elUtils.size();
382  }
383  // finding the permutation according to num_el
384  vector<int> permNode(remainEdgeVertex.size());
385  for (int i = 0; i < remainEdgeVertex.size(); ++i)
386  {
387  permNode[i] = i;
388  }
389  std::sort(permNode.begin(), permNode.end(), NodeComparator(num_el));
390  // applying the permutation to remainEdgeVertex
391  vector<NodeSharedPtr> remainEdgeVertexSort(remainEdgeVertex.size());
392  for (int i = 0; i < remainEdgeVertex.size(); ++i)
393  {
394  int j = permNode[i];
395  remainEdgeVertexSort[i] = remainEdgeVertex[j];
396  }
397 
398  retPart = CreateColoursets(remainEdgeVertexSort);
399  if(m_mesh->m_verbose)
400  {
401  cout << "Number of Edge/Vertex Coloursets: " << retPart.size() << endl;
402  }
403  for (int i = 0; i < retPart.size(); i++)
404  {
405  ret.push_back(retPart[i]);
406  }
407 
408  // face nodes
409  retPart = CreateColoursets(remainFace);
410  if(m_mesh->m_verbose)
411  {
412  cout << "Number of Face Coloursets: " << retPart.size() << endl;
413  }
414  for (int i = 0; i < retPart.size(); i++)
415  {
416  ret.push_back(retPart[i]);
417  }
418 
419  // volume nodes
420  retPart = CreateColoursets(remainVolume);
421  if(m_mesh->m_verbose)
422  {
423  cout << "Number of Volume Coloursets: " << retPart.size() << endl;
424  }
425  for (int i = 0; i < retPart.size(); i++)
426  {
427  ret.push_back(retPart[i]);
428  }
429 
430 
431  if(m_mesh->m_verbose)
432  {
433  cout << endl;
434  }
435 
436  return ret;
437 }
438 
439 vector<vector<NodeSharedPtr> > ProcessVarOpti::CreateColoursets(
440  vector<NodeSharedPtr> remain)
441 {
442  vector<vector<NodeSharedPtr> > retPart;
443 
444  // loop until all free nodes have been sorted
445  while (remain.size() > 0)
446  {
447  vector<NodeSharedPtr> layer; // one colourset
448  set<int> locked;
449  set<int> completed;
450  for (int i = 0; i < remain.size(); i++)
451  {
452  // Try to find node within all elements
453  auto it = m_nodeElMap.find(remain[i]->m_id);
454  ASSERTL0(it != m_nodeElMap.end(), "could not find node");
455 
456  // identify the vector of all associated elements of the node
457  vector<ElUtilSharedPtr> &elUtils = it->second;
458 
459  // suppose node is not locked
460  bool islocked = false;
461 
462  // loop over all associated elements of the node
463  for (int j = 0; j < elUtils.size(); j++)
464  {
465  // check all nodes of the element. if node is within the set of
466  // locked nodes then lock node and go to the next node
467  if (locked.find(elUtils[j]->GetId()) != locked.end())
468  {
469  islocked = true;
470  break;
471  }
472  }
473 
474  // if the node is not locked, insert it into the colourset and
475  // insert sorted node into the completed list. Then, loop over all
476  // other nodes of the same element and mark them as locked.
477  if (!islocked)
478  {
479  layer.push_back(remain[i]);
480  completed.insert(remain[i]->m_id);
481  for (int j = 0; j < elUtils.size(); j++)
482  {
483  locked.insert(elUtils[j]->GetId());
484  }
485  }
486  }
487 
488  // identify nodes which are not sorted, yet and create new "remain"
489  // vector
490  vector<NodeSharedPtr> tmp = remain;
491  remain.clear();
492  for (int i = 0; i < tmp.size(); i++)
493  {
494  if (completed.find(tmp[i]->m_id) == completed.end())
495  {
496  remain.push_back(tmp[i]);
497  }
498  }
499 
500  // include layer or colourset into vector of coloursets
501  retPart.push_back(layer);
502 
503  // print out progress
504  if(m_mesh->m_verbose)
505  {
507  m_res->n - remain.size(), m_res->n, "Node Coloring");
508  }
509 
510  }
511  return retPart;
512 }
513 
514 
515 void ProcessVarOpti::GetElementMap(
516  int o, map<LibUtilities::ShapeType, DerivUtilSharedPtr> derMap)
517 {
518  for (int i = 0; i < m_mesh->m_element[m_mesh->m_expDim].size(); i++)
519  {
520  ElementSharedPtr el = m_mesh->m_element[m_mesh->m_expDim][i];
521  vector<NodeSharedPtr> ns;
522  el->GetCurvedNodes(ns);
523  ElUtilSharedPtr d = std::shared_ptr<ElUtil>(new ElUtil(
524  el, derMap[el->GetShapeType()], m_res, m_mesh->m_nummode, o));
525  m_dataSet.push_back(d);
526  }
527 
528  for (int i = 0; i < m_mesh->m_element[m_mesh->m_expDim].size(); i++)
529  {
530  ElementSharedPtr el = m_mesh->m_element[m_mesh->m_expDim][i];
531  vector<NodeSharedPtr> ns;
532  el->GetCurvedNodes(ns);
533 
534  for (int j = 0; j < ns.size(); j++)
535  {
536  m_nodeElMap[ns[j]->m_id].push_back(m_dataSet[i]);
537  }
538 
539  ASSERTL0(derMap[el->GetShapeType()]->ptsStd == ns.size(),
540  "mismatch node count");
541  }
542 }
543 
544 vector<ElUtilSharedPtr> ProcessVarOpti::GetLockedElements(NekDouble thres)
545 {
546  vector<ElUtilSharedPtr> elBelowThres;
547  for (int i = 0; i < m_dataSet.size(); ++i)
548  {
549  if (m_dataSet[i]->GetScaledJac() < thres)
550  {
551  elBelowThres.push_back(m_dataSet[i]);
552  }
553  }
554 
555  std::unordered_set<int> inmesh;
556  vector<ElUtilSharedPtr> totest;
557 
558  for (int i = 0; i < elBelowThres.size(); i++)
559  {
560  auto t = inmesh.insert(elBelowThres[i]->GetId());
561 
562  vector<FaceSharedPtr> f = elBelowThres[i]->GetEl()->GetFaceList();
563  for (int j = 0; j < f.size(); j++)
564  {
565  for (int k = 0; k < f[j]->m_elLink.size(); k++)
566  {
567  if (f[j]->m_elLink[k].first->GetId() ==
568  elBelowThres[i]->GetId())
569  {
570  continue;
571  }
572 
573  t = inmesh.insert(f[j]->m_elLink[k].first->GetId());
574  if (t.second)
575  {
576  totest.push_back(
577  m_dataSet[f[j]->m_elLink[k].first->GetId()]);
578  }
579  }
580  }
581  }
582 
583  for (int i = 0; i < 6; i++)
584  {
585  vector<ElUtilSharedPtr> tmp = totest;
586  totest.clear();
587  for (int j = 0; j < tmp.size(); j++)
588  {
589  vector<FaceSharedPtr> f = tmp[j]->GetEl()->GetFaceList();
590  for (int k = 0; k < f.size(); k++)
591  {
592  for (int l = 0; l < f[k]->m_elLink.size(); l++)
593  {
594  if (f[k]->m_elLink[l].first->GetId() == tmp[j]->GetId())
595  {
596  continue;
597  }
598 
599  auto t = inmesh.insert(f[k]->m_elLink[l].first->GetId());
600  if (t.second)
601  {
602  totest.push_back(
603  m_dataSet[f[k]->m_elLink[l].first->GetId()]);
604  }
605  }
606  }
607  }
608  }
609 
610  // now need to invert the list
611  vector<ElUtilSharedPtr> ret;
612  for (int i = 0; i < m_dataSet.size(); ++i)
613  {
614  if (inmesh.find(m_dataSet[i]->GetId()) == inmesh.end())
615  {
616  ret.push_back(m_dataSet[i]);
617  }
618  }
619 
620  return ret;
621 }
622 
623 void ProcessVarOpti::RemoveLinearCurvature()
624 {
625  for(int i = 0; i < m_dataSet.size(); i++)
626  {
627  if(m_dataSet[i]->GetScaledJac() > 0.999)
628  {
629  ElementSharedPtr el = m_dataSet[i]->GetEl();
630  vector<NodeSharedPtr> ns;
631  el->SetVolumeNodes(ns);
632  }
633  }
634 
635  map<int, vector<FaceSharedPtr> > edgeToFace;
636 
637  for(auto &face : m_mesh->m_faceSet)
638  {
639  bool rm = true;
640  for(int i = 0; i < face->m_elLink.size(); i++)
641  {
642  int id = face->m_elLink[i].first->GetId();
643  if(m_dataSet[id]->GetScaledJac() <= 0.999)
644  {
645  rm = false;
646  break;
647  }
648  }
649  if(rm)
650  {
651  face->m_faceNodes.clear();
652  }
653 
654  vector<EdgeSharedPtr> es = face->m_edgeList;
655  for(int i = 0; i < es.size(); i++)
656  {
657  edgeToFace[es[i]->m_id].push_back(face);
658  }
659  }
660 
661  for(auto &edge : m_mesh->m_edgeSet)
662  {
663  auto it = edgeToFace.find(edge->m_id);
664  ASSERTL0(it != edgeToFace.end(),"not found");
665  bool rm = true;
666  for(int i = 0; i < it->second.size(); i++)
667  {
668  if(it->second[i]->m_faceNodes.size() > 0)
669  {
670  rm = false;
671  break;
672  }
673  }
674  if(rm)
675  {
676  edge->m_edgeNodes.clear();
677  }
678  }
679 }
680 
681 }
682 }
std::shared_ptr< DerivUtil > DerivUtilSharedPtr
Definition: ElUtil.h:50
#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
A class to assist in the construction of nodal simplex and hybrid elements in two and three dimension...
Definition: NodalUtil.h:83
Specialisation of the NodalUtil class to support nodal triangular elements.
Definition: NodalUtil.h:169
STL namespace.
3D Nodal Symmetric positive internal tet (Whitherden, Vincent)
Definition: PointsType.h:77
SharedMatrix GetVandermondeForDeriv(int dir)
Return the Vandermonde matrix of the derivative of the basis functions for the nodal distribution...
Definition: NodalUtil.cpp:134
std::unordered_set< NodeSharedPtr, NodeHash > NodeSet
Definition: Node.h:447
SharedMatrix GetInterpolationMatrix(Array< OneD, Array< OneD, NekDouble > > &xi)
Construct the interpolation matrix used to evaluate the basis at the points xi inside the element...
Definition: NodalUtil.cpp:201
const vector< int > & value_vector
Specialisation of the NodalUtil class to support nodal hex elements.
Definition: NodalUtil.h:352
std::shared_ptr< Element > ElementSharedPtr
Definition: Edge.h:49
3D Nodal Electrostatic Points on a Tetrahedron
Definition: PointsType.h:73
PointsManagerT & PointsManager(void)
Defines a specification for a set of points.
Definition: Points.h:59
double NekDouble
Specialisation of the NodalUtil class to support nodal prismatic elements.
Definition: NodalUtil.h:263
SharedMatrix GetVandermonde()
Return the Vandermonde matrix for the nodal distribution.
Definition: NodalUtil.cpp:101
NodeComparator(const vector< int > &val_vec)
Specialisation of the NodalUtil class to support nodal quad elements.
Definition: NodalUtil.h:311
std::shared_ptr< ElUtil > ElUtilSharedPtr
Definition: ElUtil.h:113
2D Nodal Electrostatic Points on a Triangle
Definition: PointsType.h:69
3D electrostatically spaced points on a Prism
Definition: PointsType.h:75
Specialisation of the NodalUtil class to support nodal tetrahedral elements.
Definition: NodalUtil.h:214
2D Nodal Symmetric positive internal triangle (Whitherden, Vincent)
Definition: PointsType.h:76