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"),
64  ProcessInterpPoints::create,
65  "Interpolates a field to a set of points. Requires fromfld, fromxml "
66  "to be defined, and a topts, line, plane or block of target points ");
67 
68 ProcessInterpPoints::ProcessInterpPoints(FieldSharedPtr f) : ProcessModule(f)
69 {
70  m_config["fromxml"] = ConfigOption(
71  false, "NotSet", "Xml file from which to interpolate field");
72 
73  m_config["fromfld"] = ConfigOption(
74  false, "NotSet", "Fld file from which to interpolate field");
75 
76  m_config["topts"] = ConfigOption(
77  false, "NotSet", "Pts file to which interpolate field");
78  m_config["line"] = ConfigOption(
79  false, "NotSet", "Specify a line of N points using "
80  "line=N,x0,y0,z0,z1,y1,z1");
81  m_config["plane"] = ConfigOption(
82  false, "NotSet", "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"] = ConfigOption(
85  false, "NotSet", "Specify a rectangular box of N1 x N2 x N3 points "
86  "using a box of points limited by box="
87  "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
88 
89  m_config["clamptolowervalue"] =
90  ConfigOption(false, "-10000000", "Lower bound for interpolation value");
91  m_config["clamptouppervalue"] =
92  ConfigOption(false, "10000000", "Upper bound for interpolation value");
93  m_config["defaultvalue"] =
94  ConfigOption(false, "0", "Default value if point is outside domain");
95 
96  m_config["cp"] =
97  ConfigOption(false, "NotSet",
98  "Parameters p0 and q to determine pressure coefficients");
99  m_config["realmodetoimag"] = ConfigOption(
100  false, "NotSet", "Take fields as sin mode");
101 }
102 
104 {
105 }
106 
107 void ProcessInterpPoints::Process(po::variables_map &vm)
108 {
109  m_f->SetUpExp(vm);
110 
111  CreateFieldPts(vm);
112 
113  FieldSharedPtr fromField = std::shared_ptr<Field>(new Field());
114  std::vector<std::string> files;
115  ParseUtils::GenerateVector(m_config["fromxml"].as<string>(), files);
116 
117  // set up session file for from field
118  char *argv[] = { const_cast<char *>("FieldConvert"), nullptr };
119  fromField->m_session =
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 &&
198  (fromField->m_fielddef[i]->m_basis[d1] ==
200  fromField->m_fielddef[i]->m_basis[d1] ==
202  {
203  fromField->m_fielddef[i]->m_homogeneousZIDs[0] += 2;
204  fromField->m_fielddef[i]->m_numModes[d1] = 4;
205  fromField->m_fielddef[i]->m_basis[d1] = LibUtilities::eFourier;
206  }
207  }
208 
209  //----------------------------------------------
210  // Set up Expansion information to use mode order from field
211  fromField->m_graph->SetExpansionInfo(fromField->m_fielddef);
212  int nfields = fromField->m_fielddef[0]->m_fields.size();
213  fromField->m_exp.resize(nfields);
214  fromField->m_exp[0] = fromField->SetUpFirstExpList(NumHomogeneousDir, true);
215  m_f->m_exp.resize(nfields);
216 
217  // declare auxiliary fields.
218  for (i = 1; i < nfields; ++i)
219  {
220  fromField->m_exp[i] = fromField->AppendExpList(NumHomogeneousDir);
221  }
222  // load field into expansion in fromfield.
223  set<int> sinmode;
224  if (m_config["realmodetoimag"].as<string>().compare("NotSet"))
225  {
226  vector<int> value;
228  m_config["realmodetoimag"].as<string>(), value),
229  "Failed to interpret realmodetoimag string");
230  for (int j: value)
231  {
232  sinmode.insert(j);
233  }
234  }
235  for (int j = 0; j < nfields; ++j)
236  {
237  for (i = 0; i < fromField->m_fielddef.size(); i++)
238  {
239  fromField->m_exp[j]->ExtractDataToCoeffs(
240  fromField->m_fielddef[i], fromField->m_data[i],
241  fromField->m_fielddef[0]->m_fields[j],
242  fromField->m_exp[j]->UpdateCoeffs());
243  }
244  if (NumHomogeneousDir == 1)
245  {
246  fromField->m_exp[j]->SetWaveSpace(true);
247  if (sinmode.count(j))
248  {
249  int Ncoeff = fromField->m_exp[j]->GetPlane(2)->GetNcoeffs();
250  Vmath::Smul(Ncoeff, -1.,
251  fromField->m_exp[j]->GetPlane(2)->GetCoeffs() , 1,
252  fromField->m_exp[j]->GetPlane(3)->UpdateCoeffs(), 1);
253  Vmath::Zero(Ncoeff,
254  fromField->m_exp[j]->GetPlane(2)->UpdateCoeffs(), 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  }
312  else if (m_config["line"].as<string>().compare("NotSet") != 0)
313  {
314  vector<NekDouble> values;
316  m_config["line"].as<string>(), values),
317  "Failed to interpret line string");
318 
319  ASSERTL0(values.size() > 2,
320  "line string should contain 2*Dim+1 values "
321  "N,x0,y0,z0,x1,y1,z1");
322 
323  double tmp;
324  ASSERTL0(std::modf(values[0], &tmp) == 0.0, "N is not an integer");
325  ASSERTL0(values[0] > 1, "N is not a valid number");
326 
327  int dim = (values.size() - 1) / 2;
328  int npts = values[0];
329 
330  // Information for partitioning
331  int ptsPerProc = npts / nprocs;
332  int extraPts = (rank < nprocs - 1) ? 0: npts % nprocs;
333  int locPts = ptsPerProc + extraPts;
334  int start = rank * ptsPerProc;
335  int end = start + locPts;
336 
338  Array<OneD, NekDouble> delta(dim);
339  for (int i = 0; i < dim; ++i)
340  {
341  pts[i] = Array<OneD, NekDouble>(locPts);
342  delta[i] = (values[dim + i + 1] - values[ i + 1]) / (npts - 1);
343  }
344 
345 
346 
347  for (int i = 0, cntLoc = 0; i < npts; ++i)
348  {
349  if (i >= start && i < end)
350  {
351  for (int n = 0; n < dim; ++n)
352  {
353  pts[n][cntLoc] = values[n+1] + i * delta[n];
354  }
355  ++cntLoc;
356  }
357  }
358 
359  vector<size_t> ppe;
360  ppe.push_back(npts);
361  m_f->m_fieldPts =
363  pts);
364  m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsLine);
365  m_f->m_fieldPts->SetPointsPerEdge(ppe);
366  }
367  else if (m_config["plane"].as<string>().compare("NotSet") != 0)
368  {
369  vector<NekDouble> values;
371  m_config["plane"].as<string>(), values),
372  "Failed to interpret plane string");
373 
374  ASSERTL0(values.size() > 9,
375  "plane string should contain 4 Dim+2 values "
376  "N1,N2,x0,y0,z0,x1,y1,z1,x2,y2,z2,x3,y3,z3");
377 
378  double tmp;
379  ASSERTL0(std::modf(values[0], &tmp) == 0.0, "N1 is not an integer");
380  ASSERTL0(std::modf(values[1], &tmp) == 0.0, "N2 is not an integer");
381 
382  ASSERTL0(values[0] > 1, "N1 is not a valid number");
383  ASSERTL0(values[1] > 1, "N2 is not a valid number");
384 
385  int dim = (values.size() - 2) / 4;
386 
387  Array<OneD, int> npts(2);
388  npts[0] = values[0];
389  npts[1] = values[1];
390 
391  int totpts = npts[0] * npts[1];
392 
393  // Information for partitioning
394  int ptsPerProc = totpts / nprocs;
395  int extraPts = (rank < nprocs - 1) ? 0: totpts % nprocs;
396  int locPts = ptsPerProc + extraPts;
397  int start = rank * ptsPerProc;
398  int end = start + locPts;
399 
401  Array<OneD, NekDouble> delta1(dim);
402  Array<OneD, NekDouble> delta2(dim);
403  for (int i = 0; i < dim; ++i)
404  {
405  pts[i] = Array<OneD, NekDouble>(locPts);
406  delta1[i] = (values[2+1*dim + i] - values[2+0*dim + i])/(npts[0]-1);
407  delta2[i] = (values[2+2*dim + i] - values[2+3*dim + i])/(npts[0]-1);
408  }
409 
410  for (int j = 0, cnt = 0, cntLoc = 0; j < npts[1]; ++j)
411  {
412  for (int i = 0; i < npts[0]; ++i, ++cnt)
413  {
414  if (cnt >= start && cnt < end)
415  {
416  for (int n = 0; n < dim; ++n)
417  {
418  pts[n][cntLoc] =
419  (values[2+n] + i * delta1[n]) *
420  (1.0 - j / ((NekDouble)(npts[1]-1))) +
421  (values[2 + 3*dim + n] + i * delta2[n]) *
422  ( j / ((NekDouble)(npts[1]-1)));
423  }
424  ++cntLoc;
425  }
426  }
427  }
428 
429  vector<size_t> ppe;
430  ppe.push_back(npts[0]);
431  ppe.push_back(npts[1]);
432  m_f->m_fieldPts =
434  pts);
435  m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsPlane);
436  m_f->m_fieldPts->SetPointsPerEdge(ppe);
437  }
438  else if (m_config["box"].as<string>().compare("NotSet") != 0)
439  {
440  vector<NekDouble> values;
442  m_config["box"].as<string>(), values),
443  "Failed to interpret box string");
444 
445  ASSERTL0(values.size() == 9,
446  "box string should contain 9 values "
447  "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
448 
449  int dim = 3;
450  Array<OneD, int> npts(3);
451  npts[0] = values[0];
452  npts[1] = values[1];
453  npts[2] = values[2];
454 
455  int totpts = npts[0]*npts[1]*npts[2];
456 
458  Array<OneD, NekDouble> delta(dim);
459 
460  // Information for partitioning
461  int ptsPerProc = totpts / nprocs;
462  int extraPts = (rank < nprocs - 1) ? 0: totpts % nprocs;
463  int locPts = ptsPerProc + extraPts;
464  int start = rank * ptsPerProc;
465  int end = start + locPts;
466 
467  for (int i = 0; i < dim; ++i)
468  {
469  pts[i] = Array<OneD, NekDouble>(locPts);
470  delta[i] = (values[4 + 2*i] - values[3 + 2*i]) / (npts[i] - 1);
471  }
472 
473 
474 
475  for (int k = 0, cnt = 0, cntLoc = 0; k < npts[2]; ++k)
476  {
477  for (int j = 0; j < npts[1]; ++j)
478  {
479  for (int i = 0; i < npts[0]; ++i, ++cnt)
480  {
481  if (cnt >= start && cnt < end)
482  {
483  pts[0][cntLoc] = values[3] + i * delta[0];
484  pts[1][cntLoc] = values[5] + j * delta[1];
485  pts[2][cntLoc] = values[7] + k * delta[2];
486  ++cntLoc;
487  }
488  }
489  }
490  }
491 
492  vector<size_t> ppe;
493  ppe.push_back(npts[0]);
494  ppe.push_back(npts[1]);
495  ppe.push_back(npts[2]);
496  m_f->m_fieldPts =
498  pts);
499  m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsBox);
500  m_f->m_fieldPts->SetPointsPerEdge(ppe);
501  vector<NekDouble> boxdim;
502  boxdim.assign(&values[3], &values[3] + 6);
503  m_f->m_fieldPts->SetBoxSize(boxdim);
504  }
505  else
506  {
507  ASSERTL0(false,
508  "Missing target points for ProcessInterpPoints.");
509  }
510 }
511 
513  vector<MultiRegions::ExpListSharedPtr> &field0,
515  NekDouble clamp_low,
516  NekDouble clamp_up,
517  NekDouble def_value)
518 {
519  boost::ignore_unused(def_value);
520 
521  ASSERTL0(pts->GetNFields() == field0.size(), "ptField has too few fields");
522 
523  int nfields = field0.size();
524 
525  Interpolator interp;
526  if (m_f->m_comm->GetRank() == 0)
527  {
529  this);
530  }
531  interp.Interpolate(field0, pts);
532  if (m_f->m_comm->GetRank() == 0)
533  {
534  cout << endl;
535  }
536 
537  for (int f = 0; f < nfields; ++f)
538  {
539  for (int i = 0; i < pts->GetNpoints(); ++i)
540  {
541  if (pts->GetPointVal(f, i) > clamp_up)
542  {
543  pts->SetPointVal(f, i, clamp_up);
544  }
545  else if (pts->GetPointVal(f, i) < clamp_low)
546  {
547  pts->SetPointVal(f, i, clamp_low);
548  }
549  }
550  }
551 }
552 
554 {
555  LibUtilities::PtsFieldSharedPtr pts = m_f->m_fieldPts;
556  int dim = pts->GetDim();
557  int nq1 = pts->GetNpoints();
558  int r, f;
559  int pfield = -1;
560  NekDouble p0,qinv;
561  vector<int> velid;
562 
563  vector<NekDouble> values;
564  ASSERTL0(ParseUtils::GenerateVector(m_config["cp"].as<string>(), values),
565  "Failed to interpret cp string");
566 
567  ASSERTL0(values.size() == 2,
568  "cp string should contain 2 values "
569  "p0 and q (=1/2 rho u^2)");
570 
571  p0 = values[0];
572  qinv = 1.0/values[1];
573 
574  for(int i = 0; i < pts->GetNFields(); ++i)
575  {
576  if(boost::iequals(pts->GetFieldName(i),"p"))
577  {
578  pfield = i;
579  }
580 
581  if(boost::iequals(pts->GetFieldName(i),"u")||
582  boost::iequals(pts->GetFieldName(i),"v")||
583  boost::iequals(pts->GetFieldName(i),"w"))
584  {
585  velid.push_back(i);
586  }
587  }
588 
589  if(pfield != -1)
590  {
591  if(!velid.size())
592  {
593  WARNINGL0(false,"Did not find velocity components for Cp0");
594  }
595  }
596  else
597  {
598  WARNINGL0(false,"Failed to find 'p' field to determine cp0");
599  }
600 
601  // Allocate data storage
603 
604  for (f = 0; f < 2; ++f)
605  {
606  data[f] = Array< OneD, NekDouble>(nq1, 0.0);
607  }
608 
609  for (r = 0; r < nq1; r++)
610  {
611  if(pfield != -1) // calculate cp
612  {
613  data[0][r] = qinv*(pts->GetPointVal(dim + pfield, r) - p0);
614 
615  if(velid.size()) // calculate cp0
616  {
617  NekDouble q = 0;
618  for(int i = 0; i < velid.size(); ++i)
619  {
620  q += 0.5*pts->GetPointVal(dim + velid[i], r)*
621  pts->GetPointVal(dim + velid[i], r);
622  }
623  data[1][r] = qinv*(pts->GetPointVal(dim + pfield, r)+q - p0);
624  }
625  }
626  }
627 
628  if(pfield != -1)
629  {
630  pts->AddField(data[0], "Cp");
631  m_f->m_variables.push_back("Cp");
632  if(velid.size())
633  {
634  pts->AddField(data[1], "Cp0");
635  m_f->m_variables.push_back("Cp0");
636  }
637  }
638 }
639 
641  const int goal) const
642 {
643  LibUtilities::PrintProgressbar(position, goal, "Interpolating");
644 }
645 }
646 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#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:223
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:230
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:233
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:265
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:200
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:135
static MeshGraphSharedPtr Read(const LibUtilities::SessionReaderSharedPtr pSession, LibUtilities::DomainRangeShPtr rng=LibUtilities::NullDomainRangeShPtr, bool fillGraph=true)
Definition: MeshGraph.cpp:113
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:989
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:290
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:182
std::shared_ptr< DomainRange > DomainRangeShPtr
Definition: DomainRange.h:61
std::shared_ptr< CsvIO > CsvIOSharedPtr
Definition: CsvIO.h:87
std::shared_ptr< PtsIO > PtsIOSharedPtr
Definition: PtsIO.h:101
CommFactory & GetCommFactory()
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode
Definition: BasisType.h:61
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode
Definition: BasisType.h:60
@ eFourier
Fourier Expansion .
Definition: BasisType.h:53
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:992
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:225
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436
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:892
Represents a command-line configuration option.
Definition: Module.h:134