Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ProcessQualityMetric.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: ProcessQualityMetric.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: Compute quality metric of Roca et al.
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 
36 #include <iostream>
37 #include <string>
38 using namespace std;
39 
40 #include "ProcessQualityMetric.h"
41 
45 #include <StdRegions/StdPrismExp.h>
46 #include <StdRegions/StdQuadExp.h>
47 #include <StdRegions/StdTetExp.h>
48 #include <StdRegions/StdHexExp.h>
49 #include <StdRegions/StdTriExp.h>
50 
51 namespace Nektar
52 {
53 namespace FieldUtils
54 {
55 
56 ModuleKey ProcessQualityMetric::className =
58  ModuleKey(eProcessModule, "qualitymetric"),
59  ProcessQualityMetric::create,
60  "add quality metric to field.");
61 
62 ProcessQualityMetric::ProcessQualityMetric(FieldSharedPtr f) : ProcessModule(f)
63 {
64  m_config["scaled"] =
65  ConfigOption(true, "", "use scaled jacobian instead");
66 }
67 
69 {
70 }
71 
72 void ProcessQualityMetric::Process(po::variables_map &vm)
73 {
74  if (m_f->m_verbose)
75  {
76  if (m_f->m_comm->TreatAsRankZero())
77  {
78  cout << "ProcessQualityMetric: Adding quality metric to field"
79  << endl;
80  }
81  }
82 
83  Array<OneD, NekDouble> &phys = m_f->m_exp[0]->UpdatePhys();
84  Array<OneD, NekDouble> &coeffs = m_f->m_exp[0]->UpdateCoeffs();
85 
86  for (int i = 0; i < m_f->m_exp[0]->GetExpSize(); ++i)
87  {
88  // copy Jacobian into field
89  LocalRegions::ExpansionSharedPtr Elmt = m_f->m_exp[0]->GetExp(i);
90  int offset = m_f->m_exp[0]->GetPhys_Offset(i);
91  Array<OneD, NekDouble> q = GetQ(Elmt,m_config["scaled"].m_beenSet);
92  Array<OneD, NekDouble> out = phys + offset;
93 
94  ASSERTL0(q.num_elements() == Elmt->GetTotPoints(),
95  "number of points mismatch");
96  Vmath::Vcopy(q.num_elements(), q, 1, out, 1);
97  }
98 
99  m_f->m_exp[0]->FwdTrans_IterPerExp(phys, coeffs);
100 
101  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef =
102  m_f->m_exp[0]->GetFieldDefinitions();
103  std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
104 
105  for (int i = 0; i < FieldDef.size(); ++i)
106  {
107  FieldDef[i]->m_fields.push_back("QualityMetric");
108  m_f->m_exp[0]->AppendFieldData(FieldDef[i], FieldData[i]);
109  }
110 
111  m_f->m_fielddef = FieldDef;
112  m_f->m_data = FieldData;
113 }
114 
117 {
118  vector<DNekMat> ret;
119 
120  if (geom->GetShapeType() == LibUtilities::eQuadrilateral)
121  {
122  vector<Array<OneD, NekDouble> > xy;
123  for (int i = 0; i < geom->GetNumVerts(); i++)
124  {
125  Array<OneD, NekDouble> loc(2);
126  SpatialDomains::PointGeomSharedPtr p = geom->GetVertex(i);
127  p->GetCoords(loc);
128  xy.push_back(loc);
129  }
130 
132  Array<OneD, NekDouble> u = b[0]->GetZ();
133  Array<OneD, NekDouble> v = b[1]->GetZ();
134 
135  for (int j = 0; j < b[1]->GetNumPoints(); j++)
136  {
137  for (int i = 0; i < b[0]->GetNumPoints(); i++)
138  {
139  NekDouble a1 = 0.5 * (1.0 - u[i]), a2 = 0.5 * (1.0 + u[i]);
140  NekDouble b1 = 0.5 * (1.0 - v[j]), b2 = 0.5 * (1.0 + v[j]);
141  DNekMat dxdz(2, 2, 1.0, eFULL);
142 
143  dxdz(0, 0) = 0.5 * (-b1 * xy[0][0] + b1 * xy[1][0] +
144  b2 * xy[2][0] - b2 * xy[3][0]);
145  dxdz(1, 0) = 0.5 * (-b1 * xy[0][1] + b1 * xy[1][1] +
146  b2 * xy[2][1] - b2 * xy[3][1]);
147 
148  dxdz(0, 1) = 0.5 * (-a1 * xy[0][0] - a2 * xy[1][0] +
149  a2 * xy[2][0] + a1 * xy[3][0]);
150  dxdz(1, 1) = 0.5 * (-a1 * xy[0][1] - a2 * xy[1][1] +
151  a2 * xy[2][1] + a1 * xy[3][1]);
152 
153  dxdz.Invert();
154  ret.push_back(dxdz);
155  }
156  }
157  }
158  else if (geom->GetShapeType() == LibUtilities::eTriangle)
159  {
160  vector<Array<OneD, NekDouble> > xy;
161  for (int i = 0; i < geom->GetNumVerts(); i++)
162  {
163  Array<OneD, NekDouble> loc(2);
164  SpatialDomains::PointGeomSharedPtr p = geom->GetVertex(i);
165  p->GetCoords(loc);
166  xy.push_back(loc);
167  }
168 
170  Array<OneD, NekDouble> u = b[0]->GetZ();
171  Array<OneD, NekDouble> v = b[1]->GetZ();
172 
173  for (int i = 0; i < b[0]->GetNumPoints(); i++)
174  {
175  for (int j = 0; j < b[1]->GetNumPoints(); j++)
176  {
177  DNekMat dxdz(2, 2, 1.0, eFULL);
178  dxdz(0, 0) = -xy[0][0] / 2.0 + xy[1][0] / 2.0;
179 
180  dxdz(0, 1) = -xy[0][0] / 2.0 + xy[2][0] / 2.0;
181 
182  dxdz(1, 0) = -xy[0][1] / 2.0 + xy[1][1] / 2.0;
183 
184  dxdz(1, 1) = -xy[0][1] / 2.0 + xy[2][1] / 2.0;
185 
186  dxdz.Invert();
187  ret.push_back(dxdz);
188  }
189  }
190  }
191  else if (geom->GetShapeType() == LibUtilities::eTetrahedron)
192  {
193  vector<Array<OneD, NekDouble> > xyz;
194  for (int i = 0; i < geom->GetNumVerts(); i++)
195  {
196  Array<OneD, NekDouble> loc(3);
197  SpatialDomains::PointGeomSharedPtr p = geom->GetVertex(i);
198  p->GetCoords(loc);
199  xyz.push_back(loc);
200  }
201 
203  Array<OneD, NekDouble> u = b[0]->GetZ();
204  Array<OneD, NekDouble> v = b[1]->GetZ();
205  Array<OneD, NekDouble> z = b[2]->GetZ();
206 
207  for (int i = 0; i < b[0]->GetNumPoints(); i++)
208  {
209  for (int j = 0; j < b[1]->GetNumPoints(); j++)
210  {
211  for (int k = 0; k < b[2]->GetNumPoints(); k++)
212  {
213  DNekMat dxdz(3, 3, 1.0, eFULL);
214  dxdz(0, 0) = -xyz[0][0] / 2.0 + xyz[1][0] / 2.0;
215 
216  dxdz(0, 1) = -xyz[0][0] / 2.0 + xyz[2][0] / 2.0;
217 
218  dxdz(0, 2) = -xyz[0][0] / 2.0 + xyz[3][0] / 2.0;
219 
220  dxdz(1, 0) = -xyz[0][1] / 2.0 + xyz[1][1] / 2.0;
221 
222  dxdz(1, 1) = -xyz[0][1] / 2.0 + xyz[2][1] / 2.0;
223 
224  dxdz(1, 2) = -xyz[0][1] / 2.0 + xyz[3][1] / 2.0;
225 
226  dxdz(2, 0) = -xyz[0][2] / 2.0 + xyz[1][2] / 2.0;
227 
228  dxdz(2, 1) = -xyz[0][2] / 2.0 + xyz[2][2] / 2.0;
229 
230  dxdz(2, 2) = -xyz[0][2] / 2.0 + xyz[3][2] / 2.0;
231 
232  dxdz.Invert();
233  ret.push_back(dxdz);
234  }
235  }
236  }
237  }
238  else if (geom->GetShapeType() == LibUtilities::ePrism)
239  {
240  vector<Array<OneD, NekDouble> > xyz;
241  for (int i = 0; i < geom->GetNumVerts(); i++)
242  {
243  Array<OneD, NekDouble> loc(3);
244  SpatialDomains::PointGeomSharedPtr p = geom->GetVertex(i);
245  p->GetCoords(loc);
246  xyz.push_back(loc);
247  }
248 
250  Array<OneD, NekDouble> eta1 = b[0]->GetZ();
251  Array<OneD, NekDouble> eta2 = b[1]->GetZ();
252  Array<OneD, NekDouble> eta3 = b[2]->GetZ();
253 
254  for (int k = 0; k < b[2]->GetNumPoints(); k++)
255  {
256  for (int j = 0; j < b[1]->GetNumPoints(); j++)
257  {
258  for (int i = 0; i < b[0]->GetNumPoints(); i++)
259  {
260  NekDouble xi1 = 0.5 * (1 + eta1[i]) * (1 - eta3[k]) - 1.0;
261  NekDouble a2 = 0.5 * (1 + xi1);
262  NekDouble b1 = 0.5 * (1 - eta2[j]),
263  b2 = 0.5 * (1 + eta2[j]);
264  NekDouble c1 = 0.5 * (1 - eta3[k]),
265  c2 = 0.5 * (1 + eta3[k]);
266 
267  DNekMat dxdz(3, 3, 1.0, eFULL);
268 
269  dxdz(0, 0) = 0.5 * (-b1 * xyz[0][0] + b1 * xyz[1][0] +
270  b2 * xyz[2][0] - b2 * xyz[3][0]);
271  dxdz(1, 0) = 0.5 * (-b1 * xyz[0][1] + b1 * xyz[1][1] +
272  b2 * xyz[2][1] - b2 * xyz[3][1]);
273  dxdz(2, 0) = 0.5 * (-b1 * xyz[0][2] + b1 * xyz[1][2] +
274  b2 * xyz[2][2] - b2 * xyz[3][2]);
275 
276  dxdz(0, 1) = 0.5 * ((a2 - c1) * xyz[0][0] - a2 * xyz[1][0] +
277  a2 * xyz[2][0] + (c1 - a2) * xyz[3][0] -
278  c2 * xyz[4][0] + c2 * xyz[5][0]);
279  dxdz(1, 1) = 0.5 * ((a2 - c1) * xyz[0][1] - a2 * xyz[1][1] +
280  a2 * xyz[2][1] + (c1 - a2) * xyz[3][1] -
281  c2 * xyz[4][1] + c2 * xyz[5][1]);
282  dxdz(2, 1) = 0.5 * ((a2 - c1) * xyz[0][2] - a2 * xyz[1][2] +
283  a2 * xyz[2][2] + (c1 - a2) * xyz[3][2] -
284  c2 * xyz[4][2] + c2 * xyz[5][2]);
285 
286  dxdz(0, 2) = 0.5 * (-b1 * xyz[0][0] - b2 * xyz[3][0] +
287  b1 * xyz[4][0] + b2 * xyz[5][0]);
288  dxdz(1, 2) = 0.5 * (-b1 * xyz[0][1] - b2 * xyz[3][1] +
289  b1 * xyz[4][1] + b2 * xyz[5][1]);
290  dxdz(2, 2) = 0.5 * (-b1 * xyz[0][2] - b2 * xyz[3][2] +
291  b1 * xyz[4][2] + b2 * xyz[5][2]);
292 
293  dxdz.Invert();
294  ret.push_back(dxdz);
295  }
296  }
297  }
298  }
299  else if (geom->GetShapeType() == LibUtilities::eHexahedron)
300  {
301  vector<Array<OneD, NekDouble> > xyz;
302  for (int i = 0; i < geom->GetNumVerts(); i++)
303  {
304  Array<OneD, NekDouble> loc(3);
305  SpatialDomains::PointGeomSharedPtr p = geom->GetVertex(i);
306  p->GetCoords(loc);
307  xyz.push_back(loc);
308  }
309 
311  Array<OneD, NekDouble> eta1 = b[0]->GetZ();
312  Array<OneD, NekDouble> eta2 = b[1]->GetZ();
313  Array<OneD, NekDouble> eta3 = b[2]->GetZ();
314 
315  for (int k = 0; k < b[2]->GetNumPoints(); k++)
316  {
317  for (int j = 0; j < b[1]->GetNumPoints(); j++)
318  {
319  for (int i = 0; i < b[0]->GetNumPoints(); i++)
320  {
321  NekDouble a1 = 0.5 * (1 - eta1[i]);
322  NekDouble a2 = 0.5 * (1 + eta1[i]);
323  NekDouble b1 = 0.5 * (1 - eta2[j]),
324  b2 = 0.5 * (1 + eta2[j]);
325  NekDouble c1 = 0.5 * (1 - eta3[k]),
326  c2 = 0.5 * (1 + eta3[k]);
327 
328  DNekMat dxdz(3, 3, 1.0, eFULL);
329 
330  dxdz(0, 0) =
331  -0.5 * b1 * c1 * xyz[0][0] + 0.5 * b1 * c1 * xyz[1][0] +
332  0.5 * b2 * c1 * xyz[2][0] - 0.5 * b2 * c1 * xyz[3][0] -
333  0.5 * b1 * c2 * xyz[5][0] + 0.5 * b1 * c2 * xyz[5][0] +
334  0.5 * b2 * c2 * xyz[6][0] - 0.5 * b2 * c2 * xyz[7][0];
335  dxdz(1, 0) =
336  -0.5 * b1 * c1 * xyz[0][1] + 0.5 * b1 * c1 * xyz[1][1] +
337  0.5 * b2 * c1 * xyz[2][1] - 0.5 * b2 * c1 * xyz[3][1] -
338  0.5 * b1 * c2 * xyz[5][1] + 0.5 * b1 * c2 * xyz[5][1] +
339  0.5 * b2 * c2 * xyz[6][1] - 0.5 * b2 * c2 * xyz[7][1];
340  dxdz(2, 0) =
341  -0.5 * b1 * c1 * xyz[0][2] + 0.5 * b1 * c1 * xyz[1][2] +
342  0.5 * b2 * c1 * xyz[2][2] - 0.5 * b2 * c1 * xyz[3][2] -
343  0.5 * b1 * c2 * xyz[5][2] + 0.5 * b1 * c2 * xyz[5][2] +
344  0.5 * b2 * c2 * xyz[6][2] - 0.5 * b2 * c2 * xyz[7][2];
345 
346  dxdz(0, 1) =
347  -0.5 * a1 * c1 * xyz[0][0] - 0.5 * a2 * c1 * xyz[1][0] +
348  0.5 * a2 * c1 * xyz[2][0] + 0.5 * a1 * c1 * xyz[3][0] -
349  0.5 * a1 * c2 * xyz[5][0] - 0.5 * a2 * c2 * xyz[5][0] +
350  0.5 * a2 * c2 * xyz[6][0] + 0.5 * a1 * c2 * xyz[7][0];
351  dxdz(1, 1) =
352  -0.5 * a1 * c1 * xyz[0][1] - 0.5 * a2 * c1 * xyz[1][1] +
353  0.5 * a2 * c1 * xyz[2][1] + 0.5 * a1 * c1 * xyz[3][1] -
354  0.5 * a1 * c2 * xyz[5][1] - 0.5 * a2 * c2 * xyz[5][1] +
355  0.5 * a2 * c2 * xyz[6][1] + 0.5 * a1 * c2 * xyz[7][1];
356  dxdz(2, 1) =
357  -0.5 * a1 * c1 * xyz[0][2] - 0.5 * a2 * c1 * xyz[1][2] +
358  0.5 * a2 * c1 * xyz[2][2] + 0.5 * a1 * c1 * xyz[3][2] -
359  0.5 * a1 * c2 * xyz[5][2] - 0.5 * a2 * c2 * xyz[5][2] +
360  0.5 * a2 * c2 * xyz[6][2] + 0.5 * a1 * c2 * xyz[7][2];
361 
362  dxdz(0, 0) =
363  -0.5 * b1 * a1 * xyz[0][0] - 0.5 * b1 * a2 * xyz[1][0] -
364  0.5 * b2 * a2 * xyz[2][0] - 0.5 * b2 * a1 * xyz[3][0] +
365  0.5 * b1 * a1 * xyz[5][0] + 0.5 * b1 * a2 * xyz[5][0] +
366  0.5 * b2 * a2 * xyz[6][0] + 0.5 * b2 * a1 * xyz[7][0];
367  dxdz(1, 0) =
368  -0.5 * b1 * a1 * xyz[0][1] - 0.5 * b1 * a2 * xyz[1][1] -
369  0.5 * b2 * a2 * xyz[2][1] - 0.5 * b2 * a1 * xyz[3][1] +
370  0.5 * b1 * a1 * xyz[5][1] + 0.5 * b1 * a2 * xyz[5][1] +
371  0.5 * b2 * a2 * xyz[6][1] + 0.5 * b2 * a1 * xyz[7][1];
372  dxdz(2, 0) =
373  -0.5 * b1 * a1 * xyz[0][2] - 0.5 * b1 * a2 * xyz[1][2] -
374  0.5 * b2 * a2 * xyz[2][2] - 0.5 * b2 * a1 * xyz[3][2] +
375  0.5 * b1 * a1 * xyz[5][2] + 0.5 * b1 * a2 * xyz[5][2] +
376  0.5 * b2 * a2 * xyz[6][2] + 0.5 * b2 * a1 * xyz[7][2];
377 
378  dxdz.Invert();
379  ret.push_back(dxdz);
380  }
381  }
382  }
383  }
384  else
385  {
386  ASSERTL0(false, "not coded");
387  }
388 
389  return ret;
390 }
391 
394  bool s)
395 {
396  SpatialDomains::GeometrySharedPtr geom = e->GetGeom();
397  StdRegions::StdExpansionSharedPtr chi = e->GetGeom()->GetXmap();
398  LibUtilities::PointsKeyVector p = chi->GetPointsKeys();
399  LibUtilities::PointsKeyVector pElem = e->GetPointsKeys();
400  SpatialDomains::GeomFactorsSharedPtr gfac = geom->GetGeomFactors();
401  const int expDim = chi->GetNumBases();
402  int nElemPts = 1;
403 
404  vector<LibUtilities::BasisKey> basisKeys;
405  bool needsInterp = false;
406 
407  for (int i = 0; i < expDim; ++i)
408  {
409  nElemPts *= pElem[i].GetNumPoints();
410  needsInterp =
411  needsInterp || pElem[i].GetNumPoints() < p[i].GetNumPoints() - 1;
412  }
413 
414  if (needsInterp)
415  {
416  stringstream err;
417  err << "Interpolating from higher order geometry to lower order in "
418  << "element " << geom->GetGlobalID();
419  NEKERROR(ErrorUtil::ewarning, err.str());
420  }
421 
422  for (int i = 0; i < expDim; ++i)
423  {
424  basisKeys.push_back(
425  needsInterp
426  ? chi->GetBasis(i)->GetBasisKey()
427  : LibUtilities::BasisKey(chi->GetBasisType(i),
428  chi->GetBasisNumModes(i), pElem[i]));
429  }
430 
432  switch (chi->DetShapeType())
433  {
436  basisKeys[0], basisKeys[1]);
437  break;
440  basisKeys[0], basisKeys[1]);
441  break;
444  basisKeys[0], basisKeys[1], basisKeys[2]);
445  break;
448  basisKeys[0], basisKeys[1], basisKeys[2]);
449  break;
452  basisKeys[0], basisKeys[1], basisKeys[2]);
453  break;
454  default:
455  ASSERTL0(false, "nope");
456  }
457 
458  SpatialDomains::DerivStorage deriv = gfac->GetDeriv(pElem);
459 
460  const int pts = deriv[0][0].num_elements();
461  const int nq = chiMod->GetTotPoints();
462 
463  ASSERTL0(pts == nq, "what");
464 
465  vector<DNekMat> i2rm = MappingIdealToRef(geom, chiMod);
466  Array<OneD, NekDouble> eta(nq);
467 
468  for (int k = 0; k < pts; ++k)
469  {
470  DNekMat jac(expDim, expDim, 0.0, eFULL);
471  DNekMat jacIdeal(expDim, expDim, 0.0, eFULL);
472 
473  for (int i = 0; i < expDim; ++i)
474  {
475  for (int j = 0; j < expDim; ++j)
476  {
477  jac(j, i) = deriv[i][j][k];
478  }
479  }
480 
481  jacIdeal = jac * i2rm[k];
482  NekDouble jacDet;
483 
484  if (expDim == 2)
485  {
486  jacDet = jacIdeal(0, 0) * jacIdeal(1, 1) -
487  jacIdeal(0, 1) * jacIdeal(1, 0);
488  }
489  else if (expDim == 3)
490  {
491  jacDet = jacIdeal(0, 0) * (jacIdeal(1, 1) * jacIdeal(2, 2) -
492  jacIdeal(2, 1) * jacIdeal(1, 2)) -
493  jacIdeal(0, 1) * (jacIdeal(1, 0) * jacIdeal(2, 2) -
494  jacIdeal(2, 0) * jacIdeal(1, 2)) +
495  jacIdeal(0, 2) * (jacIdeal(1, 0) * jacIdeal(2, 1) -
496  jacIdeal(2, 0) * jacIdeal(1, 1));
497  }
498  else
499  {
500  ASSERTL0(false, "silly exp dim");
501  }
502 
503  if(s)
504  {
505  eta[k] = jacDet;
506  }
507  else
508  {
509  NekDouble frob = 0.0;
510 
511  for (int i = 0; i < expDim; ++i)
512  {
513  for (int j = 0; j < expDim; ++j)
514  {
515  frob += jacIdeal(i,j) * jacIdeal(i,j);
516  }
517  }
518 
519  NekDouble sigma = 0.5*(jacDet + sqrt(jacDet*jacDet));
520  eta[k] = expDim * pow(sigma, 2.0/expDim) / frob;
521  }
522  }
523 
524  if(s)
525  {
526  NekDouble mx = -1.0 * numeric_limits<double>::max();
527  NekDouble mn = numeric_limits<double>::max();
528  for(int k = 0; k < pts; k++)
529  {
530  mx = max(mx,eta[k]);
531  mn = min(mn,eta[k]);
532  }
533  for(int k = 0; k < pts; k++)
534  {
535  eta[k] = mn/mx;
536  }
537  }
538 
539  // Project onto output stuff
540  if (needsInterp && pts != 1)
541  {
542  Array<OneD, NekDouble> tmp(nElemPts);
543 
544  if (expDim == 2)
545  {
546  LibUtilities::Interp2D(p[0], p[1], eta, pElem[0], pElem[1], tmp);
547  }
548  else if (expDim == 3)
549  {
550  LibUtilities::Interp3D(p[0], p[1], p[2], eta, pElem[0], pElem[1],
551  pElem[2], tmp);
552  }
553  else
554  {
555  ASSERTL0(false, "mesh dim makes no sense");
556  }
557 
558  eta = tmp;
559  }
560 
561  if (pts == 1)
562  {
563  Vmath::Fill(nq - 1, eta[0], &eta[1], 1);
564  }
565 
566  return eta;
567 }
568 }
569 }
map< string, ConfigOption > m_config
List of configuration values.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
vector< DNekMat > MappingIdealToRef(SpatialDomains::GeometrySharedPtr geom, StdRegions::StdExpansionSharedPtr chi)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:191
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:242
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
Represents a command-line configuration option.
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:46
virtual void Process(po::variables_map &vm)
Write mesh to output file.
STL namespace.
pair< ModuleType, string > ModuleKey
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis...
Definition: Interp.cpp:116
boost::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:740
double NekDouble
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
void Interp3D(const BasisKey &fbasis0, const BasisKey &fbasis1, const BasisKey &fbasis2, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, const BasisKey &tbasis2, Array< OneD, NekDouble > &to)
this function interpolates a 3D function evaluated at the quadrature points of the 3D basis...
Definition: Interp.cpp:186
boost::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
Abstract base class for processing modules.
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
Array< OneD, NekDouble > GetQ(LocalRegions::ExpansionSharedPtr e, bool s)
boost::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:53
Describes the specification for a Basis.
Definition: Basis.h:50
boost::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:60
ModuleFactory & GetModuleFactory()
FieldSharedPtr m_f
Field object.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215