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>
54namespace bg = boost::geometry;
55namespace bgi = boost::geometry::index;
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 ");
71 false,
"NotSet",
"Xml file from which to interpolate field");
74 false,
"NotSet",
"Fld file from which to interpolate field");
77 ConfigOption(
false,
"NotSet",
"Pts file to which interpolate field");
79 "Specify a line of N points using "
80 "line=N,x0,y0,z0,z1,y1,z1");
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");
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");
92 ConfigOption(
false,
"-10000000",
"Lower bound for interpolation value");
94 ConfigOption(
false,
"10000000",
"Upper bound for interpolation value");
96 ConfigOption(
false,
"0",
"Default value if point is outside domain");
100 "Parameters p0 and q to determine pressure coefficients");
102 ConfigOption(
false,
"NotSet",
"Take fields as sin mode");
116 std::vector<std::string> files;
120 char *argv[] = {
const_cast<char *
>(
"FieldConvert"),
nullptr};
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)
134 m_f->m_fieldPts->RemoveField(it);
138 m_f->m_fieldPts->GetPts(pts);
140 rng->m_checkShape =
false;
148 rng->m_doZrange =
true;
151 if (rng->m_zmax == rng->m_zmin)
158 rng->m_doYrange =
true;
163 rng->m_doXrange =
true;
177 fromField->m_graph->GetExpansionInfo();
181 for (
auto &expIt : expansions)
183 ElementGIDs[i++] = expIt.second->m_geomShPtr->GetGlobalID();
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)
196 int d1 = fromField->m_fielddef[i]->m_basis.size();
198 if (d1 >= 0 && (fromField->m_fielddef[i]->m_basis[d1] ==
200 fromField->m_fielddef[i]->m_basis[d1] ==
203 fromField->m_fielddef[i]->m_homogeneousZIDs[0] += 2;
204 fromField->m_fielddef[i]->m_numModes[d1] = 4;
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);
218 for (i = 1; i < nfields; ++i)
220 fromField->m_exp[i] = fromField->AppendExpList(NumHomogeneousDir);
225 if (
m_config[
"realmodetoimag"].as<string>().compare(
"NotSet"))
228 m_f->m_variables, sinmode);
230 for (
int j = 0; j < nfields; ++j)
232 for (i = 0; i < fromField->m_fielddef.size(); i++)
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());
239 if (NumHomogeneousDir == 1)
241 fromField->m_exp[j]->SetWaveSpace(
true);
242 if (sinmode.count(j))
244 int Ncoeff = fromField->m_exp[j]->GetPlane(2)->GetNcoeffs();
246 Ncoeff, -1., fromField->m_exp[j]->GetPlane(2)->GetCoeffs(),
247 1, fromField->m_exp[j]->GetPlane(3)->UpdateCoeffs(), 1);
249 fromField->m_exp[j]->GetPlane(2)->UpdateCoeffs(),
253 fromField->m_exp[j]->BwdTrans(fromField->m_exp[j]->GetCoeffs(),
254 fromField->m_exp[j]->UpdatePhys());
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]);
267 clamp_up, def_value);
269 if (!boost::iequals(
m_config[
"cp"].as<string>(),
"NotSet"))
277 boost::ignore_unused(vm);
279 int rank =
m_f->m_comm->GetSpaceComm()->GetRank();
280 int nprocs =
m_f->m_comm->GetSpaceComm()->GetSize();
283 if (
m_config[
"topts"].as<string>().compare(
"NotSet") != 0)
285 string inFile =
m_config[
"topts"].as<
string>();
287 if (boost::filesystem::path(inFile).extension() ==
".pts")
293 ptsIO->Import(inFile,
m_f->m_fieldPts);
295 else if (boost::filesystem::path(inFile).extension() ==
".csv")
301 csvIO->Import(inFile,
m_f->m_fieldPts);
305 ASSERTL0(
false,
"unknown topts file type");
308 else if (
m_config[
"line"].as<string>().compare(
"NotSet") != 0)
310 vector<NekDouble> values;
313 "Failed to interpret line string");
315 ASSERTL0(values.size() > 2,
"line string should contain 2*Dim+1 values "
316 "N,x0,y0,z0,x1,y1,z1");
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");
322 int dim = (values.size() - 1) / 2;
323 int npts = values[0];
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;
334 for (
int i = 0; i < dim; ++i)
337 delta[i] = (values[dim + i + 1] - values[i + 1]) / (npts - 1);
340 for (
int i = 0, cntLoc = 0; i < npts; ++i)
342 if (i >= start && i < end)
344 for (
int n = 0; n < dim; ++n)
346 pts[n][cntLoc] = values[n + 1] + i * delta[n];
357 m_f->m_fieldPts->SetPointsPerEdge(ppe);
359 else if (
m_config[
"plane"].as<string>().compare(
"NotSet") != 0)
361 vector<NekDouble> values;
364 "Failed to interpret plane string");
367 "plane string should contain 4 Dim+2 values "
368 "N1,N2,x0,y0,z0,x1,y1,z1,x2,y2,z2,x3,y3,z3");
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");
374 ASSERTL0(values[0] > 1,
"N1 is not a valid number");
375 ASSERTL0(values[1] > 1,
"N2 is not a valid number");
377 int dim = (values.size() - 2) / 4;
383 int totpts = npts[0] * npts[1];
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;
395 for (
int i = 0; i < dim; ++i)
398 delta1[i] = (values[2 + 1 * dim + i] - values[2 + 0 * dim + i]) /
400 delta2[i] = (values[2 + 2 * dim + i] - values[2 + 3 * dim + i]) /
404 for (
int j = 0, cnt = 0, cntLoc = 0; j < npts[1]; ++j)
406 for (
int i = 0; i < npts[0]; ++i, ++cnt)
408 if (cnt >= start && cnt < end)
410 for (
int n = 0; n < dim; ++n)
413 (values[2 + n] + i * delta1[n]) *
415 (values[2 + 3 * dim + n] + i * delta2[n]) *
424 ppe.push_back(npts[0]);
425 ppe.push_back(npts[1]);
429 m_f->m_fieldPts->SetPointsPerEdge(ppe);
431 else if (
m_config[
"box"].as<string>().compare(
"NotSet") != 0)
433 vector<NekDouble> values;
436 "Failed to interpret box string");
438 ASSERTL0(values.size() == 9,
"box string should contain 9 values "
439 "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
447 int totpts = npts[0] * npts[1] * npts[2];
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;
459 for (
int i = 0; i < dim; ++i)
462 delta[i] = (values[4 + 2 * i] - values[3 + 2 * i]) / (npts[i] - 1);
465 for (
int k = 0, cnt = 0, cntLoc = 0; k < npts[2]; ++k)
467 for (
int j = 0; j < npts[1]; ++j)
469 for (
int i = 0; i < npts[0]; ++i, ++cnt)
471 if (cnt >= start && cnt < end)
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];
483 ppe.push_back(npts[0]);
484 ppe.push_back(npts[1]);
485 ppe.push_back(npts[2]);
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);
496 ASSERTL0(
false,
"Missing target points for ProcessInterpPoints.");
501 vector<MultiRegions::ExpListSharedPtr> &field0,
505 boost::ignore_unused(def_value);
507 ASSERTL0(pts->GetNFields() == field0.size(),
"ptField has too few fields");
509 int nfields = field0.size();
512 if (
m_f->m_comm->GetRank() == 0)
519 if (
m_f->m_comm->GetRank() == 0)
524 for (
int f = 0; f < nfields; ++f)
526 for (
int i = 0; i < pts->GetNpoints(); ++i)
528 if (pts->GetPointVal(f, i) > clamp_up)
530 pts->SetPointVal(f, i, clamp_up);
532 else if (pts->GetPointVal(f, i) < clamp_low)
534 pts->SetPointVal(f, i, clamp_low);
543 int dim = pts->GetDim();
544 int nq1 = pts->GetNpoints();
550 vector<NekDouble> values;
552 "Failed to interpret cp string");
554 ASSERTL0(values.size() == 2,
"cp string should contain 2 values "
555 "p0 and q (=1/2 rho u^2)");
558 qinv = 1.0 / values[1];
560 for (
int i = 0; i < pts->GetNFields(); ++i)
562 if (boost::iequals(pts->GetFieldName(i),
"p"))
567 if (boost::iequals(pts->GetFieldName(i),
"u") ||
568 boost::iequals(pts->GetFieldName(i),
"v") ||
569 boost::iequals(pts->GetFieldName(i),
"w"))
579 WARNINGL0(
false,
"Did not find velocity components for Cp0");
584 WARNINGL0(
false,
"Failed to find 'p' field to determine cp0");
590 for (f = 0; f < 2; ++f)
595 for (r = 0; r < nq1; r++)
599 data[0][r] = qinv * (pts->GetPointVal(dim + pfield, r) - p0);
604 for (
int i = 0; i < velid.size(); ++i)
606 q += 0.5 * pts->GetPointVal(dim + velid[i], r) *
607 pts->GetPointVal(dim + velid[i], r);
610 qinv * (pts->GetPointVal(dim + pfield, r) +
q - p0);
617 pts->AddField(data[0],
"Cp");
618 m_f->m_variables.push_back(
"Cp");
621 pts->AddField(data[1],
"Cp0");
622 m_f->m_variables.push_back(
"Cp0");
628 const int goal)
const
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
#define WARNINGL0(condition, msg)
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.
std::map< std::string, ConfigOption > m_config
List of configuration values.
void InterpolateFieldToPts(std::vector< MultiRegions::ExpListSharedPtr > &field0, LibUtilities::PtsFieldSharedPtr &pts, NekDouble clamp_low, NekDouble clamp_up, NekDouble def_value)
virtual ~ProcessInterpPoints()
void CreateFieldPts(po::variables_map &vm)
void PrintProgressbar(const int position, const int goal) const
static ModuleKey className
ProcessInterpPoints(FieldSharedPtr f)
static std::shared_ptr< Module > create(FieldSharedPtr f)
Creates an instance of this class.
virtual void v_Process(po::variables_map &vm) override
Write mesh to output file.
Abstract base class for processing modules.
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.
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
static MeshGraphSharedPtr Read(const LibUtilities::SessionReaderSharedPtr pSession, LibUtilities::DomainRangeShPtr rng=LibUtilities::NullDomainRangeShPtr, bool fillGraph=true, SpatialDomains::MeshGraphSharedPtr partitionedGraph=nullptr)
std::shared_ptr< Field > FieldSharedPtr
std::pair< ModuleType, std::string > ModuleKey
ModuleFactory & GetModuleFactory()
int PrintProgressbar(const int position, const int goal, const std::string message, int lastprogress=-1)
Prints a progressbar.
static FieldMetaDataMap NullFieldMetaDataMap
std::shared_ptr< PtsField > PtsFieldSharedPtr
std::shared_ptr< DomainRange > DomainRangeShPtr
std::shared_ptr< CsvIO > CsvIOSharedPtr
std::shared_ptr< PtsIO > PtsIOSharedPtr
CommFactory & GetCommFactory()
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode .
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode .
@ eFourier
Fourier Expansion .
std::map< int, ExpansionInfoShPtr > ExpansionInfoMap
std::vector< double > q(NPUPPER *NPUPPER)
The above copyright notice and this permission notice shall be included.
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
void Zero(int n, T *x, const int incx)
Zero vector.
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
Represents a command-line configuration option.