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

#include <NodeOptiCAD.h>

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

Public Member Functions

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

CADSurfSharedPtr surf
 

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

Constructor & Destructor Documentation

◆ NodeOpti2D3D()

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

Definition at line 79 of file NodeOptiCAD.h.

83  : NodeOpti(n, e, r, d, o, 3), surf(s)
84  {
85  }
NodeOpti(NodeSharedPtr n, std::vector< ElUtilSharedPtr > e, ResidualSharedPtr r, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > d, optiType o, int dim)
Definition: NodeOpti.h:60

◆ ~NodeOpti2D3D()

Nektar::Utilities::NodeOpti2D3D::~NodeOpti2D3D ( )
inline

Definition at line 87 of file NodeOptiCAD.h.

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

87 {};

Member Function Documentation

◆ create()

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

Definition at line 92 of file NodeOptiCAD.h.

References Nektar::Utilities::NodeOpti1D3D::ProcessGradient().

95  {
96  std::vector<CADSurfSharedPtr> ss = n->GetCADSurfs();
97  return NodeOptiSharedPtr(new NodeOpti2D3D(n, e, r, d, o, ss[0]));
98  }
std::shared_ptr< NodeOpti > NodeOptiSharedPtr
Definition: NodeOpti.h:135
NodeOpti2D3D(NodeSharedPtr n, std::vector< ElUtilSharedPtr > e, ResidualSharedPtr r, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > d, optiType o, CADSurfSharedPtr s)
Definition: NodeOptiCAD.h:79

◆ Optimise()

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

Implements Nektar::Utilities::NodeOpti.

Definition at line 145 of file NodeOptiCAD.cpp.

References ASSERTL0.

146 {
147  NekDouble minJacNew;
148  NekDouble currentW = GetFunctional<3>(minJacNew);
149  NekDouble newVal = currentW;
150 
151  if (m_grad[0] * m_grad[0] + m_grad[1] * m_grad[1] + m_grad[2] * m_grad[2] >
152  gradTol())
153  {
154  // modify the gradient to be on the cad system
155  ProcessGradient();
156 
157  // needs to optimise
158  Array<OneD, NekDouble> uvc = m_node->GetCADSurfInfo(surf->GetId());
159  NekDouble xc = m_node->m_x;
160  NekDouble yc = m_node->m_y;
161  NekDouble zc = m_node->m_z;
162  Array<OneD, NekDouble> uvt(2);
163  vector<NekDouble> bd(4);
164  surf->GetBounds(bd[0], bd[1], bd[2], bd[3]);
165 
166  vector<NekDouble> sk(2);
167  NekDouble val;
168 
169  // Calculate minimum eigenvalue
170  MinEigen<2>(val);
171 
172  if (val < 1e-6)
173  {
174  // Add constant identity to Hessian matrix.
175  m_grad[2] += 1e-6 - val;
176  m_grad[4] += 1e-6 - val;
177  }
178 
179  sk[0] = -1.0 / (m_grad[2] * m_grad[4] - m_grad[3] * m_grad[3]) *
180  (m_grad[4] * m_grad[0] - m_grad[3] * m_grad[1]);
181  sk[1] = -1.0 / (m_grad[2] * m_grad[4] - m_grad[3] * m_grad[3]) *
182  (m_grad[2] * m_grad[1] - m_grad[3] * m_grad[0]);
183 
184  bool found = false;
185  NekDouble pg = (m_grad[0] * sk[0] + m_grad[1] * sk[1]);
186 
187  // normal gradient line Search
188  NekDouble alpha = 1.0;
189  while (alpha > alphaTol())
190  {
191  uvt[0] = uvc[0] + alpha * sk[0];
192  uvt[1] = uvc[1] + alpha * sk[1];
193 
194  if (uvt[0] < bd[0] || uvt[0] > bd[1] || uvt[1] < bd[2] ||
195  uvt[1] > bd[3])
196  {
197  alpha /= 2.0;
198  continue;
199  }
200 
201  surf->P(uvt, m_node->m_x, m_node->m_y, m_node->m_z);
202 
203  newVal = GetFunctional<3>(minJacNew, false);
204 
205  // Wolfe conditions
206  if (newVal <= currentW + c1() * (alpha * pg))
207  {
208  found = true;
209  break;
210  }
211 
212  alpha /= 2.0;
213  }
214 
215  if (!found)
216  {
217  // reset the node
218  surf->P(uvc, m_node->m_x, m_node->m_y, m_node->m_z);
219 
220  mtx.lock();
221  m_res->nReset[1]++;
222  mtx.unlock();
223  }
224  else
225  {
226  m_minJac = minJacNew;
227  m_node->Move(m_node->m_x, m_node->m_y, m_node->m_z, surf->GetId(),
228  uvt);
229  }
230 
231  mtx.lock();
232  m_res->val = max(sqrt((m_node->m_x - xc) * (m_node->m_x - xc) +
233  (m_node->m_y - yc) * (m_node->m_y - yc) +
234  (m_node->m_z - zc) * (m_node->m_z - zc)),
235  m_res->val);
236  mtx.unlock();
237 
238  Array<OneD, NekDouble> uva = m_node->GetCADSurfInfo(surf->GetId());
239  if (uva[0] < bd[0] || uva[0] > bd[1] || uva[1] < bd[2] ||
240  uva[1] > bd[3])
241  {
242  ASSERTL0(false, "something when very wrong and the node finished "
243  "out of its parameter plane");
244  }
245  }
246 
247  mtx.lock();
248  m_res->func += newVal;
249  mtx.unlock();
250 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
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::NodeOpti2D3D::ProcessGradient ( )
private

Definition at line 272 of file NodeOptiCAD.cpp.

273 {
274  Array<OneD, NekDouble> uvc = m_node->GetCADSurfInfo(surf->GetId());
275 
276  vector<NekDouble> grad = m_grad;
277  m_grad = vector<NekDouble>(5, 0.0);
278 
279  Array<OneD, NekDouble> d2 = surf->D2(uvc);
280  // r[0] x
281  // r[1] y
282  // r[2] z
283  // r[3] dx/du a
284  // r[4] dy/du b
285  // r[5] dz/du c
286  // r[6] dx/dv d
287  // r[7] dy/dv e
288  // r[8] dz/dv f
289  // r[9] d2x/du2
290  // r[10] d2y/du2
291  // r[11] d2z/du2
292  // r[12] d2x/dv2
293  // r[13] d2y/dv2
294  // r[14] d2z/dv2
295  // r[15] d2x/dudv
296  // r[16] d2y/dudv
297  // r[17] d2z/dudv
298  //
299  // grad[0] d/dx
300  // grad[1] d/dy
301  // grad[2] d/dz
302  // grad[3] d2/dx2
303  // grad[4] d2/dxdy
304  // grad[5] d2/dxdz
305  // grad[6] d2/dy2
306  // grad[7] d2/dydz
307  // grad[8] d2/dz2
308 
309  // Gradients
310  m_grad[0] = d2[3] * grad[0] + d2[4] * grad[1] + d2[5] * grad[2];
311  m_grad[1] = d2[6] * grad[0] + d2[7] * grad[1] + d2[8] * grad[2];
312 
313  m_grad[2] = grad[0] * d2[9] + grad[1] * d2[10] + grad[2] * d2[11] +
314  grad[3] * d2[3] * d2[3] + grad[6] * d2[4] * d2[4] +
315  grad[8] * d2[5] * d2[5] + 2.0 * grad[4] * d2[3] * d2[4] +
316  2.0 * grad[5] * d2[3] * d2[5] + 2.0 * grad[7] * d2[4] * d2[5];
317 
318  m_grad[4] = grad[0] * d2[12] + grad[1] * d2[13] + grad[2] * d2[14] +
319  grad[3] * d2[6] * d2[6] + grad[6] * d2[7] * d2[7] +
320  grad[8] * d2[8] * d2[8] + 2.0 * grad[4] * d2[6] * d2[7] +
321  2.0 * grad[5] * d2[6] * d2[8] + 2.0 * grad[7] * d2[7] * d2[8];
322 
323  m_grad[3] = grad[0] * d2[15] + grad[1] * d2[16] + grad[2] * d2[17] +
324  grad[3] * d2[3] * d2[6] + grad[6] * d2[4] * d2[7] +
325  grad[8] * d2[5] * d2[8] +
326  grad[4] * (d2[3] * d2[7] + d2[4] * d2[6]) +
327  grad[5] * (d2[3] * d2[8] + d2[5] * d2[6]) +
328  grad[7] * (d2[4] * d2[8] + d2[5] * d2[7]);
329 }
std::vector< NekDouble > m_grad
Definition: NodeOpti.h:109

Member Data Documentation

◆ m_type

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

Definition at line 91 of file NodeOptiCAD.h.

◆ surf

CADSurfSharedPtr Nektar::Utilities::NodeOpti2D3D::surf
private

Definition at line 102 of file NodeOptiCAD.h.