Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Public Member Functions | Static Public Member Functions | Static Public Attributes | List of all members
Nektar::Utilities::NodeOpti2D2D Class Reference

#include <NodeOpti.h>

Inheritance diagram for Nektar::Utilities::NodeOpti2D2D:
Inheritance graph
[legend]
Collaboration diagram for Nektar::Utilities::NodeOpti2D2D:
Collaboration graph
[legend]

Public Member Functions

 NodeOpti2D2D (NodeSharedPtr n, std::vector< ElUtilSharedPtr > e, ResidualSharedPtr r, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > d, optiType o)
 
 ~NodeOpti2D2D ()
 
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
 

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 170 of file NodeOpti.h.

Constructor & Destructor Documentation

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

Definition at line 173 of file NodeOpti.h.

Referenced by create().

177  : NodeOpti(n, e, r, d, o, 2)
178  {
179  }
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::NodeOpti2D2D::~NodeOpti2D2D ( )
inline

Definition at line 181 of file NodeOpti.h.

181 {};

Member Function Documentation

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

Definition at line 186 of file NodeOpti.h.

References NodeOpti2D2D().

189  {
190  return NodeOptiSharedPtr(new NodeOpti2D2D(n, e, r, d, o));
191  }
NodeOpti2D2D(NodeSharedPtr n, std::vector< ElUtilSharedPtr > e, ResidualSharedPtr r, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > d, optiType o)
Definition: NodeOpti.h:173
boost::shared_ptr< NodeOpti > NodeOptiSharedPtr
Definition: NodeOpti.h:135
void Nektar::Utilities::NodeOpti2D2D::Optimise ( )
virtual

Implements Nektar::Utilities::NodeOpti.

Definition at line 74 of file NodeOpti.cpp.

References Nektar::Utilities::mtx.

75 {
76  NekDouble minJacNew;
77  NekDouble currentW = GetFunctional<2>(minJacNew);
78  NekDouble newVal = currentW;
79 
80  // Gradient already zero
81  if (m_grad[0] * m_grad[0] + m_grad[1] * m_grad[1] > gradTol())
82  {
83  // needs to optimise
84  NekDouble xc = m_node->m_x;
85  NekDouble yc = m_node->m_y;
86 
87  Array<OneD, NekDouble> sk(2);
88  NekDouble val;
89 
90  // Calculate minimum eigenvalue
91  MinEigen<2>(val);
92 
93  if (val < 1e-6)
94  {
95  // Add constant identity to Hessian matrix.
96  m_grad[2] += 1e-6 - val;
97  m_grad[4] += 1e-6 - val;
98  }
99 
100  sk[0] = -1.0 / (m_grad[2] * m_grad[4] - m_grad[3] * m_grad[3]) *
101  (m_grad[4] * m_grad[0] - m_grad[3] * m_grad[1]);
102  sk[1] = -1.0 / (m_grad[2] * m_grad[4] - m_grad[3] * m_grad[3]) *
103  (m_grad[2] * m_grad[1] - m_grad[3] * m_grad[0]);
104 
105  bool found = false;
106  NekDouble pg = (m_grad[0] * sk[0] + m_grad[1] * sk[1]);
107  // normal gradient line Search
108  NekDouble alpha = 1.0;
109 
110  while (alpha > alphaTol())
111  {
112  // Update node
113  m_node->m_x = xc + alpha * sk[0];
114  m_node->m_y = yc + alpha * sk[1];
115 
116  newVal = GetFunctional<2>(minJacNew, false);
117 
118  // Wolfe conditions
119  if (newVal <= currentW + c1() * (alpha * pg))
120  {
121  found = true;
122  break;
123  }
124 
125  alpha /= 2.0;
126  }
127 
128  if (!found)
129  {
130  // reset the node
131  m_node->m_x = xc;
132  m_node->m_y = yc;
133 
134  mtx.lock();
135  m_res->nReset[2]++;
136  mtx.unlock();
137  }
138  else
139  {
140  m_minJac = minJacNew;
141 
142  mtx.lock();
143  if (alpha < 1.0)
144  {
145  m_res->alphaI++;
146  }
147  mtx.unlock();
148  }
149 
150  mtx.lock();
151  m_res->val = max(sqrt((m_node->m_x - xc) * (m_node->m_x - xc) +
152  (m_node->m_y - yc) * (m_node->m_y - yc)),
153  m_res->val);
154  mtx.unlock();
155  }
156 
157  mtx.lock();
158  m_res->func += newVal;
159  mtx.unlock();
160 }
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

Member Data Documentation

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

Definition at line 185 of file NodeOpti.h.