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