Nektar++
ProcessInterpPoints.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: ProcessInterpPoints.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: Interpolate field to a series of specified points.
32 //
33 ////////////////////////////////////////////////////////////////////////////////
34 
35 #include <iostream>
36 #include <string>
37 using namespace std;
38 
39 #include <boost/core/ignore_unused.hpp>
40 #include <boost/geometry.hpp>
41 #include <boost/lexical_cast.hpp>
42 #include <boost/math/special_functions/fpclassify.hpp>
43 
50 
51 #include "ProcessInterpPoints.h"
52 
53 namespace bg = boost::geometry;
54 namespace bgi = boost::geometry::index;
55 
56 namespace Nektar
57 {
58 namespace FieldUtils
59 {
60 
61 ModuleKey ProcessInterpPoints::className =
63  ModuleKey(eProcessModule, "interppoints"), ProcessInterpPoints::create,
64  "Interpolates a field to a set of points. Requires fromfld, fromxml "
65  "to be defined, and a topts, line, plane or block of target points ");
66 
67 ProcessInterpPoints::ProcessInterpPoints(FieldSharedPtr f) : ProcessModule(f)
68 {
69  m_config["fromxml"] = ConfigOption(
70  false, "NotSet", "Xml file from which to interpolate field");
71 
72  m_config["fromfld"] = ConfigOption(
73  false, "NotSet", "Fld file from which to interpolate field");
74 
75  m_config["topts"] =
76  ConfigOption(false, "NotSet", "Pts file to which interpolate field");
77  m_config["line"] = ConfigOption(false, "NotSet",
78  "Specify a line of N points using "
79  "line=N,x0,y0,z0,z1,y1,z1");
80  m_config["plane"] =
81  ConfigOption(false, "NotSet",
82  "Specify a plane of N1 x N2 points using "
83  "plane=N1,N2,x0,y0,z0,z1,y1,z1,x2,y2,z2,x3,y3,z3");
84  m_config["box"] =
85  ConfigOption(false, "NotSet",
86  "Specify a rectangular box of N1 x N2 x N3 points "
87  "using a box of points limited by box="
88  "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
89 
90  m_config["clamptolowervalue"] =
91  ConfigOption(false, "-10000000", "Lower bound for interpolation value");
92  m_config["clamptouppervalue"] =
93  ConfigOption(false, "10000000", "Upper bound for interpolation value");
94  m_config["defaultvalue"] =
95  ConfigOption(false, "0", "Default value if point is outside domain");
96 
97  m_config["cp"] =
98  ConfigOption(false, "NotSet",
99  "Parameters p0 and q to determine pressure coefficients");
100  m_config["realmodetoimag"] =
101  ConfigOption(false, "NotSet", "Take fields as sin mode");
102 }
103 
105 {
106 }
107 
108 void ProcessInterpPoints::Process(po::variables_map &vm)
109 {
110  m_f->SetUpExp(vm);
111 
112  CreateFieldPts(vm);
113 
114  FieldSharedPtr fromField = std::shared_ptr<Field>(new Field());
115  std::vector<std::string> files;
116  ParseUtils::GenerateVector(m_config["fromxml"].as<string>(), files);
117 
118  // set up session file for from field
119  char *argv[] = {const_cast<char *>("FieldConvert"), nullptr};
120  fromField->m_session = LibUtilities::SessionReader::CreateInstance(
121  1, argv, files,
122  LibUtilities::GetCommFactory().CreateInstance("Serial", 0, 0));
123 
124  // Set up range based on min and max of local parallel partition
127 
128  int coordim = m_f->m_fieldPts->GetDim();
129  int npts = m_f->m_fieldPts->GetNpoints();
130  std::vector<std::string> fieldNames = m_f->m_fieldPts->GetFieldNames();
131  for (auto &it : fieldNames)
132  {
133  m_f->m_fieldPts->RemoveField(it);
134  }
135 
137  m_f->m_fieldPts->GetPts(pts);
138 
139  rng->m_checkShape = false;
140  rng->m_zmin = -1;
141  rng->m_zmax = 1;
142  rng->m_ymin = -1;
143  rng->m_ymax = 1;
144  switch (coordim)
145  {
146  case 3:
147  rng->m_doZrange = true;
148  rng->m_zmin = Vmath::Vmin(npts, pts[2], 1);
149  rng->m_zmax = Vmath::Vmax(npts, pts[2], 1);
150  if (rng->m_zmax == rng->m_zmin)
151  {
152  rng->m_zmin -= 1;
153  rng->m_zmax += 1;
154  }
155  /* Falls through. */
156  case 2:
157  rng->m_doYrange = true;
158  rng->m_ymin = Vmath::Vmin(npts, pts[1], 1);
159  rng->m_ymax = Vmath::Vmax(npts, pts[1], 1);
160  /* Falls through. */
161  case 1:
162  rng->m_doXrange = true;
163  rng->m_xmin = Vmath::Vmin(npts, pts[0], 1);
164  rng->m_xmax = Vmath::Vmax(npts, pts[0], 1);
165  break;
166  default:
167  NEKERROR(ErrorUtil::efatal, "Too many values specified in range");
168  }
169 
170  // setup rng parameters.
171  fromField->m_graph =
172  SpatialDomains::MeshGraph::Read(fromField->m_session, rng);
173 
174  // Read in local from field partitions
175  const SpatialDomains::ExpansionInfoMap &expansions =
176  fromField->m_graph->GetExpansionInfo();
177  Array<OneD, int> ElementGIDs(expansions.size());
178 
179  int i = 0;
180  for (auto &expIt : expansions)
181  {
182  ElementGIDs[i++] = expIt.second->m_geomShPtr->GetGlobalID();
183  }
184  // check to see that we do have some element in the domain since
185  // possibly all points could be outside of the domain
186  ASSERTL0(i > 0, "No elements are set. Are the interpolated points "
187  "within the domain given by the xml files?");
188  string fromfld = m_config["fromfld"].as<string>();
189  m_f->FieldIOForFile(fromfld)->Import(
190  fromfld, fromField->m_fielddef, fromField->m_data,
192  int NumHomogeneousDir = fromField->m_fielddef[0]->m_numHomogeneousDir;
193  for (i = 0; i < fromField->m_fielddef.size(); ++i)
194  {
195  int d1 = fromField->m_fielddef[i]->m_basis.size();
196  d1 -= 1;
197  if (d1 >= 0 && (fromField->m_fielddef[i]->m_basis[d1] ==
199  fromField->m_fielddef[i]->m_basis[d1] ==
201  {
202  fromField->m_fielddef[i]->m_homogeneousZIDs[0] += 2;
203  fromField->m_fielddef[i]->m_numModes[d1] = 4;
204  fromField->m_fielddef[i]->m_basis[d1] = LibUtilities::eFourier;
205  }
206  }
207 
208  //----------------------------------------------
209  // Set up Expansion information to use mode order from field
210  fromField->m_graph->SetExpansionInfo(fromField->m_fielddef);
211  int nfields = fromField->m_fielddef[0]->m_fields.size();
212  fromField->m_exp.resize(nfields);
213  fromField->m_exp[0] = fromField->SetUpFirstExpList(NumHomogeneousDir, true);
214  m_f->m_exp.resize(nfields);
215 
216  // declare auxiliary fields.
217  for (i = 1; i < nfields; ++i)
218  {
219  fromField->m_exp[i] = fromField->AppendExpList(NumHomogeneousDir);
220  }
221  // load field into expansion in fromfield.
222  set<int> sinmode;
223  if (m_config["realmodetoimag"].as<string>().compare("NotSet"))
224  {
225  vector<int> value;
227  m_config["realmodetoimag"].as<string>(), value),
228  "Failed to interpret realmodetoimag string");
229  for (int j : value)
230  {
231  sinmode.insert(j);
232  }
233  }
234  for (int j = 0; j < nfields; ++j)
235  {
236  for (i = 0; i < fromField->m_fielddef.size(); i++)
237  {
238  fromField->m_exp[j]->ExtractDataToCoeffs(
239  fromField->m_fielddef[i], fromField->m_data[i],
240  fromField->m_fielddef[0]->m_fields[j],
241  fromField->m_exp[j]->UpdateCoeffs());
242  }
243  if (NumHomogeneousDir == 1)
244  {
245  fromField->m_exp[j]->SetWaveSpace(true);
246  if (sinmode.count(j))
247  {
248  int Ncoeff = fromField->m_exp[j]->GetPlane(2)->GetNcoeffs();
249  Vmath::Smul(
250  Ncoeff, -1., fromField->m_exp[j]->GetPlane(2)->GetCoeffs(),
251  1, fromField->m_exp[j]->GetPlane(3)->UpdateCoeffs(), 1);
252  Vmath::Zero(Ncoeff,
253  fromField->m_exp[j]->GetPlane(2)->UpdateCoeffs(),
254  1);
255  }
256  }
257  fromField->m_exp[j]->BwdTrans(fromField->m_exp[j]->GetCoeffs(),
258  fromField->m_exp[j]->UpdatePhys());
259  Array<OneD, NekDouble> newPts(m_f->m_fieldPts->GetNpoints());
260  m_f->m_fieldPts->AddField(newPts,
261  fromField->m_fielddef[0]->m_fields[j]);
262  m_f->m_variables.push_back(fromField->m_fielddef[0]->m_fields[j]);
263  }
264 
265  NekDouble clamp_low = m_config["clamptolowervalue"].as<NekDouble>();
266  NekDouble clamp_up = m_config["clamptouppervalue"].as<NekDouble>();
267  NekDouble def_value = m_config["defaultvalue"].as<NekDouble>();
268 
269  InterpolateFieldToPts(fromField->m_exp, m_f->m_fieldPts, clamp_low,
270  clamp_up, def_value);
271 
272  if (!boost::iequals(m_config["cp"].as<string>(), "NotSet"))
273  {
274  calcCp0();
275  }
276 }
277 
278 void ProcessInterpPoints::CreateFieldPts(po::variables_map &vm)
279 {
280  boost::ignore_unused(vm);
281 
282  int rank = m_f->m_comm->GetRank();
283  int nprocs = m_f->m_comm->GetSize();
284  // Check for command line point specification
285 
286  if (m_config["topts"].as<string>().compare("NotSet") != 0)
287  {
288  string inFile = m_config["topts"].as<string>();
289 
290  if (boost::filesystem::path(inFile).extension() == ".pts")
291  {
294  m_f->m_comm);
295 
296  ptsIO->Import(inFile, m_f->m_fieldPts);
297  }
298  else if (boost::filesystem::path(inFile).extension() == ".csv")
299  {
302  m_f->m_comm);
303 
304  csvIO->Import(inFile, m_f->m_fieldPts);
305  }
306  else
307  {
308  ASSERTL0(false, "unknown topts file type");
309  }
310  }
311  else if (m_config["line"].as<string>().compare("NotSet") != 0)
312  {
313  vector<NekDouble> values;
314  ASSERTL0(
315  ParseUtils::GenerateVector(m_config["line"].as<string>(), values),
316  "Failed to interpret line string");
317 
318  ASSERTL0(values.size() > 2, "line string should contain 2*Dim+1 values "
319  "N,x0,y0,z0,x1,y1,z1");
320 
321  double tmp;
322  ASSERTL0(std::modf(values[0], &tmp) == 0.0, "N is not an integer");
323  ASSERTL0(values[0] > 1, "N is not a valid number");
324 
325  int dim = (values.size() - 1) / 2;
326  int npts = values[0];
327 
328  // Information for partitioning
329  int ptsPerProc = npts / nprocs;
330  int extraPts = (rank < nprocs - 1) ? 0 : npts % nprocs;
331  int locPts = ptsPerProc + extraPts;
332  int start = rank * ptsPerProc;
333  int end = start + locPts;
334 
336  Array<OneD, NekDouble> delta(dim);
337  for (int i = 0; i < dim; ++i)
338  {
339  pts[i] = Array<OneD, NekDouble>(locPts);
340  delta[i] = (values[dim + i + 1] - values[i + 1]) / (npts - 1);
341  }
342 
343  for (int i = 0, cntLoc = 0; i < npts; ++i)
344  {
345  if (i >= start && i < end)
346  {
347  for (int n = 0; n < dim; ++n)
348  {
349  pts[n][cntLoc] = values[n + 1] + i * delta[n];
350  }
351  ++cntLoc;
352  }
353  }
354 
355  vector<size_t> ppe;
356  ppe.push_back(npts);
357  m_f->m_fieldPts =
359  m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsLine);
360  m_f->m_fieldPts->SetPointsPerEdge(ppe);
361  }
362  else if (m_config["plane"].as<string>().compare("NotSet") != 0)
363  {
364  vector<NekDouble> values;
365  ASSERTL0(
366  ParseUtils::GenerateVector(m_config["plane"].as<string>(), values),
367  "Failed to interpret plane string");
368 
369  ASSERTL0(values.size() > 9,
370  "plane string should contain 4 Dim+2 values "
371  "N1,N2,x0,y0,z0,x1,y1,z1,x2,y2,z2,x3,y3,z3");
372 
373  double tmp;
374  ASSERTL0(std::modf(values[0], &tmp) == 0.0, "N1 is not an integer");
375  ASSERTL0(std::modf(values[1], &tmp) == 0.0, "N2 is not an integer");
376 
377  ASSERTL0(values[0] > 1, "N1 is not a valid number");
378  ASSERTL0(values[1] > 1, "N2 is not a valid number");
379 
380  int dim = (values.size() - 2) / 4;
381 
382  Array<OneD, int> npts(2);
383  npts[0] = values[0];
384  npts[1] = values[1];
385 
386  int totpts = npts[0] * npts[1];
387 
388  // Information for partitioning
389  int ptsPerProc = totpts / nprocs;
390  int extraPts = (rank < nprocs - 1) ? 0 : totpts % nprocs;
391  int locPts = ptsPerProc + extraPts;
392  int start = rank * ptsPerProc;
393  int end = start + locPts;
394 
396  Array<OneD, NekDouble> delta1(dim);
397  Array<OneD, NekDouble> delta2(dim);
398  for (int i = 0; i < dim; ++i)
399  {
400  pts[i] = Array<OneD, NekDouble>(locPts);
401  delta1[i] = (values[2 + 1 * dim + i] - values[2 + 0 * dim + i]) /
402  (npts[0] - 1);
403  delta2[i] = (values[2 + 2 * dim + i] - values[2 + 3 * dim + i]) /
404  (npts[0] - 1);
405  }
406 
407  for (int j = 0, cnt = 0, cntLoc = 0; j < npts[1]; ++j)
408  {
409  for (int i = 0; i < npts[0]; ++i, ++cnt)
410  {
411  if (cnt >= start && cnt < end)
412  {
413  for (int n = 0; n < dim; ++n)
414  {
415  pts[n][cntLoc] =
416  (values[2 + n] + i * delta1[n]) *
417  (1.0 - j / ((NekDouble)(npts[1] - 1))) +
418  (values[2 + 3 * dim + n] + i * delta2[n]) *
419  (j / ((NekDouble)(npts[1] - 1)));
420  }
421  ++cntLoc;
422  }
423  }
424  }
425 
426  vector<size_t> ppe;
427  ppe.push_back(npts[0]);
428  ppe.push_back(npts[1]);
429  m_f->m_fieldPts =
431  m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsPlane);
432  m_f->m_fieldPts->SetPointsPerEdge(ppe);
433  }
434  else if (m_config["box"].as<string>().compare("NotSet") != 0)
435  {
436  vector<NekDouble> values;
437  ASSERTL0(
438  ParseUtils::GenerateVector(m_config["box"].as<string>(), values),
439  "Failed to interpret box string");
440 
441  ASSERTL0(values.size() == 9, "box string should contain 9 values "
442  "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
443 
444  int dim = 3;
445  Array<OneD, int> npts(3);
446  npts[0] = values[0];
447  npts[1] = values[1];
448  npts[2] = values[2];
449 
450  int totpts = npts[0] * npts[1] * npts[2];
451 
453  Array<OneD, NekDouble> delta(dim);
454 
455  // Information for partitioning
456  int ptsPerProc = totpts / nprocs;
457  int extraPts = (rank < nprocs - 1) ? 0 : totpts % nprocs;
458  int locPts = ptsPerProc + extraPts;
459  int start = rank * ptsPerProc;
460  int end = start + locPts;
461 
462  for (int i = 0; i < dim; ++i)
463  {
464  pts[i] = Array<OneD, NekDouble>(locPts);
465  delta[i] = (values[4 + 2 * i] - values[3 + 2 * i]) / (npts[i] - 1);
466  }
467 
468  for (int k = 0, cnt = 0, cntLoc = 0; k < npts[2]; ++k)
469  {
470  for (int j = 0; j < npts[1]; ++j)
471  {
472  for (int i = 0; i < npts[0]; ++i, ++cnt)
473  {
474  if (cnt >= start && cnt < end)
475  {
476  pts[0][cntLoc] = values[3] + i * delta[0];
477  pts[1][cntLoc] = values[5] + j * delta[1];
478  pts[2][cntLoc] = values[7] + k * delta[2];
479  ++cntLoc;
480  }
481  }
482  }
483  }
484 
485  vector<size_t> ppe;
486  ppe.push_back(npts[0]);
487  ppe.push_back(npts[1]);
488  ppe.push_back(npts[2]);
489  m_f->m_fieldPts =
491  m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsBox);
492  m_f->m_fieldPts->SetPointsPerEdge(ppe);
493  vector<NekDouble> boxdim;
494  boxdim.assign(&values[3], &values[3] + 6);
495  m_f->m_fieldPts->SetBoxSize(boxdim);
496  }
497  else
498  {
499  ASSERTL0(false, "Missing target points for ProcessInterpPoints.");
500  }
501 }
502 
504  vector<MultiRegions::ExpListSharedPtr> &field0,
506  NekDouble clamp_up, NekDouble def_value)
507 {
508  boost::ignore_unused(def_value);
509 
510  ASSERTL0(pts->GetNFields() == field0.size(), "ptField has too few fields");
511 
512  int nfields = field0.size();
513 
514  Interpolator interp;
515  if (m_f->m_comm->GetRank() == 0)
516  {
518  this);
519  }
520  interp.Interpolate(field0, pts);
521  if (m_f->m_comm->GetRank() == 0)
522  {
523  cout << endl;
524  }
525 
526  for (int f = 0; f < nfields; ++f)
527  {
528  for (int i = 0; i < pts->GetNpoints(); ++i)
529  {
530  if (pts->GetPointVal(f, i) > clamp_up)
531  {
532  pts->SetPointVal(f, i, clamp_up);
533  }
534  else if (pts->GetPointVal(f, i) < clamp_low)
535  {
536  pts->SetPointVal(f, i, clamp_low);
537  }
538  }
539  }
540 }
541 
543 {
544  LibUtilities::PtsFieldSharedPtr pts = m_f->m_fieldPts;
545  int dim = pts->GetDim();
546  int nq1 = pts->GetNpoints();
547  int r, f;
548  int pfield = -1;
549  NekDouble p0, qinv;
550  vector<int> velid;
551 
552  vector<NekDouble> values;
553  ASSERTL0(ParseUtils::GenerateVector(m_config["cp"].as<string>(), values),
554  "Failed to interpret cp string");
555 
556  ASSERTL0(values.size() == 2, "cp string should contain 2 values "
557  "p0 and q (=1/2 rho u^2)");
558 
559  p0 = values[0];
560  qinv = 1.0 / values[1];
561 
562  for (int i = 0; i < pts->GetNFields(); ++i)
563  {
564  if (boost::iequals(pts->GetFieldName(i), "p"))
565  {
566  pfield = i;
567  }
568 
569  if (boost::iequals(pts->GetFieldName(i), "u") ||
570  boost::iequals(pts->GetFieldName(i), "v") ||
571  boost::iequals(pts->GetFieldName(i), "w"))
572  {
573  velid.push_back(i);
574  }
575  }
576 
577  if (pfield != -1)
578  {
579  if (!velid.size())
580  {
581  WARNINGL0(false, "Did not find velocity components for Cp0");
582  }
583  }
584  else
585  {
586  WARNINGL0(false, "Failed to find 'p' field to determine cp0");
587  }
588 
589  // Allocate data storage
591 
592  for (f = 0; f < 2; ++f)
593  {
594  data[f] = Array<OneD, NekDouble>(nq1, 0.0);
595  }
596 
597  for (r = 0; r < nq1; r++)
598  {
599  if (pfield != -1) // calculate cp
600  {
601  data[0][r] = qinv * (pts->GetPointVal(dim + pfield, r) - p0);
602 
603  if (velid.size()) // calculate cp0
604  {
605  NekDouble q = 0;
606  for (int i = 0; i < velid.size(); ++i)
607  {
608  q += 0.5 * pts->GetPointVal(dim + velid[i], r) *
609  pts->GetPointVal(dim + velid[i], r);
610  }
611  data[1][r] =
612  qinv * (pts->GetPointVal(dim + pfield, r) + q - p0);
613  }
614  }
615  }
616 
617  if (pfield != -1)
618  {
619  pts->AddField(data[0], "Cp");
620  m_f->m_variables.push_back("Cp");
621  if (velid.size())
622  {
623  pts->AddField(data[1], "Cp0");
624  m_f->m_variables.push_back("Cp0");
625  }
626  }
627 }
628 
630  const int goal) const
631 {
632  LibUtilities::PrintProgressbar(position, goal, "Interpolating");
633 }
634 } // namespace FieldUtils
635 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:222
A class that contains algorithms for interpolation between pts fields, expansions and different meshe...
FIELD_UTILS_EXPORT void Interpolate(const std::vector< MultiRegions::ExpListSharedPtr > expInField, std::vector< MultiRegions::ExpListSharedPtr > &expOutField, NekDouble def_value=0.0)
Interpolate from an expansion to an expansion.
FieldSharedPtr m_f
Field object.
Definition: Module.h:225
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:228
void InterpolateFieldToPts(std::vector< MultiRegions::ExpListSharedPtr > &field0, LibUtilities::PtsFieldSharedPtr &pts, NekDouble clamp_low, NekDouble clamp_up, NekDouble def_value)
virtual void Process(po::variables_map &vm)
Write mesh to output file.
void PrintProgressbar(const int position, const int goal) const
Abstract base class for processing modules.
Definition: Module.h:260
void SetProgressCallback(FuncPointerT func, ObjectPointerT obj)
sets a callback funtion which gets called every time the interpolation progresses
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
static SessionReaderSharedPtr CreateInstance(int argc, char *argv[])
Creates an instance of the SessionReader class.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:129
static MeshGraphSharedPtr Read(const LibUtilities::SessionReaderSharedPtr pSession, LibUtilities::DomainRangeShPtr rng=LibUtilities::NullDomainRangeShPtr, bool fillGraph=true)
Definition: MeshGraph.cpp:110
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:989
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:285
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:49
int PrintProgressbar(const int position, const int goal, const std::string message, int lastprogress=-1)
Prints a progressbar.
Definition: Progressbar.hpp:67
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:53
std::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:190
std::shared_ptr< DomainRange > DomainRangeShPtr
Definition: DomainRange.h:66
std::shared_ptr< CsvIO > CsvIOSharedPtr
Definition: CsvIO.h:83
std::shared_ptr< PtsIO > PtsIOSharedPtr
Definition: PtsIO.h:97
CommFactory & GetCommFactory()
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode .
Definition: BasisType.h:70
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode .
Definition: BasisType.h:68
@ eFourier
Fourier Expansion .
Definition: BasisType.h:57
std::map< int, ExpansionInfoShPtr > ExpansionInfoMap
Definition: MeshGraph.h:143
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
Definition: Vmath.cpp:1050
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:248
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
Definition: Vmath.cpp:945
Represents a command-line configuration option.
Definition: Module.h:131