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

#include <ProcessVarOpti.h>

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

Public Member Functions

 ProcessVarOpti (MeshSharedPtr m)
 
virtual ~ProcessVarOpti ()
 
virtual void Process ()
 
- 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 (MeshSharedPtr m)
 Creates an instance of this class. More...
 

Static Public Attributes

static ModuleKey className
 

Private Types

typedef std::map< int, std::vector< ElUtilSharedPtr > > NodeElMap
 

Private Member Functions

void Analytics ()
 
std::map< LibUtilities::ShapeType, DerivUtilSharedPtrBuildDerivUtil (int o)
 
void GetElementMap (int o, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > derMap)
 
std::vector< ElUtilSharedPtrGetLockedElements (NekDouble thres)
 
std::vector< std::vector< NodeSharedPtr > > CreateColoursets (std::vector< NodeSharedPtr > remain)
 
std::vector< std::vector< NodeSharedPtr > > GetColouredNodes (std::vector< ElUtilSharedPtr > elLock)
 
void RemoveLinearCurvature ()
 

Private Attributes

NodeElMap m_nodeElMap
 
std::vector< ElUtilSharedPtrm_dataSet
 
ResidualSharedPtr m_res
 
optiType m_opti
 

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 82 of file ProcessVarOpti.h.

Member Typedef Documentation

◆ NodeElMap

typedef std::map<int, std::vector<ElUtilSharedPtr> > Nektar::Utilities::ProcessVarOpti::NodeElMap
private

Definition at line 100 of file ProcessVarOpti.h.

Constructor & Destructor Documentation

◆ ProcessVarOpti()

Nektar::Utilities::ProcessVarOpti::ProcessVarOpti ( MeshSharedPtr  m)

Definition at line 72 of file ProcessVarOpti.cpp.

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

72  : ProcessModule(m)
73 {
74  // clang-format off
75  m_config["linearelastic"] =
76  ConfigOption(true, "", "Optimise for linear elasticity");
77  m_config["winslow"] =
78  ConfigOption(true, "", "Optimise for winslow");
79  m_config["roca"] =
80  ConfigOption(true, "", "Optimise for roca method");
81  m_config["hyperelastic"] =
82  ConfigOption(true, "", "Optimise for hyper elasticity");
83  m_config["numthreads"] =
84  ConfigOption(false, "1", "Number of threads");
85  m_config["restol"] =
86  ConfigOption(false, "1e-6", "Tolerance criterion");
87  m_config["maxiter"] =
88  ConfigOption(false, "500", "Maximum number of iterations");
89  m_config["nq"] =
90  ConfigOption(false, "-1", "Order of mesh");
91  m_config["region"] =
92  ConfigOption(false, "0.0", "create regions based on target");
93  m_config["resfile"] =
94  ConfigOption(false, "", "writes residual values to file");
95  m_config["histfile"] =
96  ConfigOption(false, "", "histogram of scaled jac");
97  m_config["overint"] =
98  ConfigOption(false, "6", "over integration order");
99  m_config["analytics"] =
100  ConfigOption(false, "", "basic analytics module");
101  // clang-format on
102 }
Represents a command-line configuration option.
NEKMESHUTILS_EXPORT ProcessModule(MeshSharedPtr p_m)
std::map< std::string, ConfigOption > m_config
List of configuration values.

◆ ~ProcessVarOpti()

Nektar::Utilities::ProcessVarOpti::~ProcessVarOpti ( )
virtual

Definition at line 104 of file ProcessVarOpti.cpp.

105 {
106 }

Member Function Documentation

◆ Analytics()

void Nektar::Utilities::ProcessVarOpti::Analytics ( )
private

Definition at line 451 of file ProcessVarOpti.cpp.

References BuildDerivUtil(), Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), GetElementMap(), Nektar::Utilities::GetNodeOptiFactory(), m_dataSet, Nektar::NekMeshUtils::Module::m_mesh, m_nodeElMap, m_opti, and m_res.

Referenced by Process().

452 {
453  // Grab the first element from the list
454  ElementSharedPtr elmt = m_mesh->m_element[m_mesh->m_expDim][0];
455 
456  // Get curved nodes
457  vector<NodeSharedPtr> nodes;
458  elmt->GetCurvedNodes(nodes);
459 
460  // We're going to investigate only the first node (corner node)
461  NodeSharedPtr node = nodes[4];
462 
463  // Loop over overintegration orders
464  const int nPoints = 200;
465  const int overInt = 40;
466  const NekDouble originX = -1.0;
467  const NekDouble originY = -1.0;
468  const NekDouble length = 2.0;
469  const NekDouble dx = length / (nPoints - 1);
470 
471  cout << "# overint = " << overInt << endl;
472  cout << "# Columns: x, y, over-integration orders (0 -> " << overInt - 1
473  << "), "
474  << " min(scaledJac)" << endl;
475 
476  // Loop over square defined by (originX, originY), length
477  for (int k = 0; k < nPoints; ++k)
478  {
479  node->m_y = originY + k * dx;
480  for (int j = 0; j < nPoints; ++j)
481  {
482  node->m_x = originX + j * dx;
483  cout << node->m_x << " " << node->m_y << " ";
484 
485  NekDouble minJacNew;
486 
487  for (int i = 0; i < overInt; ++i)
488  {
489  // Clear any existing node to element mapping.
490  m_dataSet.clear();
491  m_nodeElMap.clear();
492 
493  // Build deriv utils and element map.
494  map<LibUtilities::ShapeType, DerivUtilSharedPtr> derivUtils =
495  BuildDerivUtil(i);
496 
497  // Reconstruct element map
498  GetElementMap(i, derivUtils);
499 
500  for (int j = 0; j < m_dataSet.size(); j++)
501  {
502  m_dataSet[j]->Evaluate();
503  m_dataSet[j]->InitialMinJac();
504  }
505 
506  // Create NodeOpti object.
507  NodeOptiSharedPtr nodeOpti =
509  m_mesh->m_spaceDim * 11, node,
510  m_nodeElMap.find(node->m_id)->second, m_res, derivUtils,
511  m_opti);
512 
513  minJacNew = 0.0;
514 
515  // Evaluate functional.
516  nodeOpti->CalcMinJac();
517  cout << nodeOpti->GetFunctional<2>(minJacNew) << " ";
518  // NekDouble eigen;
519  // nodeOpti->GetFunctional<2>(minJacNew);
520  // nodeOpti->MinEigen<2>(eigen);
521  // cout << eigen << " ";
522  }
523 
524  cout << minJacNew << endl;
525  }
526  }
527 }
NodeOptiFactory & GetNodeOptiFactory()
Definition: NodeOpti.cpp:51
void GetElementMap(int o, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > derMap)
std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > BuildDerivUtil(int o)
std::vector< ElUtilSharedPtr > m_dataSet
std::shared_ptr< NodeOpti > NodeOptiSharedPtr
Definition: NodeOpti.h:135
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
std::shared_ptr< Element > ElementSharedPtr
Definition: Edge.h:49
double NekDouble

◆ BuildDerivUtil()

map< LibUtilities::ShapeType, DerivUtilSharedPtr > Nektar::Utilities::ProcessVarOpti::BuildDerivUtil ( int  o)
private

Definition at line 48 of file PreProcessing.cpp.

References ASSERTL0, Nektar::LibUtilities::eHexahedron, Nektar::LibUtilities::eNodalPrismElec, Nektar::LibUtilities::eNodalPrismSPI, Nektar::LibUtilities::eNodalQuadElec, Nektar::LibUtilities::eNodalTetElec, Nektar::LibUtilities::eNodalTetSPI, Nektar::LibUtilities::eNodalTriElec, Nektar::LibUtilities::eNodalTriSPI, Nektar::LibUtilities::ePrism, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eTetrahedron, Nektar::LibUtilities::eTriangle, Nektar::LibUtilities::NodalUtil::GetInterpolationMatrix(), Nektar::LibUtilities::NodalUtil::GetVandermonde(), Nektar::LibUtilities::NodalUtil::GetVandermondeForDeriv(), and Nektar::LibUtilities::PointsManager().

Referenced by Analytics(), and Process().

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 
155  Array<OneD, NekDouble> qds =
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 }
std::shared_ptr< DerivUtil > DerivUtilSharedPtr
Definition: ElUtil.h:50
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
3D Nodal Symmetric positive internal tet (Whitherden, Vincent)
Definition: PointsType.h:77
3D Nodal Electrostatic Points on a Tetrahedron
Definition: PointsType.h:73
PointsManagerT & PointsManager(void)
2D Nodal Electrostatic Points on a Triangle
Definition: PointsType.h:69
3D electrostatically spaced points on a Prism
Definition: PointsType.h:75
2D Nodal Symmetric positive internal triangle (Whitherden, Vincent)
Definition: PointsType.h:76

◆ create()

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

Creates an instance of this class.

Definition at line 86 of file ProcessVarOpti.h.

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

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

◆ CreateColoursets()

vector< vector< NodeSharedPtr > > Nektar::Utilities::ProcessVarOpti::CreateColoursets ( std::vector< NodeSharedPtr remain)
private

Definition at line 439 of file PreProcessing.cpp.

References ASSERTL0, and Nektar::LibUtilities::PrintProgressbar().

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 }
#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

◆ GetColouredNodes()

vector< vector< NodeSharedPtr > > Nektar::Utilities::ProcessVarOpti::GetColouredNodes ( std::vector< ElUtilSharedPtr elLock)
private

Definition at line 184 of file PreProcessing.cpp.

References ASSERTL0.

Referenced by Process().

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 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::vector< std::vector< NodeSharedPtr > > CreateColoursets(std::vector< NodeSharedPtr > remain)
std::unordered_set< NodeSharedPtr, NodeHash > NodeSet
Definition: Node.h:447

◆ GetElementMap()

void Nektar::Utilities::ProcessVarOpti::GetElementMap ( int  o,
std::map< LibUtilities::ShapeType, DerivUtilSharedPtr derMap 
)
private

Definition at line 515 of file PreProcessing.cpp.

References ASSERTL0.

Referenced by Analytics(), and Process().

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 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::vector< ElUtilSharedPtr > m_dataSet
std::shared_ptr< Element > ElementSharedPtr
Definition: Edge.h:49
std::shared_ptr< ElUtil > ElUtilSharedPtr
Definition: ElUtil.h:113

◆ GetLockedElements()

vector< ElUtilSharedPtr > Nektar::Utilities::ProcessVarOpti::GetLockedElements ( NekDouble  thres)
private

Definition at line 544 of file PreProcessing.cpp.

Referenced by Process().

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 }
std::vector< ElUtilSharedPtr > m_dataSet

◆ Process()

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

Implements Nektar::NekMeshUtils::Module.

Definition at line 108 of file ProcessVarOpti.cpp.

References Analytics(), ASSERTL0, BuildDerivUtil(), Nektar::Thread::ThreadMaster::CreateInstance(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::Utilities::eHypEl, Nektar::Utilities::eLinEl, Nektar::Utilities::eRoca, Nektar::Utilities::eWins, GetColouredNodes(), GetElementMap(), GetLockedElements(), Nektar::Utilities::GetNodeOptiFactory(), Nektar::NekMeshUtils::Module::m_config, m_dataSet, Nektar::NekMeshUtils::Module::m_mesh, m_nodeElMap, m_opti, m_res, CellMLToNektar.pycml::name, CellMLToNektar.cellml_metadata::p, Nektar::Thread::ThreadMaster::SessionJob, Nektar::Thread::ThreadMaster::SetThreadingType(), Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and Nektar::LibUtilities::Timer::TimePerTest().

109 {
110  if (m_mesh->m_verbose)
111  {
112  cout << "ProcessVarOpti: Optimising... " << endl;
113  }
114 
115  if (m_config["linearelastic"].beenSet)
116  {
117  m_opti = eLinEl;
118  }
119  else if (m_config["winslow"].beenSet)
120  {
121  m_opti = eWins;
122  }
123  else if (m_config["roca"].beenSet)
124  {
125  m_opti = eRoca;
126  }
127  else if (m_config["hyperelastic"].beenSet)
128  {
129  m_opti = eHypEl;
130  }
131  else
132  {
133  ASSERTL0(false, "not opti type set");
134  }
135 
136  const int maxIter = m_config["maxiter"].as<int>();
137  const NekDouble restol = m_config["restol"].as<NekDouble>();
138 
139  // m_mesh->m_nummode = m_config["nq"].as<int>();
140 
141  bool fd = false;
142 
143  if (m_config["nq"].beenSet)
144  {
145  m_mesh->m_nummode = m_config["nq"].as<int>();
146  fd = true;
147  }
148 
149  if (!fd)
150  {
151  for (auto &edge : m_mesh->m_edgeSet)
152  {
153  if (edge->m_edgeNodes.size() > 0)
154  {
155  m_mesh->m_nummode = edge->m_edgeNodes.size() + 2;
156  fd = true;
157  break;
158  }
159  }
160  }
161  ASSERTL0(fd, "failed to find order of mesh");
162 
163  // Safety feature: limit over-integration order for high-order triangles
164  // over order 5.
165  int intOrder = m_config["overint"].as<int>();
166  intOrder = m_mesh->m_nummode + intOrder <= 11 ?
167  intOrder : 11 - m_mesh->m_nummode;
168 
169  if (m_mesh->m_verbose)
170  {
171  cout << "Identified mesh order as: " << m_mesh->m_nummode - 1 << endl;
172  }
173 
174  if (m_mesh->m_expDim == 2 && m_mesh->m_spaceDim == 3)
175  {
176  ASSERTL0(false, "cannot deal with manifolds");
177  }
178 
179  m_res = std::shared_ptr<Residual>(new Residual);
180  m_res->val = 1.0;
181 
182 
183  m_mesh->MakeOrder(m_mesh->m_nummode - 1,
185 
186  if (m_config["analytics"].beenSet)
187  {
188  Analytics();
189  return;
190  }
191 
192  map<LibUtilities::ShapeType, DerivUtilSharedPtr> derivUtils =
193  BuildDerivUtil(intOrder);
194 
195  GetElementMap(intOrder, derivUtils);
196 
197  m_res->startInv = 0;
198  m_res->worstJac = numeric_limits<double>::max();
199  for (int i = 0; i < m_dataSet.size(); i++)
200  {
201  m_dataSet[i]->Evaluate();
202  m_dataSet[i]->InitialMinJac();
203  }
204 
205  vector<ElUtilSharedPtr> elLock;
206 
207  if (m_config["region"].beenSet)
208  {
209  elLock = GetLockedElements(m_config["region"].as<NekDouble>());
210  }
211 
212 
213  vector<vector<NodeSharedPtr> > freenodes = GetColouredNodes(elLock);
214  vector<vector<NodeOptiSharedPtr> > optiNodes;
215 
216  // turn the free nodes into optimisable objects with all required data
217  set<int> check;
218  for (int i = 0; i < freenodes.size(); i++)
219  {
220  vector<NodeOptiSharedPtr> ns;
221  for (int j = 0; j < freenodes[i].size(); j++)
222  {
223  auto it = m_nodeElMap.find(freenodes[i][j]->m_id);
224  ASSERTL0(it != m_nodeElMap.end(), "could not find");
225 
226  int optiKind = m_mesh->m_spaceDim;
227 
228  if (freenodes[i][j]->GetNumCadCurve() == 1)
229  {
230  optiKind += 10;
231  }
232  else if (freenodes[i][j]->GetNumCADSurf() == 1)
233  {
234  optiKind += 20;
235  }
236  else
237  {
238  optiKind += 10 * m_mesh->m_expDim;
239  }
240 
241  auto c = check.find(freenodes[i][j]->m_id);
242  ASSERTL0(c == check.end(), "duplicate node");
243  check.insert(freenodes[i][j]->m_id);
244 
245  ns.push_back(GetNodeOptiFactory().CreateInstance(
246  optiKind, freenodes[i][j], it->second, m_res, derivUtils,
247  m_opti));
248  }
249  optiNodes.push_back(ns);
250  }
251 
252  int nset = optiNodes.size();
253  int p = 0;
254  int mn = numeric_limits<int>::max();
255  int mx = 0;
256  for (int i = 0; i < nset; i++)
257  {
258  p += optiNodes[i].size();
259  mn = min(mn, int(optiNodes[i].size()));
260  mx = max(mx, int(optiNodes[i].size()));
261  }
262 
263  if (m_config["histfile"].beenSet)
264  {
265  ofstream histFile;
266  string name = m_config["histfile"].as<string>() + "_start.txt";
267  histFile.open(name.c_str());
268 
269  for (int i = 0; i < m_dataSet.size(); i++)
270  {
271  histFile << m_dataSet[i]->GetScaledJac() << endl;
272  }
273  histFile.close();
274  }
275 
276  if(m_mesh->m_verbose)
277  {
278  cout << scientific << endl;
279  cout << "N elements:\t\t" << m_mesh->m_element[m_mesh->m_expDim].size() - elLock.size() << endl
280  << "N elements invalid:\t" << m_res->startInv << endl
281  << "Worst jacobian:\t\t" << m_res->worstJac << endl
282  << "N free nodes:\t\t" << m_res->n << endl
283  << "N Dof:\t\t\t" << m_res->nDoF << endl
284  << "N color sets:\t\t" << nset << endl
285  << "Avg set colors:\t\t" << p/nset << endl
286  << "Min set:\t\t" << mn << endl
287  << "Max set:\t\t" << mx << endl
288  << "Residual tolerance:\t" << restol << endl;
289  }
290 
291  int nThreads = m_config["numthreads"].as<int>();
292 
293  int ctr = 0;
294  Thread::ThreadMaster tms;
295  tms.SetThreadingType("ThreadManagerBoost");
297  tms.CreateInstance(Thread::ThreadMaster::SessionJob, nThreads);
298 
299  LibUtilities::Timer t;
300  t.Start();
301 
302  ofstream resFile;
303  if (m_config["resfile"].beenSet)
304  {
305  resFile.open(m_config["resfile"].as<string>().c_str());
306  }
307 
308  for (int i = 0; i < optiNodes.size(); i++)
309  {
310  vector<Thread::ThreadJob *> jobs(optiNodes[i].size());
311  for (int j = 0; j < optiNodes[i].size(); j++)
312  {
313  optiNodes[i][j]->CalcMinJac();
314  }
315  }
316 
317  while (m_res->val > restol)
318  {
319  ctr++;
320  m_res->val = 0.0;
321  m_res->func = 0.0;
322  m_res->nReset[0] = 0;
323  m_res->nReset[1] = 0;
324  m_res->nReset[2] = 0;
325  m_res->alphaI = 0;
326  for (int i = 0; i < optiNodes.size(); i++)
327  {
328  vector<Thread::ThreadJob *> jobs(optiNodes[i].size());
329  for (int j = 0; j < optiNodes[i].size(); j++)
330  {
331  jobs[j] = optiNodes[i][j]->GetJob();
332  }
333 
334  tm->SetNumWorkers(0);
335  tm->QueueJobs(jobs);
336  tm->SetNumWorkers(nThreads);
337  tm->Wait();
338  }
339 
340  m_res->startInv = 0;
341  m_res->worstJac = numeric_limits<double>::max();
342 
343  vector<Thread::ThreadJob *> elJobs(m_dataSet.size());
344  for (int i = 0; i < m_dataSet.size(); i++)
345  {
346  elJobs[i] = m_dataSet[i]->GetJob();
347  }
348 
349  tm->SetNumWorkers(0);
350  tm->QueueJobs(elJobs);
351  tm->SetNumWorkers(nThreads);
352  tm->Wait();
353 
354  if (m_config["resfile"].beenSet)
355  {
356  resFile << m_res->val << " " << m_res->worstJac << " "
357  << m_res->func << endl;
358  }
359 
360  if(m_mesh->m_verbose)
361  {
362  cout << ctr
363  << "\tResidual: " << m_res->val
364  << "\tMin Jac: " << m_res->worstJac
365  << "\tInvalid: " << m_res->startInv
366  << "\tReset nodes: " << m_res->nReset[0] << "/" << m_res->nReset[1]
367  << "/" << m_res->nReset[2] << "\tFunctional: " << m_res->func
368  << endl;
369  }
370 
371  if(ctr >= maxIter)
372  {
373  break;
374  }
375  }
376 
377  if (m_config["histfile"].beenSet)
378  {
379  ofstream histFile;
380  string name = m_config["histfile"].as<string>() + "_end.txt";
381  histFile.open(name.c_str());
382 
383  for (int i = 0; i < m_dataSet.size(); i++)
384  {
385  histFile << m_dataSet[i]->GetScaledJac() << endl;
386  }
387  histFile.close();
388  }
389  if (m_config["resfile"].beenSet)
390  {
391  resFile.close();
392  }
393 
394  t.Stop();
395 
396  //RemoveLinearCurvature();
397 
398  if(m_mesh->m_verbose)
399  {
400  cout << "Time to compute: " << t.TimePerTest(1) << endl;
401  cout << "Invalid at end: " << m_res->startInv << endl;
402  cout << "Worst at end: " << m_res->worstJac << endl;
403  }
404 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
NodeOptiFactory & GetNodeOptiFactory()
Definition: NodeOpti.cpp:51
std::shared_ptr< ThreadManager > ThreadManagerSharedPtr
Definition: Thread.h:72
std::vector< ElUtilSharedPtr > GetLockedElements(NekDouble thres)
void GetElementMap(int o, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > derMap)
std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > BuildDerivUtil(int o)
std::vector< ElUtilSharedPtr > m_dataSet
std::vector< std::vector< NodeSharedPtr > > GetColouredNodes(std::vector< ElUtilSharedPtr > elLock)
double NekDouble
std::map< std::string, ConfigOption > m_config
List of configuration values.
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51

◆ RemoveLinearCurvature()

void Nektar::Utilities::ProcessVarOpti::RemoveLinearCurvature ( )
private

Definition at line 623 of file PreProcessing.cpp.

References ASSERTL0.

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 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::vector< ElUtilSharedPtr > m_dataSet
std::shared_ptr< Element > ElementSharedPtr
Definition: Edge.h:49

Member Data Documentation

◆ className

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

Definition at line 90 of file ProcessVarOpti.h.

◆ m_dataSet

std::vector<ElUtilSharedPtr> Nektar::Utilities::ProcessVarOpti::m_dataSet
private

Definition at line 114 of file ProcessVarOpti.h.

Referenced by Analytics(), and Process().

◆ m_nodeElMap

NodeElMap Nektar::Utilities::ProcessVarOpti::m_nodeElMap
private

Definition at line 113 of file ProcessVarOpti.h.

Referenced by Analytics(), and Process().

◆ m_opti

optiType Nektar::Utilities::ProcessVarOpti::m_opti
private

Definition at line 117 of file ProcessVarOpti.h.

Referenced by Analytics(), and Process().

◆ m_res

ResidualSharedPtr Nektar::Utilities::ProcessVarOpti::m_res
private

Definition at line 116 of file ProcessVarOpti.h.

Referenced by Analytics(), and Process().