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 
38 #include <boost/core/ignore_unused.hpp>
39 #include <boost/geometry.hpp>
40 #include <boost/lexical_cast.hpp>
41 #include <boost/math/special_functions/fpclassify.hpp>
42 
49 
50 #include "ProcessInterpPoints.h"
51 
52 using namespace std;
53 
54 namespace bg = boost::geometry;
55 namespace bgi = boost::geometry::index;
56 
57 namespace Nektar
58 {
59 namespace FieldUtils
60 {
61 
62 ModuleKey ProcessInterpPoints::className =
64  ModuleKey(eProcessModule, "interppoints"), 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"] =
77  ConfigOption(false, "NotSet", "Pts file to which interpolate field");
78  m_config["line"] = ConfigOption(false, "NotSet",
79  "Specify a line of N points using "
80  "line=N,x0,y0,z0,z1,y1,z1");
81  m_config["plane"] =
82  ConfigOption(false, "NotSet",
83  "Specify a plane of N1 x N2 points using "
84  "plane=N1,N2,x0,y0,z0,z1,y1,z1,x2,y2,z2,x3,y3,z3");
85  m_config["box"] =
86  ConfigOption(false, "NotSet",
87  "Specify a rectangular box of N1 x N2 x N3 points "
88  "using a box of points limited by box="
89  "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
90 
91  m_config["clamptolowervalue"] =
92  ConfigOption(false, "-10000000", "Lower bound for interpolation value");
93  m_config["clamptouppervalue"] =
94  ConfigOption(false, "10000000", "Upper bound for interpolation value");
95  m_config["defaultvalue"] =
96  ConfigOption(false, "0", "Default value if point is outside domain");
97 
98  m_config["cp"] =
99  ConfigOption(false, "NotSet",
100  "Parameters p0 and q to determine pressure coefficients");
101  m_config["realmodetoimag"] =
102  ConfigOption(false, "NotSet", "Take fields as sin mode");
103 }
104 
106 {
107 }
108 
109 void ProcessInterpPoints::v_Process(po::variables_map &vm)
110 {
111  m_f->SetUpExp(vm);
112 
113  CreateFieldPts(vm);
114 
115  FieldSharedPtr fromField = std::shared_ptr<Field>(new Field());
116  std::vector<std::string> files;
117  ParseUtils::GenerateVector(m_config["fromxml"].as<string>(), files);
118 
119  // set up session file for from field
120  char *argv[] = {const_cast<char *>("FieldConvert"), nullptr};
121  fromField->m_session = LibUtilities::SessionReader::CreateInstance(
122  1, argv, files,
123  LibUtilities::GetCommFactory().CreateInstance("Serial", 0, 0));
124 
125  // Set up range based on min and max of local parallel partition
128 
129  int coordim = m_f->m_fieldPts->GetDim();
130  int npts = m_f->m_fieldPts->GetNpoints();
131  std::vector<std::string> fieldNames = m_f->m_fieldPts->GetFieldNames();
132  for (auto &it : fieldNames)
133  {
134  m_f->m_fieldPts->RemoveField(it);
135  }
136 
138  m_f->m_fieldPts->GetPts(pts);
139 
140  rng->m_checkShape = false;
141  rng->m_zmin = -1;
142  rng->m_zmax = 1;
143  rng->m_ymin = -1;
144  rng->m_ymax = 1;
145  switch (coordim)
146  {
147  case 3:
148  rng->m_doZrange = true;
149  rng->m_zmin = Vmath::Vmin(npts, pts[2], 1);
150  rng->m_zmax = Vmath::Vmax(npts, pts[2], 1);
151  if (rng->m_zmax == rng->m_zmin)
152  {
153  rng->m_zmin -= 1;
154  rng->m_zmax += 1;
155  }
156  /* Falls through. */
157  case 2:
158  rng->m_doYrange = true;
159  rng->m_ymin = Vmath::Vmin(npts, pts[1], 1);
160  rng->m_ymax = Vmath::Vmax(npts, pts[1], 1);
161  /* Falls through. */
162  case 1:
163  rng->m_doXrange = true;
164  rng->m_xmin = Vmath::Vmin(npts, pts[0], 1);
165  rng->m_xmax = Vmath::Vmax(npts, pts[0], 1);
166  break;
167  default:
168  NEKERROR(ErrorUtil::efatal, "Too many values specified in range");
169  }
170 
171  // setup rng parameters.
172  fromField->m_graph =
173  SpatialDomains::MeshGraph::Read(fromField->m_session, rng);
174 
175  // Read in local from field partitions
176  const SpatialDomains::ExpansionInfoMap &expansions =
177  fromField->m_graph->GetExpansionInfo();
178  Array<OneD, int> ElementGIDs(expansions.size());
179 
180  int i = 0;
181  for (auto &expIt : expansions)
182  {
183  ElementGIDs[i++] = expIt.second->m_geomShPtr->GetGlobalID();
184  }
185  // check to see that we do have some element in the domain since
186  // possibly all points could be outside of the domain
187  ASSERTL0(i > 0, "No elements are set. Are the interpolated points "
188  "within the domain given by the xml files?");
189  string fromfld = m_config["fromfld"].as<string>();
190  m_f->FieldIOForFile(fromfld)->Import(
191  fromfld, fromField->m_fielddef, fromField->m_data,
193  int NumHomogeneousDir = fromField->m_fielddef[0]->m_numHomogeneousDir;
194  for (i = 0; i < fromField->m_fielddef.size(); ++i)
195  {
196  int d1 = fromField->m_fielddef[i]->m_basis.size();
197  d1 -= 1;
198  if (d1 >= 0 && (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 
223  // load field into expansion in fromfield.
224  set<int> sinmode;
225  if (m_config["realmodetoimag"].as<string>().compare("NotSet"))
226  {
227  ParseUtils::GenerateVariableSet(m_config["realmodetoimag"].as<string>(),
228  m_f->m_variables, sinmode);
229  }
230  for (int j = 0; j < nfields; ++j)
231  {
232  for (i = 0; i < fromField->m_fielddef.size(); i++)
233  {
234  fromField->m_exp[j]->ExtractDataToCoeffs(
235  fromField->m_fielddef[i], fromField->m_data[i],
236  fromField->m_fielddef[0]->m_fields[j],
237  fromField->m_exp[j]->UpdateCoeffs());
238  }
239  if (NumHomogeneousDir == 1)
240  {
241  fromField->m_exp[j]->SetWaveSpace(true);
242  if (sinmode.count(j))
243  {
244  int Ncoeff = fromField->m_exp[j]->GetPlane(2)->GetNcoeffs();
245  Vmath::Smul(
246  Ncoeff, -1., fromField->m_exp[j]->GetPlane(2)->GetCoeffs(),
247  1, fromField->m_exp[j]->GetPlane(3)->UpdateCoeffs(), 1);
248  Vmath::Zero(Ncoeff,
249  fromField->m_exp[j]->GetPlane(2)->UpdateCoeffs(),
250  1);
251  }
252  }
253  fromField->m_exp[j]->BwdTrans(fromField->m_exp[j]->GetCoeffs(),
254  fromField->m_exp[j]->UpdatePhys());
255 
256  Array<OneD, NekDouble> newPts(m_f->m_fieldPts->GetNpoints());
257  m_f->m_fieldPts->AddField(newPts,
258  fromField->m_fielddef[0]->m_fields[j]);
259  m_f->m_variables.push_back(fromField->m_fielddef[0]->m_fields[j]);
260  }
261 
262  NekDouble clamp_low = m_config["clamptolowervalue"].as<NekDouble>();
263  NekDouble clamp_up = m_config["clamptouppervalue"].as<NekDouble>();
264  NekDouble def_value = m_config["defaultvalue"].as<NekDouble>();
265 
266  InterpolateFieldToPts(fromField->m_exp, m_f->m_fieldPts, clamp_low,
267  clamp_up, def_value);
268 
269  if (!boost::iequals(m_config["cp"].as<string>(), "NotSet"))
270  {
271  calcCp0();
272  }
273 }
274 
275 void ProcessInterpPoints::CreateFieldPts(po::variables_map &vm)
276 {
277  boost::ignore_unused(vm);
278 
279  int rank = m_f->m_comm->GetSpaceComm()->GetRank();
280  int nprocs = m_f->m_comm->GetSpaceComm()->GetSize();
281  // Check for command line point specification
282 
283  if (m_config["topts"].as<string>().compare("NotSet") != 0)
284  {
285  string inFile = m_config["topts"].as<string>();
286 
287  if (boost::filesystem::path(inFile).extension() == ".pts")
288  {
291  m_f->m_comm);
292 
293  ptsIO->Import(inFile, m_f->m_fieldPts);
294  }
295  else if (boost::filesystem::path(inFile).extension() == ".csv")
296  {
299  m_f->m_comm);
300 
301  csvIO->Import(inFile, m_f->m_fieldPts);
302  }
303  else
304  {
305  ASSERTL0(false, "unknown topts file type");
306  }
307  }
308  else if (m_config["line"].as<string>().compare("NotSet") != 0)
309  {
310  vector<NekDouble> values;
311  ASSERTL0(
312  ParseUtils::GenerateVector(m_config["line"].as<string>(), values),
313  "Failed to interpret line string");
314 
315  ASSERTL0(values.size() > 2, "line string should contain 2*Dim+1 values "
316  "N,x0,y0,z0,x1,y1,z1");
317 
318  double tmp;
319  ASSERTL0(std::modf(values[0], &tmp) == 0.0, "N is not an integer");
320  ASSERTL0(values[0] > 1, "N is not a valid number");
321 
322  int dim = (values.size() - 1) / 2;
323  int npts = values[0];
324 
325  // Information for partitioning
326  int ptsPerProc = npts / nprocs;
327  int extraPts = (rank < nprocs - 1) ? 0 : npts % nprocs;
328  int locPts = ptsPerProc + extraPts;
329  int start = rank * ptsPerProc;
330  int end = start + locPts;
331 
333  Array<OneD, NekDouble> delta(dim);
334  for (int i = 0; i < dim; ++i)
335  {
336  pts[i] = Array<OneD, NekDouble>(locPts);
337  delta[i] = (values[dim + i + 1] - values[i + 1]) / (npts - 1);
338  }
339 
340  for (int i = 0, cntLoc = 0; i < npts; ++i)
341  {
342  if (i >= start && i < end)
343  {
344  for (int n = 0; n < dim; ++n)
345  {
346  pts[n][cntLoc] = values[n + 1] + i * delta[n];
347  }
348  ++cntLoc;
349  }
350  }
351 
352  vector<size_t> ppe;
353  ppe.push_back(npts);
354  m_f->m_fieldPts =
356  m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsLine);
357  m_f->m_fieldPts->SetPointsPerEdge(ppe);
358  }
359  else if (m_config["plane"].as<string>().compare("NotSet") != 0)
360  {
361  vector<NekDouble> values;
362  ASSERTL0(
363  ParseUtils::GenerateVector(m_config["plane"].as<string>(), values),
364  "Failed to interpret plane string");
365 
366  ASSERTL0(values.size() > 9,
367  "plane string should contain 4 Dim+2 values "
368  "N1,N2,x0,y0,z0,x1,y1,z1,x2,y2,z2,x3,y3,z3");
369 
370  double tmp;
371  ASSERTL0(std::modf(values[0], &tmp) == 0.0, "N1 is not an integer");
372  ASSERTL0(std::modf(values[1], &tmp) == 0.0, "N2 is not an integer");
373 
374  ASSERTL0(values[0] > 1, "N1 is not a valid number");
375  ASSERTL0(values[1] > 1, "N2 is not a valid number");
376 
377  int dim = (values.size() - 2) / 4;
378 
379  Array<OneD, int> npts(2);
380  npts[0] = values[0];
381  npts[1] = values[1];
382 
383  int totpts = npts[0] * npts[1];
384 
385  // Information for partitioning
386  int ptsPerProc = totpts / nprocs;
387  int extraPts = (rank < nprocs - 1) ? 0 : totpts % nprocs;
388  int locPts = ptsPerProc + extraPts;
389  int start = rank * ptsPerProc;
390  int end = start + locPts;
391 
393  Array<OneD, NekDouble> delta1(dim);
394  Array<OneD, NekDouble> delta2(dim);
395  for (int i = 0; i < dim; ++i)
396  {
397  pts[i] = Array<OneD, NekDouble>(locPts);
398  delta1[i] = (values[2 + 1 * dim + i] - values[2 + 0 * dim + i]) /
399  (npts[0] - 1);
400  delta2[i] = (values[2 + 2 * dim + i] - values[2 + 3 * dim + i]) /
401  (npts[0] - 1);
402  }
403 
404  for (int j = 0, cnt = 0, cntLoc = 0; j < npts[1]; ++j)
405  {
406  for (int i = 0; i < npts[0]; ++i, ++cnt)
407  {
408  if (cnt >= start && cnt < end)
409  {
410  for (int n = 0; n < dim; ++n)
411  {
412  pts[n][cntLoc] =
413  (values[2 + n] + i * delta1[n]) *
414  (1.0 - j / ((NekDouble)(npts[1] - 1))) +
415  (values[2 + 3 * dim + n] + i * delta2[n]) *
416  (j / ((NekDouble)(npts[1] - 1)));
417  }
418  ++cntLoc;
419  }
420  }
421  }
422 
423  vector<size_t> ppe;
424  ppe.push_back(npts[0]);
425  ppe.push_back(npts[1]);
426  m_f->m_fieldPts =
428  m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsPlane);
429  m_f->m_fieldPts->SetPointsPerEdge(ppe);
430  }
431  else if (m_config["box"].as<string>().compare("NotSet") != 0)
432  {
433  vector<NekDouble> values;
434  ASSERTL0(
435  ParseUtils::GenerateVector(m_config["box"].as<string>(), values),
436  "Failed to interpret box string");
437 
438  ASSERTL0(values.size() == 9, "box string should contain 9 values "
439  "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
440 
441  int dim = 3;
442  Array<OneD, int> npts(3);
443  npts[0] = values[0];
444  npts[1] = values[1];
445  npts[2] = values[2];
446 
447  int totpts = npts[0] * npts[1] * npts[2];
448 
450  Array<OneD, NekDouble> delta(dim);
451 
452  // Information for partitioning
453  int ptsPerProc = totpts / nprocs;
454  int extraPts = (rank < nprocs - 1) ? 0 : totpts % nprocs;
455  int locPts = ptsPerProc + extraPts;
456  int start = rank * ptsPerProc;
457  int end = start + locPts;
458 
459  for (int i = 0; i < dim; ++i)
460  {
461  pts[i] = Array<OneD, NekDouble>(locPts);
462  delta[i] = (values[4 + 2 * i] - values[3 + 2 * i]) / (npts[i] - 1);
463  }
464 
465  for (int k = 0, cnt = 0, cntLoc = 0; k < npts[2]; ++k)
466  {
467  for (int j = 0; j < npts[1]; ++j)
468  {
469  for (int i = 0; i < npts[0]; ++i, ++cnt)
470  {
471  if (cnt >= start && cnt < end)
472  {
473  pts[0][cntLoc] = values[3] + i * delta[0];
474  pts[1][cntLoc] = values[5] + j * delta[1];
475  pts[2][cntLoc] = values[7] + k * delta[2];
476  ++cntLoc;
477  }
478  }
479  }
480  }
481 
482  vector<size_t> ppe;
483  ppe.push_back(npts[0]);
484  ppe.push_back(npts[1]);
485  ppe.push_back(npts[2]);
486  m_f->m_fieldPts =
488  m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsBox);
489  m_f->m_fieldPts->SetPointsPerEdge(ppe);
490  vector<NekDouble> boxdim;
491  boxdim.assign(&values[3], &values[3] + 6);
492  m_f->m_fieldPts->SetBoxSize(boxdim);
493  }
494  else
495  {
496  ASSERTL0(false, "Missing target points for ProcessInterpPoints.");
497  }
498 }
499 
501  vector<MultiRegions::ExpListSharedPtr> &field0,
503  NekDouble clamp_up, NekDouble def_value)
504 {
505  boost::ignore_unused(def_value);
506 
507  ASSERTL0(pts->GetNFields() == field0.size(), "ptField has too few fields");
508 
509  int nfields = field0.size();
510 
512  if (m_f->m_comm->GetRank() == 0)
513  {
515  this);
516  }
517  interp.Interpolate(field0, pts);
518 
519  if (m_f->m_comm->GetRank() == 0)
520  {
521  cout << endl;
522  }
523 
524  for (int f = 0; f < nfields; ++f)
525  {
526  for (int i = 0; i < pts->GetNpoints(); ++i)
527  {
528  if (pts->GetPointVal(f, i) > clamp_up)
529  {
530  pts->SetPointVal(f, i, clamp_up);
531  }
532  else if (pts->GetPointVal(f, i) < clamp_low)
533  {
534  pts->SetPointVal(f, i, clamp_low);
535  }
536  }
537  }
538 }
539 
541 {
542  LibUtilities::PtsFieldSharedPtr pts = m_f->m_fieldPts;
543  int dim = pts->GetDim();
544  int nq1 = pts->GetNpoints();
545  int r, f;
546  int pfield = -1;
547  NekDouble p0, qinv;
548  vector<int> velid;
549 
550  vector<NekDouble> values;
551  ASSERTL0(ParseUtils::GenerateVector(m_config["cp"].as<string>(), values),
552  "Failed to interpret cp string");
553 
554  ASSERTL0(values.size() == 2, "cp string should contain 2 values "
555  "p0 and q (=1/2 rho u^2)");
556 
557  p0 = values[0];
558  qinv = 1.0 / values[1];
559 
560  for (int i = 0; i < pts->GetNFields(); ++i)
561  {
562  if (boost::iequals(pts->GetFieldName(i), "p"))
563  {
564  pfield = i;
565  }
566 
567  if (boost::iequals(pts->GetFieldName(i), "u") ||
568  boost::iequals(pts->GetFieldName(i), "v") ||
569  boost::iequals(pts->GetFieldName(i), "w"))
570  {
571  velid.push_back(i);
572  }
573  }
574 
575  if (pfield != -1)
576  {
577  if (!velid.size())
578  {
579  WARNINGL0(false, "Did not find velocity components for Cp0");
580  }
581  }
582  else
583  {
584  WARNINGL0(false, "Failed to find 'p' field to determine cp0");
585  }
586 
587  // Allocate data storage
589 
590  for (f = 0; f < 2; ++f)
591  {
592  data[f] = Array<OneD, NekDouble>(nq1, 0.0);
593  }
594 
595  for (r = 0; r < nq1; r++)
596  {
597  if (pfield != -1) // calculate cp
598  {
599  data[0][r] = qinv * (pts->GetPointVal(dim + pfield, r) - p0);
600 
601  if (velid.size()) // calculate cp0
602  {
603  NekDouble q = 0;
604  for (int i = 0; i < velid.size(); ++i)
605  {
606  q += 0.5 * pts->GetPointVal(dim + velid[i], r) *
607  pts->GetPointVal(dim + velid[i], r);
608  }
609  data[1][r] =
610  qinv * (pts->GetPointVal(dim + pfield, r) + q - p0);
611  }
612  }
613  }
614 
615  if (pfield != -1)
616  {
617  pts->AddField(data[0], "Cp");
618  m_f->m_variables.push_back("Cp");
619  if (velid.size())
620  {
621  pts->AddField(data[1], "Cp0");
622  m_f->m_variables.push_back("Cp0");
623  }
624  }
625 }
626 
628  const int goal) const
629 {
630  LibUtilities::PrintProgressbar(position, goal, "Interpolating");
631 }
632 } // namespace FieldUtils
633 } // 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 T expInField, T &expOutField, NekDouble def_value=0.0)
Interpolate from an expansion to an expansion.
FieldSharedPtr m_f
Field object.
Definition: Module.h:234
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:263
void InterpolateFieldToPts(std::vector< MultiRegions::ExpListSharedPtr > &field0, LibUtilities::PtsFieldSharedPtr &pts, NekDouble clamp_low, NekDouble clamp_up, NekDouble def_value)
void PrintProgressbar(const int position, const int goal) const
virtual void v_Process(po::variables_map &vm) override
Write mesh to output file.
Abstract base class for processing modules.
Definition: Module.h:292
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 GenerateVariableSet(const std::string &str, const std::vector< std::string > &variables, std::set< int > &out)
Generate a set of variable locations.
Definition: ParseUtils.cpp:166
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:131
static MeshGraphSharedPtr Read(const LibUtilities::SessionReaderSharedPtr pSession, LibUtilities::DomainRangeShPtr rng=LibUtilities::NullDomainRangeShPtr, bool fillGraph=true)
Definition: MeshGraph.cpp:111
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:991
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:317
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:80
std::shared_ptr< PtsIO > PtsIOSharedPtr
Definition: PtsIO.h:96
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:2
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