Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 // 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: Interpolate field to a series of specified points.
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 
36 #include <iostream>
37 #include <string>
38 using namespace std;
39 
40 #include "ProcessInterpPoints.h"
41 
46 #include <boost/lexical_cast.hpp>
47 #include <boost/math/special_functions/fpclassify.hpp>
48 
49 namespace Nektar
50 {
51 namespace FieldUtils
52 {
53 
54 ModuleKey ProcessInterpPoints::className =
56  ModuleKey(eProcessModule, "interppoints"),
57  ProcessInterpPoints::create,
58  "Interpolates a set of points to another, requires fromfld and "
59  "fromxml to be defined, a line, plane or block of points can be "
60  "defined");
61 
62 ProcessInterpPoints::ProcessInterpPoints(FieldSharedPtr f) : ProcessModule(f)
63 {
64 
65  m_config["fromxml"] = ConfigOption(
66  false, "NotSet", "Xml file from which to interpolate field");
67 
68  ASSERTL0(m_config["fromxml"].as<string>().compare("NotSet") != 0,
69  "Need to specify fromxml=file.xml");
70 
71  m_config["fromfld"] = ConfigOption(
72  false, "NotSet", "Fld file from which to interpolate field");
73 
74  ASSERTL0(m_config["fromfld"].as<string>().compare("NotSet") != 0,
75  "Need to specify fromfld=file.fld ");
76 
77  m_config["clamptolowervalue"] =
78  ConfigOption(false, "-10000000", "Lower bound for interpolation value");
79  m_config["clamptouppervalue"] =
80  ConfigOption(false, "10000000", "Upper bound for interpolation value");
81  m_config["defaultvalue"] =
82  ConfigOption(false, "0", "Default value if point is outside domain");
83  m_config["line"] =
84  ConfigOption(false, "NotSet", "Specify a line of N points using "
85  "line=N,x0,y0,z0,z1,y1,z1");
86  m_config["plane"] = ConfigOption(
87  false, "NotSet", "Specify a plane of N1 x N2 points using "
88  "plane=N1,N2,x0,y0,z0,z1,y1,z1,x2,y2,z2,x3,"
89  "y3,z3");
90  m_config["box"] = ConfigOption(
91  false, "NotSet", "Specify a rectangular box of N1 x N2 x N3 points "
92  "using a box of points limited by box="
93  "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
94 
95  m_config["cp"] =
96  ConfigOption(false, "NotSet",
97  "Parameters p0 and q to determine pressure coefficients "
98  "(box only currently)");
99 }
100 
102 {
103 }
104 
105 void ProcessInterpPoints::Process(po::variables_map &vm)
106 {
107  if (m_f->m_verbose)
108  {
109  if (m_f->m_comm->TreatAsRankZero())
110  {
111  cout << "ProcessInterpPoints: interpolating to points..." << endl;
112  }
113  }
114 
115  int rank = m_f->m_comm->GetRank();
116  int nprocs = m_f->m_comm->GetSize();
117 
118  // Check for command line point specification if no .pts file
119  // specified
120  if (m_f->m_fieldPts == LibUtilities::NullPtsField)
121  {
122  if (m_config["line"].as<string>().compare("NotSet") != 0)
123  {
124  string help = m_config["line"].as<string>();
125  vector<NekDouble> values;
127  m_config["line"].as<string>().c_str(), values),
128  "Failed to interpret line string");
129 
130  ASSERTL0(values.size() > 2,
131  "line string should contain 2 Dim+1 values "
132  "N,x0,y0,z0,x1,y1,z1");
133 
134  double tmp;
135  ASSERTL0(std::modf(values[0], &tmp) == 0.0, "N is not an integer");
136  ASSERTL0(values[0] > 1, "N is not a valid number");
137 
138  int dim = (values.size() - 1) / 2;
139  int npts = values[0];
141 
142  for (int i = 0; i < dim; ++i)
143  {
144  pts[i] = Array<OneD, NekDouble>(npts);
145  }
146 
147  for (int i = 0; i < npts; ++i)
148  {
149  pts[0][i] =
150  values[1] +
151  i / ((NekDouble)(npts - 1)) * (values[dim + 1] - values[1]);
152  if (dim > 1)
153  {
154  pts[1][i] = values[2] +
155  i / ((NekDouble)(npts - 1)) *
156  (values[dim + 2] - values[2]);
157 
158  if (dim > 2)
159  {
160  pts[2][i] = values[3] +
161  i / ((NekDouble)(npts - 1)) *
162  (values[dim + 3] - values[3]);
163  }
164  }
165  }
166 
167  vector<int> ppe;
168  ppe.push_back(npts);
169  m_f->m_fieldPts =
171  pts);
172  m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsLine);
173  m_f->m_fieldPts->SetPointsPerEdge(ppe);
174  }
175  else if (m_config["plane"].as<string>().compare("NotSet") != 0)
176  {
177  string help = m_config["plane"].as<string>();
178  vector<NekDouble> values;
180  m_config["plane"].as<string>().c_str(), values),
181  "Failed to interpret plane string");
182 
183  ASSERTL0(values.size() > 9,
184  "plane string should contain 4 Dim+2 values "
185  "N1,N2,x0,y0,z0,x1,y1,z1,x2,y2,z2,x3,y3,z3");
186 
187  double tmp;
188  ASSERTL0(std::modf(values[0], &tmp) == 0.0, "N1 is not an integer");
189  ASSERTL0(std::modf(values[1], &tmp) == 0.0, "N2 is not an integer");
190 
191  ASSERTL0(values[0] > 1, "N1 is not a valid number");
192  ASSERTL0(values[1] > 1, "N2 is not a valid number");
193 
194  int dim = (values.size() - 2) / 4;
195 
196  int npts1 = values[0];
197  int npts2 = values[1];
198 
200 
201  int totpts = npts1 * npts2;
202  int nlocpts = totpts / nprocs;
203 
204  if (rank < nprocs - 1)
205  {
206  for (int i = 0; i < dim; ++i)
207  {
208  pts[i] = Array<OneD, NekDouble>(nlocpts);
209  }
210 
211  int cnt = 0;
212  int cntloc = 0;
213 
214  for (int j = 0; j < npts2; ++j)
215  {
216  for (int i = 0; i < npts1; ++i)
217  {
218 
219  if ((cnt >= rank * nlocpts) &&
220  (cnt < (rank + 1) * nlocpts))
221  {
222  pts[0][cntloc] =
223  (values[2] +
224  i / ((NekDouble)(npts1 - 1)) *
225  (values[dim + 2] - values[2])) *
226  (1.0 - j / ((NekDouble)(npts2 - 1))) +
227  (values[3 * dim + 2] +
228  i / ((NekDouble)(npts1 - 1)) *
229  (values[2 * dim + 2] -
230  values[3 * dim + 2])) *
231  (j / ((NekDouble)(npts2 - 1)));
232 
233  pts[1][cntloc] =
234  (values[3] +
235  i / ((NekDouble)(npts1 - 1)) *
236  (values[dim + 3] - values[3])) *
237  (1.0 - j / ((NekDouble)(npts2 - 1))) +
238  (values[3 * dim + 3] +
239  i / ((NekDouble)(npts1 - 1)) *
240  (values[2 * dim + 3] -
241  values[3 * dim + 3])) *
242  (j / ((NekDouble)(npts2 - 1)));
243 
244  if (dim > 2)
245  {
246  pts[2][cntloc] =
247  (values[4] +
248  i / ((NekDouble)(npts1 - 1)) *
249  (values[dim + 4] - values[4])) *
250  (1.0 - j / ((NekDouble)(npts2 - 1))) +
251  (values[3 * dim + 4] +
252  i / ((NekDouble)(npts1 - 1)) *
253  (values[2 * dim + 4] -
254  values[3 * dim + 4])) *
255  (j / ((NekDouble)(npts2 - 1)));
256  }
257  cntloc++;
258  }
259  cnt++;
260  }
261  }
262  }
263  else
264  {
265  totpts = totpts - rank * nlocpts;
266 
267  for (int i = 0; i < dim; ++i)
268  {
269  pts[i] = Array<OneD, NekDouble>(totpts);
270  }
271 
272  int cnt = 0;
273  int cntloc = 0;
274 
275  for (int j = 0; j < npts2; ++j)
276  {
277  for (int i = 0; i < npts1; ++i)
278  {
279 
280  if (cnt >= rank * nlocpts)
281  {
282  pts[0][cntloc] =
283  (values[2] +
284  i / ((NekDouble)(npts1 - 1)) *
285  (values[dim + 2] - values[2])) *
286  (1.0 - j / ((NekDouble)(npts2 - 1))) +
287  (values[3 * dim + 2] +
288  i / ((NekDouble)(npts1 - 1)) *
289  (values[2 * dim + 2] -
290  values[3 * dim + 2])) *
291  (j / ((NekDouble)(npts2 - 1)));
292 
293  pts[1][cntloc] =
294  (values[3] +
295  i / ((NekDouble)(npts1 - 1)) *
296  (values[dim + 3] - values[3])) *
297  (1.0 - j / ((NekDouble)(npts2 - 1))) +
298  (values[3 * dim + 3] +
299  i / ((NekDouble)(npts1 - 1)) *
300  (values[2 * dim + 3] -
301  values[3 * dim + 3])) *
302  (j / ((NekDouble)(npts2 - 1)));
303 
304  if (dim > 2)
305  {
306  pts[2][cntloc] =
307  (values[4] +
308  i / ((NekDouble)(npts1 - 1)) *
309  (values[dim + 4] - values[4])) *
310  (1.0 - j / ((NekDouble)(npts2 - 1))) +
311  (values[3 * dim + 4] +
312  i / ((NekDouble)(npts1 - 1)) *
313  (values[2 * dim + 4] -
314  values[3 * dim + 4])) *
315  (j / ((NekDouble)(npts2 - 1)));
316  }
317  cntloc++;
318  }
319  cnt++;
320  }
321  }
322  }
323 
324  vector<int> ppe;
325  ppe.push_back(npts1);
326  ppe.push_back(npts2);
327  m_f->m_fieldPts =
329  pts);
330  m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsPlane);
331  m_f->m_fieldPts->SetPointsPerEdge(ppe);
332  }
333  else if (m_config["box"].as<string>().compare("NotSet") != 0)
334  {
335  string help = m_config["box"].as<string>();
336  vector<NekDouble> values;
338  m_config["box"].as<string>().c_str(), values),
339  "Failed to interpret box string");
340 
341  ASSERTL0(values.size() == 9,
342  "box string should contain 9 values "
343  "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
344 
345  int dim = 3;
346 
347  int npts1 = values[0];
348  int npts2 = values[1];
349  int npts3 = values[2];
350 
352 
353  int totpts = npts1 * npts2 * npts3;
354  int nlocpts = totpts / nprocs;
355 
356  if (rank < nprocs - 1) // for rank 0 to nproc-1
357  {
358  totpts = nlocpts;
359 
360  for (int i = 0; i < dim; ++i)
361  {
362  pts[i] = Array<OneD, NekDouble>(totpts);
363  }
364 
365  int cnt = 0;
366  int cntloc = 0;
367 
368  for (int k = 0; k < npts3; ++k)
369  {
370  for (int j = 0; j < npts2; ++j)
371  {
372  for (int i = 0; i < npts1; ++i)
373  {
374  if ((cnt >= rank * nlocpts) &&
375  (cnt < (rank + 1) * nlocpts))
376  {
377  pts[0][cntloc] = values[3] +
378  i / ((NekDouble)(npts1 - 1)) *
379  (values[4] - values[3]);
380  pts[1][cntloc] = values[5] +
381  j / ((NekDouble)(npts2 - 1)) *
382  (values[6] - values[5]);
383  pts[2][cntloc] = values[7] +
384  k / ((NekDouble)(npts3 - 1)) *
385  (values[8] - values[7]);
386  cntloc++;
387  }
388  cnt++;
389  }
390  }
391  }
392  }
393  else // give last rank all remaining points
394  {
395  totpts = totpts - rank * nlocpts;
396 
397  for (int i = 0; i < dim; ++i)
398  {
399  pts[i] = Array<OneD, NekDouble>(totpts);
400  }
401 
402  int cnt = 0;
403  int cntloc = 0;
404 
405  for (int k = 0; k < npts3; ++k)
406  {
407  for (int j = 0; j < npts2; ++j)
408  {
409  for (int i = 0; i < npts1; ++i)
410  {
411  if (cnt >= rank * nlocpts)
412  {
413  pts[0][cntloc] = values[3] +
414  i / ((NekDouble)(npts1 - 1)) *
415  (values[4] - values[3]);
416  pts[1][cntloc] = values[5] +
417  j / ((NekDouble)(npts2 - 1)) *
418  (values[6] - values[5]);
419  pts[2][cntloc] = values[7] +
420  k / ((NekDouble)(npts3 - 1)) *
421  (values[8] - values[7]);
422  cntloc++;
423  }
424  cnt++;
425  }
426  }
427  }
428  }
429 
430  vector<int> ppe;
431  ppe.push_back(npts1);
432  ppe.push_back(npts2);
433  ppe.push_back(npts3);
434  m_f->m_fieldPts =
436  pts);
437  m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsBox);
438  m_f->m_fieldPts->SetPointsPerEdge(ppe);
439  vector<NekDouble> boxdim;
440  boxdim.assign(&values[3], &values[3] + 6);
441  m_f->m_fieldPts->SetBoxSize(boxdim);
442  }
443  }
444 
445  FieldSharedPtr fromField = boost::shared_ptr<Field>(new Field());
446 
447  std::vector<std::string> files;
448  ParseUtils::GenerateOrderedStringVector(m_config["fromxml"].as<string>().c_str(), files);
449  // set up session file for from field
450  fromField->m_session =
452 
453  // Set up range based on min and max of local parallel partition
456 
457  int coordim = m_f->m_fieldPts->GetDim();
458  int npts = m_f->m_fieldPts->GetNpoints();
460  m_f->m_fieldPts->GetPts(pts);
461 
462  rng->m_checkShape = false;
463  rng->m_zmin = -1;
464  rng->m_zmax = 1;
465  rng->m_ymin = -1;
466  rng->m_ymax = 1;
467  switch (coordim)
468  {
469  case 3:
470  rng->m_doZrange = true;
471  rng->m_zmin = Vmath::Vmin(npts, pts[2], 1);
472  rng->m_zmax = Vmath::Vmax(npts, pts[2], 1);
473  if (rng->m_zmax == rng->m_zmin)
474  {
475  rng->m_zmin -= 1;
476  rng->m_zmax += 1;
477  }
478  case 2:
479  rng->m_doYrange = true;
480  rng->m_ymin = Vmath::Vmin(npts, pts[1], 1);
481  rng->m_ymax = Vmath::Vmax(npts, pts[1], 1);
482  case 1:
483  rng->m_doXrange = true;
484  rng->m_xmin = Vmath::Vmin(npts, pts[0], 1);
485  rng->m_xmax = Vmath::Vmax(npts, pts[0], 1);
486  break;
487  default:
488  ASSERTL0(false, "too many values specfied in range");
489  }
490 
491  // setup rng parameters.
492  fromField->m_graph =
493  SpatialDomains::MeshGraph::Read(fromField->m_session, rng);
494 
495  // Read in local from field partitions
496  const SpatialDomains::ExpansionMap &expansions =
497  fromField->m_graph->GetExpansions();
498 
499  Array<OneD, int> ElementGIDs(expansions.size());
500  SpatialDomains::ExpansionMap::const_iterator expIt;
501 
502  int i = 0;
503  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
504  {
505  ElementGIDs[i++] = expIt->second->m_geomShPtr->GetGlobalID();
506  }
507 
508  // check to see that we do have some elmement in teh domain since
509  // possibly all points could be outside of the domain
510  ASSERTL0(i > 0, "No elements are set. Are the interpolated points "
511  "wihtin the domain given by the xml files?");
512 
513  string fromfld = m_config["fromfld"].as<string>();
514  m_f->FieldIOForFile(fromfld)->Import(
515  fromfld, fromField->m_fielddef, fromField->m_data,
517 
518  int NumHomogeneousDir = fromField->m_fielddef[0]->m_numHomogeneousDir;
519 
520  //----------------------------------------------
521  // Set up Expansion information to use mode order from field
522  fromField->m_graph->SetExpansions(fromField->m_fielddef);
523 
524  int nfields = fromField->m_fielddef[0]->m_fields.size();
525 
526  fromField->m_exp.resize(nfields);
527  fromField->m_exp[0] = fromField->SetUpFirstExpList(NumHomogeneousDir, true);
528 
529  m_f->m_exp.resize(nfields);
530 
531  // declare auxiliary fields.
532  for (i = 1; i < nfields; ++i)
533  {
534  fromField->m_exp[i] = fromField->AppendExpList(NumHomogeneousDir);
535  }
536 
537  // load field into expansion in fromfield.
538  for (int j = 0; j < nfields; ++j)
539  {
540  for (i = 0; i < fromField->m_fielddef.size(); i++)
541  {
542  fromField->m_exp[j]->ExtractDataToCoeffs(
543  fromField->m_fielddef[i], fromField->m_data[i],
544  fromField->m_fielddef[0]->m_fields[j],
545  fromField->m_exp[j]->UpdateCoeffs());
546  }
547  fromField->m_exp[j]->BwdTrans(fromField->m_exp[j]->GetCoeffs(),
548  fromField->m_exp[j]->UpdatePhys());
549 
550  Array<OneD, NekDouble> newPts(m_f->m_fieldPts->GetNpoints());
551  m_f->m_fieldPts->AddField(newPts,
552  fromField->m_fielddef[0]->m_fields[j]);
553  }
554 
555  NekDouble clamp_low = m_config["clamptolowervalue"].as<NekDouble>();
556  NekDouble clamp_up = m_config["clamptouppervalue"].as<NekDouble>();
557  NekDouble def_value = m_config["defaultvalue"].as<NekDouble>();
558 
559  InterpolateFieldToPts(fromField->m_exp, m_f->m_fieldPts, clamp_low,
560  clamp_up, def_value);
561 
562  if (!boost::iequals(m_config["cp"].as<string>(), "NotSet"))
563  {
564  calcCp0();
565  }
566 }
567 
569  vector<MultiRegions::ExpListSharedPtr> &field0,
571  NekDouble clamp_low,
572  NekDouble clamp_up,
573  NekDouble def_value)
574 {
575  ASSERTL0(pts->GetNFields() >= field0.size(), "ptField has too few fields");
576 
577  int nfields = field0.size();
578 
579  Interpolator interp;
580  if (m_f->m_comm->GetRank() == 0)
581  {
583  this);
584  }
585  interp.Interpolate(field0, pts);
586  if (m_f->m_comm->GetRank() == 0)
587  {
588  cout << endl;
589  }
590 
591  for (int f = 0; f < nfields; ++f)
592  {
593  for (int i = 0; i < pts->GetNpoints(); ++i)
594  {
595  if (pts->GetPointVal(f, i) > clamp_up)
596  {
597  pts->SetPointVal(f, i, clamp_up);
598  }
599  else if (pts->GetPointVal(f, i) < clamp_low)
600  {
601  pts->SetPointVal(f, i, clamp_low);
602  }
603  }
604  }
605 }
606 
608 {
609  LibUtilities::PtsFieldSharedPtr pts = m_f->m_fieldPts;
610  int dim = pts->GetDim();
611  int nq1 = pts->GetNpoints();
612  int r, f;
613  int pfield = -1;
614  NekDouble p0,qinv;
615  vector<int> velid;
616 
617  vector<NekDouble> values;
619  m_config["cp"].as<string>().c_str(),values),
620  "Failed to interpret cp string");
621 
622  ASSERTL0(values.size() == 2,
623  "cp string should contain 2 values "
624  "p0 and q (=1/2 rho u^2)");
625 
626  p0 = values[0];
627  qinv = 1.0/values[1];
628 
629  for(int i = 0; i < pts->GetNFields(); ++i)
630  {
631  if(boost::iequals(pts->GetFieldName(i),"p"))
632  {
633  pfield = i;
634  }
635 
636  if(boost::iequals(pts->GetFieldName(i),"u")||
637  boost::iequals(pts->GetFieldName(i),"v")||
638  boost::iequals(pts->GetFieldName(i),"w"))
639  {
640  velid.push_back(i);
641  }
642  }
643 
644  if(pfield != -1)
645  {
646  if(!velid.size())
647  {
648  WARNINGL0(false,"Did not find velocity components for Cp0");
649  }
650  }
651  else
652  {
653  WARNINGL0(false,"Failed to find 'p' field to determine cp0");
654  }
655 
656  // Allocate data storage
658 
659  for (f = 0; f < 2; ++f)
660  {
661  data[f] = Array< OneD, NekDouble>(nq1, 0.0);
662  }
663 
664  for (r = 0; r < nq1; r++)
665  {
666  if(pfield != -1) // calculate cp
667  {
668  data[0][r] = qinv*(pts->GetPointVal(dim + pfield, r) - p0);
669 
670  if(velid.size()) // calculate cp0
671  {
672  NekDouble q = 0;
673  for(int i = 0; i < velid.size(); ++i)
674  {
675  q += 0.5*pts->GetPointVal(dim + velid[i], r)*
676  pts->GetPointVal(dim + velid[i], r);
677  }
678  data[1][r] = qinv*(pts->GetPointVal(dim + pfield, r)+q - p0);
679  }
680  }
681  }
682 
683  if(pfield != -1)
684  {
685  pts->AddField(data[0], "Cp");
686  if(velid.size())
687  {
688  pts->AddField(data[1], "Cp0");
689  }
690  }
691 }
692 
694  const int goal) const
695 {
696  LibUtilities::PrintProgressbar(position, goal, "Interpolating");
697 }
698 }
699 }
map< string, ConfigOption > m_config
List of configuration values.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
static bool GenerateOrderedStringVector(const char *const str, std::vector< std::string > &vec)
Definition: ParseUtils.hpp:143
A class that contains algorithms for interpolation between pts fields, expansions and different meshe...
Definition: Interpolator.h:78
void SetProgressCallback(FuncPointerT func, ObjectPointerT obj)
sets a callback funtion which gets called every time the interpolation progresses ...
Definition: Interpolator.h:159
int PrintProgressbar(const int position, const int goal, const string message, int lastprogress=-1)
Prints a progressbar.
Definition: Progressbar.hpp:69
static boost::shared_ptr< MeshGraph > Read(const LibUtilities::SessionReaderSharedPtr &pSession, DomainRangeShPtr &rng=NullDomainRangeShPtr)
Definition: MeshGraph.cpp:124
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
Represents a command-line configuration option.
FIELD_UTILS_EXPORT void Interpolate(const LibUtilities::PtsFieldSharedPtr ptsInField, LibUtilities::PtsFieldSharedPtr &ptsOutField)
Interpolate from a pts field to a pts field.
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:779
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:871
STL namespace.
pair< ModuleType, string > ModuleKey
boost::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:178
static SessionReaderSharedPtr CreateInstance(int argc, char *argv[])
Creates an instance of the SessionReader class.
virtual void Process(po::variables_map &vm)
Write mesh to output file.
boost::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:740
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:204
static std::string npts
Definition: InputFld.cpp:43
boost::shared_ptr< DomainRange > DomainRangeShPtr
Definition: MeshGraph.h:157
double NekDouble
static PtsFieldSharedPtr NullPtsField
Definition: PtsField.h:179
void PrintProgressbar(const int position, const int goal) const
static bool GenerateUnOrderedVector(const char *const str, std::vector< NekDouble > &vec)
Definition: ParseUtils.hpp:128
Abstract base class for processing modules.
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:55
void InterpolateFieldToPts(vector< MultiRegions::ExpListSharedPtr > &field0, LibUtilities::PtsFieldSharedPtr &pts, NekDouble clamp_low, NekDouble clamp_up, NekDouble def_value)
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:174
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