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

#include <NodeOpti.h>

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

Public Member Functions

 NodeOpti3D3D (NodeSharedPtr n, std::vector< ElUtilSharedPtr > e, ResidualSharedPtr r, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > d, optiType o)
 
 ~NodeOpti3D3D ()
 
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
 
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 144 of file NodeOpti.h.

Constructor & Destructor Documentation

◆ NodeOpti3D3D()

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

Definition at line 147 of file NodeOpti.h.

151  : NodeOpti(n, e, r, d, o, 3)
152  {
153  }
NodeOpti(NodeSharedPtr n, std::vector< ElUtilSharedPtr > e, ResidualSharedPtr r, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > d, optiType o, int dim)
Definition: NodeOpti.h:60

◆ ~NodeOpti3D3D()

Nektar::Utilities::NodeOpti3D3D::~NodeOpti3D3D ( )
inline

Definition at line 155 of file NodeOpti.h.

References Nektar::Utilities::NodeOpti::Optimise().

155 {};

Member Function Documentation

◆ create()

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

Definition at line 160 of file NodeOpti.h.

163  {
164  return NodeOptiSharedPtr(new NodeOpti3D3D(n, e, r, d, o));
165  }
std::shared_ptr< NodeOpti > NodeOptiSharedPtr
Definition: NodeOpti.h:135
NodeOpti3D3D(NodeSharedPtr n, std::vector< ElUtilSharedPtr > e, ResidualSharedPtr r, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > d, optiType o)
Definition: NodeOpti.h:147

◆ Optimise()

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

Implements Nektar::Utilities::NodeOpti.

Definition at line 163 of file NodeOpti.cpp.

References Nektar::Utilities::mtx.

164 {
165  NekDouble minJacNew;
166  NekDouble currentW = GetFunctional<3>(minJacNew);
167  NekDouble newVal = currentW;
168 
169  if (m_grad[0] * m_grad[0] + m_grad[1] * m_grad[1] + m_grad[2] * m_grad[2] >
170  gradTol())
171  {
172  // needs to optimise
173  NekDouble xc = m_node->m_x;
174  NekDouble yc = m_node->m_y;
175  NekDouble zc = m_node->m_z;
176 
177  vector<NekDouble> sk(3);
178  NekDouble val;
179 
180  // Calculate minimum eigenvalue
181  MinEigen<3>(val);
182 
183  if (val < 1e-6)
184  {
185  // Add constant identity to Hessian matrix.
186  m_grad[3] += 1e-6 - val;
187  m_grad[6] += 1e-6 - val;
188  m_grad[8] += 1e-6 - val;
189  }
190 
191  // calculate sk
192  NekDouble det =
193  m_grad[3] * (m_grad[6] * m_grad[8] - m_grad[7] * m_grad[7]) -
194  m_grad[4] * (m_grad[4] * m_grad[8] - m_grad[5] * m_grad[7]) +
195  m_grad[5] * (m_grad[4] * m_grad[7] - m_grad[5] * m_grad[6]);
196 
197  sk[0] = m_grad[0] * (m_grad[6] * m_grad[8] - m_grad[7] * m_grad[7]) +
198  m_grad[1] * (m_grad[5] * m_grad[7] - m_grad[4] * m_grad[8]) +
199  m_grad[2] * (m_grad[4] * m_grad[7] - m_grad[3] * m_grad[7]);
200  sk[1] = m_grad[0] * (m_grad[7] * m_grad[5] - m_grad[4] * m_grad[5]) +
201  m_grad[1] * (m_grad[3] * m_grad[8] - m_grad[5] * m_grad[5]) +
202  m_grad[2] * (m_grad[4] * m_grad[5] - m_grad[3] * m_grad[7]);
203  sk[2] = m_grad[0] * (m_grad[4] * m_grad[7] - m_grad[6] * m_grad[5]) +
204  m_grad[1] * (m_grad[4] * m_grad[5] - m_grad[3] * m_grad[7]) +
205  m_grad[2] * (m_grad[3] * m_grad[6] - m_grad[4] * m_grad[4]);
206 
207  sk[0] /= det * -1.0;
208  sk[1] /= det * -1.0;
209  sk[2] /= det * -1.0;
210 
211  bool found = false;
212 
213  NekDouble pg =
214  (m_grad[0] * sk[0] + m_grad[1] * sk[1] + m_grad[2] * sk[2]);
215 
216  // normal gradient line Search
217  NekDouble alpha = 1.0;
218 
219  while (alpha > alphaTol())
220  {
221  // Update node
222  m_node->m_x = xc + alpha * sk[0];
223  m_node->m_y = yc + alpha * sk[1];
224  m_node->m_z = zc + alpha * sk[2];
225 
226  newVal = GetFunctional<3>(minJacNew, false);
227  // dont need the hessian again this function updates G to be the new
228  // location
229  //
230  // Wolfe conditions
231  if (newVal <= currentW + c1() * alpha * pg)
232  {
233  found = true;
234  break;
235  }
236 
237  alpha /= 2.0;
238  }
239 
240  if (!found)
241  {
242  m_node->m_x = xc;
243  m_node->m_y = yc;
244  m_node->m_z = zc;
245 
246  mtx.lock();
247  m_res->nReset[2]++;
248  mtx.unlock();
249  }
250  else
251  {
252  m_minJac = minJacNew;
253 
254  mtx.lock();
255  if (alpha < 1.0)
256  {
257  m_res->alphaI++;
258  }
259  mtx.unlock();
260  }
261 
262  mtx.lock();
263  m_res->val = max(sqrt((m_node->m_x - xc) * (m_node->m_x - xc) +
264  (m_node->m_y - yc) * (m_node->m_y - yc) +
265  (m_node->m_z - zc) * (m_node->m_z - zc)),
266  m_res->val);
267  mtx.unlock();
268  }
269  mtx.lock();
270  m_res->func += newVal;
271  mtx.unlock();
272 }
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

Member Data Documentation

◆ m_type

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

Definition at line 159 of file NodeOpti.h.