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/geometry.hpp>
39#include <boost/lexical_cast.hpp>
40
48
49#include "ProcessInterpPoints.h"
50
51using namespace std;
52
53namespace bg = boost::geometry;
54namespace bgi = boost::geometry::index;
55
56namespace Nektar::FieldUtils
57{
58
62 "Interpolates a field to a set of points. Requires fromfld, fromxml "
63 "to be defined, and a topts, line, plane or block of target points ");
64
66{
67 m_config["fromxml"] = ConfigOption(
68 false, "NotSet", "Xml file from which to interpolate field");
69
70 m_config["fromfld"] = ConfigOption(
71 false, "NotSet", "Fld file from which to interpolate field");
72
73 m_config["topts"] =
74 ConfigOption(false, "NotSet", "Pts file to which interpolate field");
75 m_config["line"] = ConfigOption(false, "NotSet",
76 "Specify a line of N points using "
77 "line=N,x0,y0,z0,z1,y1,z1");
78 m_config["plane"] =
79 ConfigOption(false, "NotSet",
80 "Specify a plane of N1 x N2 points using "
81 "plane=N1,N2,x0,y0,z0,z1,y1,z1,x2,y2,z2,x3,y3,z3");
82 m_config["box"] =
83 ConfigOption(false, "NotSet",
84 "Specify a rectangular box of N1 x N2 x N3 points "
85 "using a box of points limited by box="
86 "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
87
88 m_config["clamptolowervalue"] =
89 ConfigOption(false, "-10000000", "Lower bound for interpolation value");
90 m_config["clamptouppervalue"] =
91 ConfigOption(false, "10000000", "Upper bound for interpolation value");
92 m_config["defaultvalue"] =
93 ConfigOption(false, "0", "Default value if point is outside domain");
94
95 m_config["cp"] =
96 ConfigOption(false, "NotSet",
97 "Parameters p0 and q to determine pressure coefficients");
98 m_config["realmodetoimag"] =
99 ConfigOption(false, "NotSet", "Take fields as sin mode");
100 m_config["distTolerance"] =
101 ConfigOption(false, "1e-5", "The maximum acceptable distance.");
102}
103
105{
106}
107
108void ProcessInterpPoints::v_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::MeshGraphIO::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
222 // load field into expansion in fromfield.
223 set<int> sinmode;
224 if (m_config["realmodetoimag"].as<string>().compare("NotSet"))
225 {
226 ParseUtils::GenerateVariableSet(m_config["realmodetoimag"].as<string>(),
227 m_f->m_variables, sinmode);
228 }
229 for (int j = 0; j < nfields; ++j)
230 {
231 for (i = 0; i < fromField->m_fielddef.size(); i++)
232 {
233 fromField->m_exp[j]->ExtractDataToCoeffs(
234 fromField->m_fielddef[i], fromField->m_data[i],
235 fromField->m_fielddef[0]->m_fields[j],
236 fromField->m_exp[j]->UpdateCoeffs());
237 }
238 if (NumHomogeneousDir == 1)
239 {
240 fromField->m_exp[j]->SetWaveSpace(true);
241 if (sinmode.count(j))
242 {
243 int Ncoeff = fromField->m_exp[j]->GetPlane(2)->GetNcoeffs();
245 Ncoeff, -1., fromField->m_exp[j]->GetPlane(2)->GetCoeffs(),
246 1, fromField->m_exp[j]->GetPlane(3)->UpdateCoeffs(), 1);
247 Vmath::Zero(Ncoeff,
248 fromField->m_exp[j]->GetPlane(2)->UpdateCoeffs(),
249 1);
250 }
251 }
252 fromField->m_exp[j]->BwdTrans(fromField->m_exp[j]->GetCoeffs(),
253 fromField->m_exp[j]->UpdatePhys());
254
255 Array<OneD, NekDouble> newPts(m_f->m_fieldPts->GetNpoints());
256 m_f->m_fieldPts->AddField(newPts,
257 fromField->m_fielddef[0]->m_fields[j]);
258 m_f->m_variables.push_back(fromField->m_fielddef[0]->m_fields[j]);
259 }
260
261 NekDouble clamp_low = m_config["clamptolowervalue"].as<NekDouble>();
262 NekDouble clamp_up = m_config["clamptouppervalue"].as<NekDouble>();
263 NekDouble def_value = m_config["defaultvalue"].as<NekDouble>();
264 NekDouble tolerance = m_config["distTolerance"].as<NekDouble>();
265
266 InterpolateFieldToPts(fromField->m_exp, m_f->m_fieldPts, clamp_low,
267 clamp_up, def_value, tolerance);
268
269 if (!boost::iequals(m_config["cp"].as<string>(), "NotSet"))
270 {
271 calcCp0();
272 }
273}
274
275void ProcessInterpPoints::CreateFieldPts([[maybe_unused]] po::variables_map &vm)
276{
277 int rank = m_f->m_comm->GetSpaceComm()->GetRank();
278 int nprocs = m_f->m_comm->GetSpaceComm()->GetSize();
279 // Check for command line point specification
280
281 if (m_config["topts"].as<string>().compare("NotSet") != 0)
282 {
283 string inFile = m_config["topts"].as<string>();
284
285 if (fs::path(inFile).extension() == ".pts")
286 {
289 m_f->m_comm);
290
291 ptsIO->Import(inFile, m_f->m_fieldPts);
292 }
293 else if (fs::path(inFile).extension() == ".csv")
294 {
297 m_f->m_comm);
298
299 csvIO->Import(inFile, m_f->m_fieldPts);
300 }
301 else
302 {
303 ASSERTL0(false, "unknown topts file type");
304 }
305 }
306 else if (m_config["line"].as<string>().compare("NotSet") != 0)
307 {
308 vector<NekDouble> values;
309 ASSERTL0(
310 ParseUtils::GenerateVector(m_config["line"].as<string>(), values),
311 "Failed to interpret line string");
312
313 ASSERTL0(values.size() > 2, "line string should contain 2*Dim+1 values "
314 "N,x0,y0,z0,x1,y1,z1");
315
316 double tmp;
317 ASSERTL0(std::modf(values[0], &tmp) == 0.0, "N is not an integer");
318 ASSERTL0(values[0] > 1, "N is not a valid number");
319
320 int dim = (values.size() - 1) / 2;
321 int npts = values[0];
322
323 // Information for partitioning
324 int ptsPerProc = npts / nprocs;
325 int extraPts = (rank < nprocs - 1) ? 0 : npts % nprocs;
326 int locPts = ptsPerProc + extraPts;
327 int start = rank * ptsPerProc;
328 int end = start + locPts;
329
331 Array<OneD, NekDouble> delta(dim);
332 for (int i = 0; i < dim; ++i)
333 {
334 pts[i] = Array<OneD, NekDouble>(locPts);
335 delta[i] = (values[dim + i + 1] - values[i + 1]) / (npts - 1);
336 }
337
338 for (int i = 0, cntLoc = 0; i < npts; ++i)
339 {
340 if (i >= start && i < end)
341 {
342 for (int n = 0; n < dim; ++n)
343 {
344 pts[n][cntLoc] = values[n + 1] + i * delta[n];
345 }
346 ++cntLoc;
347 }
348 }
349
350 vector<size_t> ppe;
351 ppe.push_back(npts);
352 m_f->m_fieldPts =
354 m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsLine);
355 m_f->m_fieldPts->SetPointsPerEdge(ppe);
356 }
357 else if (m_config["plane"].as<string>().compare("NotSet") != 0)
358 {
359 vector<NekDouble> values;
360 ASSERTL0(
361 ParseUtils::GenerateVector(m_config["plane"].as<string>(), values),
362 "Failed to interpret plane string");
363
364 ASSERTL0(values.size() > 9,
365 "plane string should contain 4 Dim+2 values "
366 "N1,N2,x0,y0,z0,x1,y1,z1,x2,y2,z2,x3,y3,z3");
367
368 double tmp;
369 ASSERTL0(std::modf(values[0], &tmp) == 0.0, "N1 is not an integer");
370 ASSERTL0(std::modf(values[1], &tmp) == 0.0, "N2 is not an integer");
371
372 ASSERTL0(values[0] > 1, "N1 is not a valid number");
373 ASSERTL0(values[1] > 1, "N2 is not a valid number");
374
375 int dim = (values.size() - 2) / 4;
376
377 Array<OneD, int> npts(2);
378 npts[0] = values[0];
379 npts[1] = values[1];
380
381 int totpts = npts[0] * npts[1];
382
383 // Information for partitioning
384 int ptsPerProc = totpts / nprocs;
385 int extraPts = (rank < nprocs - 1) ? 0 : totpts % nprocs;
386 int locPts = ptsPerProc + extraPts;
387 int start = rank * ptsPerProc;
388 int end = start + locPts;
389
391 Array<OneD, NekDouble> delta1(dim);
392 Array<OneD, NekDouble> delta2(dim);
393 for (int i = 0; i < dim; ++i)
394 {
395 pts[i] = Array<OneD, NekDouble>(locPts);
396 delta1[i] = (values[2 + 1 * dim + i] - values[2 + 0 * dim + i]) /
397 (npts[0] - 1);
398 delta2[i] = (values[2 + 2 * dim + i] - values[2 + 3 * dim + i]) /
399 (npts[0] - 1);
400 }
401
402 for (int j = 0, cnt = 0, cntLoc = 0; j < npts[1]; ++j)
403 {
404 for (int i = 0; i < npts[0]; ++i, ++cnt)
405 {
406 if (cnt >= start && cnt < end)
407 {
408 for (int n = 0; n < dim; ++n)
409 {
410 pts[n][cntLoc] =
411 (values[2 + n] + i * delta1[n]) *
412 (1.0 - j / ((NekDouble)(npts[1] - 1))) +
413 (values[2 + 3 * dim + n] + i * delta2[n]) *
414 (j / ((NekDouble)(npts[1] - 1)));
415 }
416 ++cntLoc;
417 }
418 }
419 }
420
421 vector<size_t> ppe;
422 ppe.push_back(npts[0]);
423 ppe.push_back(npts[1]);
424 m_f->m_fieldPts =
426 m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsPlane);
427 m_f->m_fieldPts->SetPointsPerEdge(ppe);
428 }
429 else if (m_config["box"].as<string>().compare("NotSet") != 0)
430 {
431 vector<NekDouble> values;
432 ASSERTL0(
433 ParseUtils::GenerateVector(m_config["box"].as<string>(), values),
434 "Failed to interpret box string");
435
436 ASSERTL0(values.size() == 9, "box string should contain 9 values "
437 "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
438
439 int dim = 3;
440 Array<OneD, int> npts(3);
441 npts[0] = values[0];
442 npts[1] = values[1];
443 npts[2] = values[2];
444
445 int totpts = npts[0] * npts[1] * npts[2];
446
448 Array<OneD, NekDouble> delta(dim);
449
450 // Information for partitioning
451 int ptsPerProc = totpts / nprocs;
452 int extraPts = (rank < nprocs - 1) ? 0 : totpts % nprocs;
453 int locPts = ptsPerProc + extraPts;
454 int start = rank * ptsPerProc;
455 int end = start + locPts;
456
457 for (int i = 0; i < dim; ++i)
458 {
459 pts[i] = Array<OneD, NekDouble>(locPts);
460 delta[i] = (values[4 + 2 * i] - values[3 + 2 * i]) / (npts[i] - 1);
461 }
462
463 for (int k = 0, cnt = 0, cntLoc = 0; k < npts[2]; ++k)
464 {
465 for (int j = 0; j < npts[1]; ++j)
466 {
467 for (int i = 0; i < npts[0]; ++i, ++cnt)
468 {
469 if (cnt >= start && cnt < end)
470 {
471 pts[0][cntLoc] = values[3] + i * delta[0];
472 pts[1][cntLoc] = values[5] + j * delta[1];
473 pts[2][cntLoc] = values[7] + k * delta[2];
474 ++cntLoc;
475 }
476 }
477 }
478 }
479
480 vector<size_t> ppe;
481 ppe.push_back(npts[0]);
482 ppe.push_back(npts[1]);
483 ppe.push_back(npts[2]);
484 m_f->m_fieldPts =
486 m_f->m_fieldPts->SetPtsType(LibUtilities::ePtsBox);
487 m_f->m_fieldPts->SetPointsPerEdge(ppe);
488 vector<NekDouble> boxdim;
489 boxdim.assign(&values[3], &values[3] + 6);
490 m_f->m_fieldPts->SetBoxSize(boxdim);
491 }
492 else
493 {
494 ASSERTL0(false, "Missing target points for ProcessInterpPoints.");
495 }
496}
497
499 vector<MultiRegions::ExpListSharedPtr> &field0,
501 NekDouble clamp_up, NekDouble def_value, NekDouble tolerance)
502{
503 ASSERTL0(pts->GetNFields() == field0.size(), "ptField has too few fields");
504
505 int nfields = field0.size();
506
508 if (m_f->m_comm->GetRank() == 0)
509 {
511 this);
512 }
513 interp.Interpolate(field0, pts, def_value, tolerance);
514
515 if (m_f->m_comm->GetRank() == 0)
516 {
517 cout << endl;
518 }
519
520 for (int f = 0; f < nfields; ++f)
521 {
522 for (int i = 0; i < pts->GetNpoints(); ++i)
523 {
524 if (pts->GetPointVal(f, i) > clamp_up)
525 {
526 pts->SetPointVal(f, i, clamp_up);
527 }
528 else if (pts->GetPointVal(f, i) < clamp_low)
529 {
530 pts->SetPointVal(f, i, clamp_low);
531 }
532 }
533 }
534}
535
537{
538 LibUtilities::PtsFieldSharedPtr pts = m_f->m_fieldPts;
539 int dim = pts->GetDim();
540 int nq1 = pts->GetNpoints();
541 int r, f;
542 int pfield = -1;
543 NekDouble p0, qinv;
544 vector<int> velid;
545
546 vector<NekDouble> values;
547 ASSERTL0(ParseUtils::GenerateVector(m_config["cp"].as<string>(), values),
548 "Failed to interpret cp string");
549
550 ASSERTL0(values.size() == 2, "cp string should contain 2 values "
551 "p0 and q (=1/2 rho u^2)");
552
553 p0 = values[0];
554 qinv = 1.0 / values[1];
555
556 for (int i = 0; i < pts->GetNFields(); ++i)
557 {
558 if (boost::iequals(pts->GetFieldName(i), "p"))
559 {
560 pfield = i;
561 }
562
563 if (boost::iequals(pts->GetFieldName(i), "u") ||
564 boost::iequals(pts->GetFieldName(i), "v") ||
565 boost::iequals(pts->GetFieldName(i), "w"))
566 {
567 velid.push_back(i);
568 }
569 }
570
571 if (pfield != -1)
572 {
573 if (!velid.size())
574 {
575 WARNINGL0(false, "Did not find velocity components for Cp0");
576 }
577 }
578 else
579 {
580 WARNINGL0(false, "Failed to find 'p' field to determine cp0");
581 }
582
583 // Allocate data storage
585
586 for (f = 0; f < 2; ++f)
587 {
588 data[f] = Array<OneD, NekDouble>(nq1, 0.0);
589 }
590
591 for (r = 0; r < nq1; r++)
592 {
593 if (pfield != -1) // calculate cp
594 {
595 data[0][r] = qinv * (pts->GetPointVal(dim + pfield, r) - p0);
596
597 if (velid.size()) // calculate cp0
598 {
599 NekDouble q = 0;
600 for (int i = 0; i < velid.size(); ++i)
601 {
602 q += 0.5 * pts->GetPointVal(dim + velid[i], r) *
603 pts->GetPointVal(dim + velid[i], r);
604 }
605 data[1][r] =
606 qinv * (pts->GetPointVal(dim + pfield, r) + q - p0);
607 }
608 }
609 }
610
611 if (pfield != -1)
612 {
613 pts->AddField(data[0], "Cp");
614 m_f->m_variables.push_back("Cp");
615 if (velid.size())
616 {
617 pts->AddField(data[1], "Cp0");
618 m_f->m_variables.push_back("Cp0");
619 }
620 }
621}
622
624 const int goal) const
625{
626 LibUtilities::PrintProgressbar(position, goal, "Interpolating");
627}
628} // namespace Nektar::FieldUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:215
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., NekDouble tolerance=NekConstants::kFindDistanceMin)
Interpolate from an expansion to an expansion.
FieldSharedPtr m_f
Field object.
Definition: Module.h:239
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:272
void InterpolateFieldToPts(std::vector< MultiRegions::ExpListSharedPtr > &field0, LibUtilities::PtsFieldSharedPtr &pts, NekDouble clamp_low, NekDouble clamp_up, NekDouble def_value, NekDouble tolerance)
void PrintProgressbar(const int position, const int goal) const
static std::shared_ptr< Module > create(FieldSharedPtr f)
Creates an instance of this class.
void v_Process(po::variables_map &vm) override
Write mesh to output file.
Abstract base class for processing modules.
Definition: Module.h:301
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.
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:168
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:130
static MeshGraphSharedPtr Read(const LibUtilities::SessionReaderSharedPtr pSession, LibUtilities::DomainRangeShPtr rng=LibUtilities::NullDomainRangeShPtr, bool fillGraph=true, SpatialDomains::MeshGraphSharedPtr partitionedGraph=nullptr)
Definition: MeshGraphIO.cpp:53
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:1026
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:180
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:47
int PrintProgressbar(const int position, const int goal, const std::string message, int lastprogress=-1)
Prints a progressbar.
Definition: Progressbar.hpp:65
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:51
std::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:184
std::shared_ptr< DomainRange > DomainRangeShPtr
Definition: DomainRange.h:64
std::shared_ptr< CsvIO > CsvIOSharedPtr
Definition: CsvIO.h:74
std::shared_ptr< PtsIO > PtsIOSharedPtr
Definition: PtsIO.h:90
CommFactory & GetCommFactory()
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode .
Definition: BasisType.h:68
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode .
Definition: BasisType.h:66
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55
std::map< int, ExpansionInfoShPtr > ExpansionInfoMap
Definition: MeshGraph.h:141
std::vector< double > q(NPUPPER *NPUPPER)
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.hpp:725
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.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
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.hpp:644
STL namespace.
Represents a command-line configuration option.
Definition: Module.h:129