Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::Utilities::NodeOpti1D3D:
Collaboration graph
[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
 
boost::mutex mtx
 
std::map
< LibUtilities::ShapeType,
std::vector< ElUtilSharedPtr > > 
m_data
 
Array< OneD, NekDoublem_grad
 
std::vector< NekDoublem_tmpStore
 
boost::unordered_map
< LibUtilities::ShapeType,
DerivArray
m_derivs
 
NekDouble m_minJac
 
ResidualSharedPtr m_res
 
std::map
< LibUtilities::ShapeType,
DerivUtilSharedPtr
m_derivUtils
 
optiType m_opti
 

Detailed Description

Definition at line 48 of file NodeOptiCAD.h.

Constructor & Destructor Documentation

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 51 of file NodeOptiCAD.h.

Referenced by create().

55  : NodeOpti(n, e, r, d, o, 3), curve(c)
56  {
57  }
NodeOpti(NodeSharedPtr n, std::vector< ElUtilSharedPtr > e, ResidualSharedPtr r, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > d, optiType o, int dim)
Definition: NodeOpti.h:59
Nektar::Utilities::NodeOpti1D3D::~NodeOpti1D3D ( )
inline

Definition at line 59 of file NodeOptiCAD.h.

59 {};

Member Function Documentation

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 64 of file NodeOptiCAD.h.

References NodeOpti1D3D().

67  {
68  std::vector<std::pair<int, CADCurveSharedPtr> > cs = n->GetCADCurves();
69  return NodeOptiSharedPtr(new NodeOpti1D3D(n, e, r, d, o, cs[0].second));
70  }
NodeOpti1D3D(NodeSharedPtr n, std::vector< ElUtilSharedPtr > e, ResidualSharedPtr r, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > d, optiType o, CADCurveSharedPtr c)
Definition: NodeOptiCAD.h:51
boost::shared_ptr< NodeOpti > NodeOptiSharedPtr
Definition: NodeOpti.h:135
void Nektar::Utilities::NodeOpti1D3D::Optimise ( )
virtual

Implements Nektar::Utilities::NodeOpti.

Definition at line 52 of file NodeOptiCAD.cpp.

References CellMLToNektar.cellml_metadata::p.

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

Definition at line 263 of file NodeOptiCAD.cpp.

264 {
265  NekDouble tc = m_node->GetCADCurveInfo(curve->GetId());
266  Array<OneD, NekDouble> grad = m_grad;
267  m_grad = Array<OneD, NekDouble>(2, 0.0);
268 
269  // Grab first and second order CAD derivatives
270  Array<OneD, NekDouble> d2 = curve->D2(tc);
271 
272  // Multiply gradient by derivative of CAD
273  m_grad[0] = grad[0] * d2[3] + grad[1] * d2[4] + grad[2] * d2[5];
274 
275  // Second order: product rule of above, so multiply gradient by second
276  // order CAD derivatives and Hessian by gradient of CAD
277  m_grad[1] = grad[0] * d2[6] + grad[1] * d2[7] + grad[2] * d2[8] +
278  d2[3] * (grad[3] * d2[3] + grad[4] * d2[4] + grad[5] * d2[5]) +
279  d2[4] * (grad[4] * d2[3] + grad[6] * d2[4] + grad[7] * d2[5]) +
280  d2[5] * (grad[5] * d2[3] + grad[7] * d2[4] + grad[8] * d2[5]);
281 }
Array< OneD, NekDouble > m_grad
Definition: NodeOpti.h:109
double NekDouble

Member Data Documentation

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

Definition at line 74 of file NodeOptiCAD.h.

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

Definition at line 63 of file NodeOptiCAD.h.