Nektar++
NodeOpti.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: NodeOpti.cpp
4 //
5 // For more information, please see: http://www.nektar.info/
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Calculate jacobians of elements.
32 //
33 ////////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/multi_array.hpp>
36 
37 #include <limits>
38 
39 #include "NodeOpti.h"
40 #include "Evaluator.hxx"
41 #include "Hessian.hxx"
42 
43 using namespace std;
44 using namespace Nektar::NekMeshUtils;
45 
46 namespace Nektar
47 {
48 namespace Utilities
49 {
50 
52 {
53  static NodeOptiFactory asd;
54  return asd;
55 }
56 
57 void NodeOpti::CalcMinJac()
58 {
59  m_minJac = numeric_limits<double>::max();
60  for (auto &typeIt : m_data)
61  {
62  for (int i = 0; i < typeIt.second.size(); i++)
63  {
64  m_minJac = min(m_minJac, typeIt.second[i]->GetMinJac());
65  }
66  }
67 }
68 
69 int NodeOpti2D2D::m_type = GetNodeOptiFactory().RegisterCreatorFunction(
70  22, NodeOpti2D2D::create, "2D2D");
71 
72 void NodeOpti2D2D::Optimise()
73 {
74  NekDouble minJacNew;
75  NekDouble currentW = GetFunctional<2>(minJacNew);
76  NekDouble newVal = currentW;
77 
78  // Gradient already zero
79  if (m_grad[0] * m_grad[0] + m_grad[1] * m_grad[1] > gradTol())
80  {
81  // needs to optimise
82  NekDouble xc = m_node->m_x;
83  NekDouble yc = m_node->m_y;
84 
85  vector<NekDouble> sk(2);
86  NekDouble val;
87 
88  // Calculate minimum eigenvalue
89  MinEigen<2>(val);
90 
91  if (val < 1e-6)
92  {
93  // Add constant identity to Hessian matrix.
94  m_grad[2] += 1e-6 - val;
95  m_grad[4] += 1e-6 - val;
96  }
97 
98  sk[0] = -1.0 / (m_grad[2] * m_grad[4] - m_grad[3] * m_grad[3]) *
99  (m_grad[4] * m_grad[0] - m_grad[3] * m_grad[1]);
100  sk[1] = -1.0 / (m_grad[2] * m_grad[4] - m_grad[3] * m_grad[3]) *
101  (m_grad[2] * m_grad[1] - m_grad[3] * m_grad[0]);
102 
103  bool found = false;
104  NekDouble pg = (m_grad[0] * sk[0] + m_grad[1] * sk[1]);
105  // normal gradient line Search
106  NekDouble alpha = 1.0;
107 
108  while (alpha > alphaTol())
109  {
110  // Update node
111  m_node->m_x = xc + alpha * sk[0];
112  m_node->m_y = yc + alpha * sk[1];
113 
114  newVal = GetFunctional<2>(minJacNew, false);
115 
116  // Wolfe conditions
117  if (newVal <= currentW + c1() * (alpha * pg))
118  {
119  found = true;
120  break;
121  }
122 
123  alpha /= 2.0;
124  }
125 
126  if (!found)
127  {
128  // reset the node
129  m_node->m_x = xc;
130  m_node->m_y = yc;
131 
132  mtx.lock();
133  m_res->nReset[2]++;
134  mtx.unlock();
135  }
136  else
137  {
138  m_minJac = minJacNew;
139 
140  mtx.lock();
141  if (alpha < 1.0)
142  {
143  m_res->alphaI++;
144  }
145  mtx.unlock();
146  }
147 
148  mtx.lock();
149  m_res->val = max(sqrt((m_node->m_x - xc) * (m_node->m_x - xc) +
150  (m_node->m_y - yc) * (m_node->m_y - yc)),
151  m_res->val);
152  mtx.unlock();
153  }
154 
155  mtx.lock();
156  m_res->func += newVal;
157  mtx.unlock();
158 }
159 
160 int NodeOpti3D3D::m_type = GetNodeOptiFactory().RegisterCreatorFunction(
161  33, NodeOpti3D3D::create, "3D3D");
162 
163 void NodeOpti3D3D::Optimise()
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 }
273 
274 NodeOptiJob *NodeOpti::GetJob()
275 {
276  return new NodeOptiJob(this);
277 }
278 }
279 }
NodeOptiFactory & GetNodeOptiFactory()
Definition: NodeOpti.cpp:51
STL namespace.
double NekDouble
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:199
Provides a generic Factory class.
Definition: NekFactory.hpp:103