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