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

#include <NodeOptiCAD.h>

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

Public Member Functions

 NodeOpti1D3D (NodeSharedPtr n, std::vector< ElUtilSharedPtr > e, ResidualSharedPtr r, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > d, optiType o, CADCurveSharedPtr c)
 
 ~NodeOpti1D3D ()
 
void Optimise ()
 
- Public Member Functions inherited from Nektar::Utilities::NodeOpti
 NodeOpti (NodeSharedPtr n, std::vector< ElUtilSharedPtr > e, ResidualSharedPtr r, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > d, optiType o, int dim)
 
virtual ~NodeOpti ()
 
void CalcMinJac ()
 
NodeOptiJobGetJob ()
 
template<int DIM>
NekDouble GetFunctional (NekDouble &minJacNew, bool gradient=true)
 Evaluate functional for elements connected to a node. More...
 
template<int DIM>
void MinEigen (NekDouble &val)
 Calculates minimum eigenvalue of Hessian matrix. More...
 
template<>
void MinEigen (NekDouble &val)
 
template<>
void MinEigen (NekDouble &val)
 

Static Public Member Functions

static NodeOptiSharedPtr create (NodeSharedPtr n, std::vector< ElUtilSharedPtr > e, ResidualSharedPtr r, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > d, optiType o)
 

Static Public Attributes

static int m_type
 

Private Member Functions

void ProcessGradient ()
 

Private Attributes

CADCurveSharedPtr curve
 

Additional Inherited Members

- Protected Member Functions inherited from Nektar::Utilities::NodeOpti
template<int DIM>
int IsIndefinite ()
 Returns 1 if Hessian matrix is indefinite and 0 otherwise. More...
 
template<>
int IsIndefinite ()
 
template<>
int IsIndefinite ()
 
- Static Protected Member Functions inherited from Nektar::Utilities::NodeOpti
static NekDouble c1 ()
 
static NekDouble gradTol ()
 
static NekDouble alphaTol ()
 
- Protected Attributes inherited from Nektar::Utilities::NodeOpti
NodeSharedPtr m_node
 
std::mutex mtx
 
std::map< LibUtilities::ShapeType, std::vector< ElUtilSharedPtr > > m_data
 
std::vector< NekDoublem_grad
 
std::vector< NekDoublem_tmpStore
 
std::unordered_map< LibUtilities::ShapeType, DerivArray, EnumHashm_derivs
 
NekDouble m_minJac
 
ResidualSharedPtr m_res
 
std::map< LibUtilities::ShapeType, DerivUtilSharedPtrm_derivUtils
 
optiType m_opti
 

Detailed Description

Definition at line 47 of file NodeOptiCAD.h.

Constructor & Destructor Documentation

◆ NodeOpti1D3D()

Nektar::Utilities::NodeOpti1D3D::NodeOpti1D3D ( NodeSharedPtr  n,
std::vector< ElUtilSharedPtr e,
ResidualSharedPtr  r,
std::map< LibUtilities::ShapeType, DerivUtilSharedPtr d,
optiType  o,
CADCurveSharedPtr  c 
)
inline

Definition at line 50 of file NodeOptiCAD.h.

Referenced by create().

54  : NodeOpti(n, e, r, d, o, 3), curve(c)
55  {
56  }
NodeOpti(NodeSharedPtr n, std::vector< ElUtilSharedPtr > e, ResidualSharedPtr r, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > d, optiType o, int dim)
Definition: NodeOpti.h:60

◆ ~NodeOpti1D3D()

Nektar::Utilities::NodeOpti1D3D::~NodeOpti1D3D ( )
inline

Definition at line 58 of file NodeOptiCAD.h.

References Optimise().

58 {};

Member Function Documentation

◆ create()

static NodeOptiSharedPtr Nektar::Utilities::NodeOpti1D3D::create ( NodeSharedPtr  n,
std::vector< ElUtilSharedPtr e,
ResidualSharedPtr  r,
std::map< LibUtilities::ShapeType, DerivUtilSharedPtr d,
optiType  o 
)
inlinestatic

Definition at line 63 of file NodeOptiCAD.h.

References NodeOpti1D3D(), and ProcessGradient().

66  {
67  std::vector<CADCurveSharedPtr> cs = n->GetCADCurves();
68  return NodeOptiSharedPtr(new NodeOpti1D3D(n, e, r, d, o, cs[0]));
69  }
NodeOpti1D3D(NodeSharedPtr n, std::vector< ElUtilSharedPtr > e, ResidualSharedPtr r, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > d, optiType o, CADCurveSharedPtr c)
Definition: NodeOptiCAD.h:50
std::shared_ptr< NodeOpti > NodeOptiSharedPtr
Definition: NodeOpti.h:135

◆ Optimise()

void Nektar::Utilities::NodeOpti1D3D::Optimise ( )
virtual

Implements Nektar::Utilities::NodeOpti.

Definition at line 51 of file NodeOptiCAD.cpp.

References Nektar::Utilities::GetNodeOptiFactory(), and Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::RegisterCreatorFunction().

Referenced by ~NodeOpti1D3D(), and Nektar::Utilities::NodeOpti2D3D::~NodeOpti2D3D().

52 {
53  NekDouble minJacNew;
54  NekDouble currentW = GetFunctional<3>(minJacNew);
55  NekDouble newVal = currentW;
56 
57  if (m_grad[0] * m_grad[0] + m_grad[1] * m_grad[1] + m_grad[2] * m_grad[2] >
58  gradTol())
59  {
60  // modify the gradient to be on the cad system
62 
63  // needs to optimise
64  NekDouble tc = m_node->GetCADCurveInfo(curve->GetId());
65  NekDouble xc = m_node->m_x;
66  NekDouble yc = m_node->m_y;
67  NekDouble zc = m_node->m_z;
68  NekDouble nt;
69 
70  vector<NekDouble> sk(1);
71 
72  if (m_grad[1] < 1e-6)
73  {
74  m_grad[1] = 1e-6 - m_grad[1];
75  }
76 
77  sk[0] = m_grad[0] / m_grad[1] * -1.0;
78 
79  vector<NekDouble> bd(2);
80  curve->GetBounds(bd[0], bd[1]);
81 
82  bool found = false;
83 
84  NekDouble pg = m_grad[0] * sk[0];
85 
86  // normal gradient line Search
87  NekDouble alpha = 1.0;
88 
89  while (alpha > alphaTol())
90  {
91  // Update node
92  nt = tc + alpha * sk[0];
93  if (nt < bd[0] || nt > bd[1])
94  {
95  alpha /= 2.0;
96  continue;
97  }
98 
99  curve->P(nt, m_node->m_x, m_node->m_y, m_node->m_z);
100 
101  newVal = GetFunctional<3>(minJacNew, false);
102 
103  // Wolfe conditions
104  if (newVal <= currentW + c1() * alpha * pg)
105  {
106  found = true;
107  break;
108  }
109 
110  alpha /= 2.0;
111  }
112 
113  if (!found)
114  {
115  // reset the node
116  nt = tc;
117  curve->P(nt, m_node->m_x, m_node->m_y, m_node->m_z);
118 
119  mtx.lock();
120  m_res->nReset[0]++;
121  mtx.unlock();
122  }
123  else
124  {
125  m_minJac = minJacNew;
126  m_node->Move(m_node->m_x, m_node->m_y, m_node->m_z, curve->GetId(),
127  nt);
128  }
129  mtx.lock();
130  m_res->val = max(sqrt((m_node->m_x - xc) * (m_node->m_x - xc) +
131  (m_node->m_y - yc) * (m_node->m_y - yc) +
132  (m_node->m_z - zc) * (m_node->m_z - zc)),
133  m_res->val);
134  mtx.unlock();
135  }
136 
137  mtx.lock();
138  m_res->func += newVal;
139  mtx.unlock();
140 }
static NekDouble alphaTol()
Definition: NodeOpti.h:129
static NekDouble c1()
Definition: NodeOpti.h:121
double NekDouble
ResidualSharedPtr m_res
Definition: NodeOpti.h:117
std::vector< NekDouble > m_grad
Definition: NodeOpti.h:109
static NekDouble gradTol()
Definition: NodeOpti.h:125

◆ ProcessGradient()

void Nektar::Utilities::NodeOpti1D3D::ProcessGradient ( )
private

Definition at line 252 of file NodeOptiCAD.cpp.

Referenced by create(), and Nektar::Utilities::NodeOpti2D3D::create().

253 {
254  NekDouble tc = m_node->GetCADCurveInfo(curve->GetId());
255  vector<NekDouble> grad = m_grad;
256  m_grad = vector<NekDouble>(2, 0.0);
257 
258  // Grab first and second order CAD derivatives
259  Array<OneD, NekDouble> d2 = curve->D2(tc);
260 
261  // Multiply gradient by derivative of CAD
262  m_grad[0] = grad[0] * d2[3] + grad[1] * d2[4] + grad[2] * d2[5];
263 
264  // Second order: product rule of above, so multiply gradient by second
265  // order CAD derivatives and Hessian by gradient of CAD
266  m_grad[1] = grad[0] * d2[6] + grad[1] * d2[7] + grad[2] * d2[8] +
267  d2[3] * (grad[3] * d2[3] + grad[4] * d2[4] + grad[5] * d2[5]) +
268  d2[4] * (grad[4] * d2[3] + grad[6] * d2[4] + grad[7] * d2[5]) +
269  d2[5] * (grad[5] * d2[3] + grad[7] * d2[4] + grad[8] * d2[5]);
270 }
double NekDouble
std::vector< NekDouble > m_grad
Definition: NodeOpti.h:109

Member Data Documentation

◆ curve

CADCurveSharedPtr Nektar::Utilities::NodeOpti1D3D::curve
private

Definition at line 73 of file NodeOptiCAD.h.

◆ m_type

int Nektar::Utilities::NodeOpti1D3D::m_type
static
Initial value:

Definition at line 62 of file NodeOptiCAD.h.