Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ProcessVarOpti.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: ProcessJac.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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Calculate Jacobians of elements.
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 
38 
39 #include "ElUtil.h"
40 #include "NodeOpti.h"
41 #include "ProcessVarOpti.h"
43 
44 #include <boost/thread/mutex.hpp>
45 
46 #include <StdRegions/StdPrismExp.h>
47 #include <StdRegions/StdQuadExp.h>
48 #include <StdRegions/StdTetExp.h>
49 #include <StdRegions/StdTriExp.h>
50 
53 
54 using namespace std;
55 using namespace Nektar::NekMeshUtils;
56 
57 namespace Nektar
58 {
59 namespace Utilities
60 {
61 
62 ModuleKey ProcessVarOpti::className =
64  ModuleKey(eProcessModule, "varopti"), ProcessVarOpti::create,
65  "Optimise mesh locations.");
66 
67 ProcessVarOpti::ProcessVarOpti(MeshSharedPtr m) : ProcessModule(m)
68 {
69  // clang-format off
70  m_config["linearelastic"] =
71  ConfigOption(true, "", "Optimise for linear elasticity");
72  m_config["winslow"] =
73  ConfigOption(true, "", "Optimise for winslow");
74  m_config["roca"] =
75  ConfigOption(true, "", "Optimise for roca method");
76  m_config["hyperelastic"] =
77  ConfigOption(true, "", "Optimise for hyper elasticity");
78  m_config["numthreads"] =
79  ConfigOption(false, "1", "Number of threads");
80  m_config["restol"] =
81  ConfigOption(false, "1e-6", "Tolerance criterion");
82  m_config["maxiter"] =
83  ConfigOption(false, "500", "Maximum number of iterations");
84  m_config["nq"] =
85  ConfigOption(false, "-1", "Order of mesh");
86  m_config["region"] =
87  ConfigOption(false, "0.0", "create regions based on target");
88  m_config["resfile"] =
89  ConfigOption(false, "", "writes residual values to file");
90  m_config["histfile"] =
91  ConfigOption(false, "", "histogram of scaled jac");
92  m_config["overint"] =
93  ConfigOption(false, "6", "over integration order");
94  m_config["analytics"] =
95  ConfigOption(false, "", "basic analytics module");
96  // clang-format on
97 }
98 
100 {
101 }
102 
104 {
105  if (m_mesh->m_verbose)
106  {
107  cout << "ProcessVarOpti: Optimising... " << endl;
108  }
109 
110  if (m_config["linearelastic"].beenSet)
111  {
112  m_opti = eLinEl;
113  }
114  else if (m_config["winslow"].beenSet)
115  {
116  m_opti = eWins;
117  }
118  else if (m_config["roca"].beenSet)
119  {
120  m_opti = eRoca;
121  }
122  else if (m_config["hyperelastic"].beenSet)
123  {
124  m_opti = eHypEl;
125  }
126  else
127  {
128  ASSERTL0(false, "not opti type set");
129  }
130 
131  const int maxIter = m_config["maxiter"].as<int>();
132  const NekDouble restol = m_config["restol"].as<NekDouble>();
133 
134  // m_mesh->m_nummode = m_config["nq"].as<int>();
135 
136  EdgeSet::iterator eit;
137  bool fd = false;
138 
139  if (m_config["nq"].beenSet)
140  {
141  m_mesh->m_nummode = m_config["nq"].as<int>();
142  fd = true;
143  }
144 
145  if (!fd)
146  {
147  for (eit = m_mesh->m_edgeSet.begin(); eit != m_mesh->m_edgeSet.end();
148  eit++)
149  {
150  if ((*eit)->m_edgeNodes.size() > 0)
151  {
152  m_mesh->m_nummode = (*eit)->m_edgeNodes.size() + 2;
153  fd = true;
154  break;
155  }
156  }
157  }
158  ASSERTL0(fd, "failed to find order of mesh");
159 
160  int intOrder = m_config["overint"].as<NekDouble>();
161 
162  if (m_mesh->m_verbose)
163  {
164  cout << "Identified mesh order as: " << m_mesh->m_nummode - 1 << endl;
165  }
166 
167  if (m_mesh->m_expDim == 2 && m_mesh->m_spaceDim == 3)
168  {
169  ASSERTL0(false, "cannot deal with manifolds");
170  }
171 
172  m_res = boost::shared_ptr<Residual>(new Residual);
173  m_res->val = 1.0;
174 
175 
176  m_mesh->MakeOrder(m_mesh->m_nummode - 1,
178 
179  if (m_config["analytics"].beenSet)
180  {
181  Analytics();
182  return;
183  }
184 
185  map<LibUtilities::ShapeType, DerivUtilSharedPtr> derivUtils =
186  BuildDerivUtil(intOrder);
187 
188  GetElementMap(intOrder, derivUtils);
189 
190  m_res->startInv = 0;
191  m_res->worstJac = numeric_limits<double>::max();
192  for (int i = 0; i < m_dataSet.size(); i++)
193  {
194  m_dataSet[i]->Evaluate();
195  m_dataSet[i]->InitialMinJac();
196  }
197 
198  vector<ElUtilSharedPtr> elLock;
199 
200  if (m_config["region"].beenSet)
201  {
202  elLock = GetLockedElements(m_config["region"].as<NekDouble>());
203  }
204 
205 
206  vector<vector<NodeSharedPtr> > freenodes = GetColouredNodes(elLock);
207  vector<vector<NodeOptiSharedPtr> > optiNodes;
208 
209  // turn the free nodes into optimisable objects with all required data
210  set<int> check;
211  for (int i = 0; i < freenodes.size(); i++)
212  {
213  vector<NodeOptiSharedPtr> ns;
214  for (int j = 0; j < freenodes[i].size(); j++)
215  {
216  NodeElMap::iterator it = m_nodeElMap.find(freenodes[i][j]->m_id);
217  ASSERTL0(it != m_nodeElMap.end(), "could not find");
218 
219  int optiKind = m_mesh->m_spaceDim;
220 
221  if (freenodes[i][j]->GetNumCadCurve() == 1)
222  {
223  optiKind += 10;
224  }
225  else if (freenodes[i][j]->GetNumCADSurf() == 1)
226  {
227  optiKind += 20;
228  }
229  else
230  {
231  optiKind += 10 * m_mesh->m_expDim;
232  }
233 
234  set<int>::iterator c = check.find(freenodes[i][j]->m_id);
235  ASSERTL0(c == check.end(), "duplicate node");
236  check.insert(freenodes[i][j]->m_id);
237 
238  ns.push_back(GetNodeOptiFactory().CreateInstance(
239  optiKind, freenodes[i][j], it->second, m_res, derivUtils,
240  m_opti));
241  }
242  optiNodes.push_back(ns);
243  }
244 
245  int nset = optiNodes.size();
246  int p = 0;
247  int mn = numeric_limits<int>::max();
248  int mx = 0;
249  for (int i = 0; i < nset; i++)
250  {
251  p += optiNodes[i].size();
252  mn = min(mn, int(optiNodes[i].size()));
253  mx = max(mx, int(optiNodes[i].size()));
254  }
255 
256  if (m_config["histfile"].beenSet)
257  {
258  ofstream histFile;
259  string name = m_config["histfile"].as<string>() + "_start.txt";
260  histFile.open(name.c_str());
261 
262  for (int i = 0; i < m_dataSet.size(); i++)
263  {
264  histFile << m_dataSet[i]->GetScaledJac() << endl;
265  }
266  histFile.close();
267  }
268 
269  if(m_mesh->m_verbose)
270  {
271  cout << scientific << endl;
272  cout << "N elements:\t\t" << m_mesh->m_element[m_mesh->m_expDim].size() - elLock.size() << endl
273  << "N elements invalid:\t" << m_res->startInv << endl
274  << "Worst jacobian:\t\t" << m_res->worstJac << endl
275  << "N free nodes:\t\t" << m_res->n << endl
276  << "N Dof:\t\t\t" << m_res->nDoF << endl
277  << "N color sets:\t\t" << nset << endl
278  << "Avg set colors:\t\t" << p/nset << endl
279  << "Min set:\t\t" << mn << endl
280  << "Max set:\t\t" << mx << endl
281  << "Residual tolerance:\t" << restol << endl;
282  }
283 
284  int nThreads = m_config["numthreads"].as<int>();
285 
286  int ctr = 0;
288  tms.SetThreadingType("ThreadManagerBoost");
291 
292  Timer t;
293  t.Start();
294 
295  ofstream resFile;
296  if (m_config["resfile"].beenSet)
297  {
298  resFile.open(m_config["resfile"].as<string>().c_str());
299  }
300 
301  for (int i = 0; i < optiNodes.size(); i++)
302  {
303  vector<Thread::ThreadJob *> jobs(optiNodes[i].size());
304  for (int j = 0; j < optiNodes[i].size(); j++)
305  {
306  optiNodes[i][j]->CalcMinJac();
307  }
308  }
309 
310  while (m_res->val > restol)
311  {
312  ctr++;
313  m_res->val = 0.0;
314  m_res->func = 0.0;
315  m_res->nReset[0] = 0;
316  m_res->nReset[1] = 0;
317  m_res->nReset[2] = 0;
318  m_res->alphaI = 0;
319  for (int i = 0; i < optiNodes.size(); i++)
320  {
321  vector<Thread::ThreadJob *> jobs(optiNodes[i].size());
322  for (int j = 0; j < optiNodes[i].size(); j++)
323  {
324  jobs[j] = optiNodes[i][j]->GetJob();
325  }
326 
327  tm->SetNumWorkers(0);
328  tm->QueueJobs(jobs);
329  tm->SetNumWorkers(nThreads);
330  tm->Wait();
331  }
332 
333  m_res->startInv = 0;
334  m_res->worstJac = numeric_limits<double>::max();
335 
336  vector<Thread::ThreadJob *> elJobs(m_dataSet.size());
337  for (int i = 0; i < m_dataSet.size(); i++)
338  {
339  elJobs[i] = m_dataSet[i]->GetJob();
340  }
341 
342  tm->SetNumWorkers(0);
343  tm->QueueJobs(elJobs);
344  tm->SetNumWorkers(nThreads);
345  tm->Wait();
346 
347  if (m_config["resfile"].beenSet)
348  {
349  resFile << m_res->val << " " << m_res->worstJac << " "
350  << m_res->func << endl;
351  }
352 
353  if(m_mesh->m_verbose)
354  {
355  cout << ctr
356  << "\tResidual: " << m_res->val
357  << "\tMin Jac: " << m_res->worstJac
358  << "\tInvalid: " << m_res->startInv
359  << "\tReset nodes: " << m_res->nReset[0] << "/" << m_res->nReset[1]
360  << "/" << m_res->nReset[2] << "\tFunctional: " << m_res->func
361  << endl;
362  }
363 
364  if(ctr >= maxIter)
365  {
366  break;
367  }
368  }
369 
370  if (m_config["histfile"].beenSet)
371  {
372  ofstream histFile;
373  string name = m_config["histfile"].as<string>() + "_end.txt";
374  histFile.open(name.c_str());
375 
376  for (int i = 0; i < m_dataSet.size(); i++)
377  {
378  histFile << m_dataSet[i]->GetScaledJac() << endl;
379  }
380  histFile.close();
381  }
382  if (m_config["resfile"].beenSet)
383  {
384  resFile.close();
385  }
386 
387  t.Stop();
388 
390 
391  if(m_mesh->m_verbose)
392  {
393  cout << "Time to compute: " << t.TimePerTest(1) << endl;
394  cout << "Invalid at end: " << m_res->startInv << endl;
395  cout << "Worst at end: " << m_res->worstJac << endl;
396  }
397 }
398 
400 {
401 public:
404  : NodalUtilTriangle(degree, r, s)
405  {
406  }
407 
409  {
410  }
411 
412 protected:
413  virtual NekVector<NekDouble> v_OrthoBasis(const int mode)
414  {
415  // Monomial basis.
416  std::pair<int, int> modes = m_ordering[mode];
418 
419  for (int i = 0; i < m_numPoints; ++i)
420  {
421  ret(i) =
422  pow(m_xi[0][i], modes.first) * pow(m_xi[1][i], modes.second);
423  }
424 
425  return ret;
426  }
427 
429  const int mode)
430  {
431  ASSERTL0(false, "not supported");
432  return NekVector<NekDouble>();
433  }
434 
435  virtual boost::shared_ptr<NodalUtil> v_CreateUtil(
437  {
439  m_degree, xi[0], xi[1]);
440  }
441 };
442 
444 {
445  // Grab the first element from the list
446  ElementSharedPtr elmt = m_mesh->m_element[m_mesh->m_expDim][0];
447 
448  // Get curved nodes
449  vector<NodeSharedPtr> nodes;
450  elmt->GetCurvedNodes(nodes);
451 
452  // We're going to investigate only the first node (corner node)
453  NodeSharedPtr node = nodes[4];
454 
455  // Loop over overintegration orders
456  const int nPoints = 200;
457  const int overInt = 40;
458  const NekDouble originX = -1.0;
459  const NekDouble originY = -1.0;
460  const NekDouble length = 2.0;
461  const NekDouble dx = length / (nPoints - 1);
462 
463  cout << "# overint = " << overInt << endl;
464  cout << "# Columns: x, y, over-integration orders (0 -> " << overInt - 1
465  << "), "
466  << " min(scaledJac)" << endl;
467 
468  // Loop over square defined by (originX, originY), length
469  for (int k = 0; k < nPoints; ++k)
470  {
471  node->m_y = originY + k * dx;
472  for (int j = 0; j < nPoints; ++j)
473  {
474  node->m_x = originX + j * dx;
475  cout << node->m_x << " " << node->m_y << " ";
476 
477  NekDouble minJacNew;
478 
479  for (int i = 0; i < overInt; ++i)
480  {
481  // Clear any existing node to element mapping.
482  m_dataSet.clear();
483  m_nodeElMap.clear();
484 
485  // Build deriv utils and element map.
486  map<LibUtilities::ShapeType, DerivUtilSharedPtr> derivUtils =
487  BuildDerivUtil(i);
488 
489  // Reconstruct element map
490  GetElementMap(i, derivUtils);
491 
492  for (int j = 0; j < m_dataSet.size(); j++)
493  {
494  m_dataSet[j]->Evaluate();
495  m_dataSet[j]->InitialMinJac();
496  }
497 
498  // Create NodeOpti object.
499  NodeOptiSharedPtr nodeOpti =
501  m_mesh->m_spaceDim * 11, node,
502  m_nodeElMap.find(node->m_id)->second, m_res, derivUtils,
503  m_opti);
504 
505  minJacNew = 0.0;
506 
507  // Evaluate functional.
508  nodeOpti->CalcMinJac();
509  cout << nodeOpti->GetFunctional<2>(minJacNew) << " ";
510  // NekDouble eigen;
511  // nodeOpti->GetFunctional<2>(minJacNew);
512  // nodeOpti->MinEigen<2>(eigen);
513  // cout << eigen << " ";
514  }
515 
516  cout << minJacNew << endl;
517  }
518  }
519 }
520 }
521 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
NodeOptiFactory & GetNodeOptiFactory()
Definition: NodeOpti.cpp:52
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
Specialisation of the NodalUtil class to support nodal triangular elements.
Definition: NodalUtil.h:170
virtual NekVector< NekDouble > v_OrthoBasis(const int mode)
Return the value of the modal functions for the triangular element at the nodal points m_xi for a giv...
std::vector< ElUtilSharedPtr > GetLockedElements(NekDouble thres)
void GetElementMap(int o, std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > derMap)
std::map< LibUtilities::ShapeType, DerivUtilSharedPtr > BuildDerivUtil(int o)
std::vector< ElUtilSharedPtr > m_dataSet
STL namespace.
pair< ModuleType, string > ModuleKey
virtual boost::shared_ptr< NodalUtil > v_CreateUtil(Array< OneD, Array< OneD, NekDouble > > &xi)
Construct a NodalUtil object of the appropriate element type for a given set of points.
std::vector< std::vector< NodeSharedPtr > > GetColouredNodes(std::vector< ElUtilSharedPtr > elLock)
NodalUtilTriMonomial(int degree, Array< OneD, NekDouble > r, Array< OneD, NekDouble > s)
void Stop()
Definition: Timer.cpp:62
std::vector< std::pair< int, int > > m_ordering
Mapping from the indexing of the basis to a continuous ordering.
Definition: NodalUtil.h:184
int m_degree
Degree of the nodal element.
Definition: NodalUtil.h:108
Represents a command-line configuration option.
boost::shared_ptr< Node > NodeSharedPtr
Definition: Node.h:50
double NekDouble
std::map< std::string, ConfigOption > m_config
List of configuration values.
boost::shared_ptr< NodeOpti > NodeOptiSharedPtr
Definition: NodeOpti.h:107
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
ThreadManagerSharedPtr CreateInstance(const ThreadManagerName t, unsigned int nThr)
Creates an instance of a ThreadManager (which one is determined by a previous call to SetThreadingTyp...
Definition: Thread.cpp:200
boost::shared_ptr< Mesh > MeshSharedPtr
Shared pointer to a mesh.
Definition: Mesh.h:147
Abstract base class for processing modules.
virtual NekVector< NekDouble > v_OrthoBasisDeriv(const int dir, const int mode)
Return the value of the derivative of the modal functions for the triangular element at the nodal poi...
Array< OneD, Array< OneD, NekDouble > > m_xi
Coordinates of the nodal points defining the basis.
Definition: NodalUtil.h:112
void SetThreadingType(const std::string &p_type)
Sets what ThreadManagers will be created in CreateInstance.
Definition: Thread.cpp:162
int m_numPoints
Total number of nodal points.
Definition: NodalUtil.h:110
NodalUtilTriangle(int degree, Array< OneD, NekDouble > r, Array< OneD, NekDouble > s)
Construct the nodal utility class for a triangle.
Definition: NodalUtil.cpp:239
boost::shared_ptr< Element > ElementSharedPtr
Definition: Edge.h:49
void Start()
Definition: Timer.cpp:51
std::pair< ModuleType, std::string > ModuleKey
NekDouble TimePerTest(unsigned int n)
Returns amount of seconds per iteration in a test with n iterations.
Definition: Timer.cpp:108
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:52
ModuleFactory & GetModuleFactory()
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215
boost::shared_ptr< ThreadManager > ThreadManagerSharedPtr
Definition: Thread.h:74