39 #include <boost/core/ignore_unused.hpp> 40 #include <boost/geometry.hpp> 41 #include <boost/lexical_cast.hpp> 42 #include <boost/math/special_functions/fpclassify.hpp> 53 namespace bg = boost::geometry;
54 namespace bgi = boost::geometry::index;
61 ModuleKey ProcessInterpPoints::className =
64 ProcessInterpPoints::create,
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 false,
"NotSet",
"Pts file to which interpolate field");
79 false,
"NotSet",
"Specify a line of N points using " 80 "line=N,x0,y0,z0,z1,y1,z1");
82 false,
"NotSet",
"Specify a plane of N1 x N2 points using " 83 "plane=N1,N2,x0,y0,z0,z1,y1,z1,x2,y2,z2,x3,y3,z3");
85 false,
"NotSet",
"Specify a rectangular box of N1 x N2 x N3 points " 86 "using a box of points limited by box=" 87 "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
90 ConfigOption(
false,
"-10000000",
"Lower bound for interpolation value");
92 ConfigOption(
false,
"10000000",
"Upper bound for interpolation value");
94 ConfigOption(
false,
"0",
"Default value if point is outside domain");
98 "Parameters p0 and q to determine pressure coefficients");
110 std::vector<std::string> files;
114 char *argv[] = {
const_cast<char *
>(
"FieldConvert"),
nullptr };
115 fromField->m_session =
122 int coordim =
m_f->m_fieldPts->GetDim();
123 int npts =
m_f->m_fieldPts->GetNpoints();
124 std::vector<std::string> fieldNames =
m_f->m_fieldPts->GetFieldNames();
125 for (
auto &it : fieldNames)
127 m_f->m_fieldPts->RemoveField(it);
131 m_f->m_fieldPts->GetPts(pts);
133 rng->m_checkShape =
false;
141 rng->m_doZrange =
true;
144 if (rng->m_zmax == rng->m_zmin)
151 rng->m_doYrange =
true;
156 rng->m_doXrange =
true;
168 fromField->m_graph->GetExpansions();
172 for (
auto &expIt : expansions)
174 ElementGIDs[i++] = expIt.second->m_geomShPtr->GetGlobalID();
178 ASSERTL0(i > 0,
"No elements are set. Are the interpolated points " 179 "wihtin the domain given by the xml files?");
180 string fromfld =
m_config[
"fromfld"].as<
string>();
181 m_f->FieldIOForFile(fromfld)->Import(
182 fromfld, fromField->m_fielddef, fromField->m_data,
184 int NumHomogeneousDir = fromField->m_fielddef[0]->m_numHomogeneousDir;
187 fromField->m_graph->SetExpansions(fromField->m_fielddef);
188 int nfields = fromField->m_fielddef[0]->m_fields.size();
189 fromField->m_exp.resize(nfields);
190 fromField->m_exp[0] = fromField->SetUpFirstExpList(NumHomogeneousDir,
true);
191 m_f->m_exp.resize(nfields);
193 for (i = 1; i < nfields; ++i)
195 fromField->m_exp[i] = fromField->AppendExpList(NumHomogeneousDir);
198 for (
int j = 0; j < nfields; ++j)
200 for (i = 0; i < fromField->m_fielddef.size(); i++)
202 fromField->m_exp[j]->ExtractDataToCoeffs(
203 fromField->m_fielddef[i], fromField->m_data[i],
204 fromField->m_fielddef[0]->m_fields[j],
205 fromField->m_exp[j]->UpdateCoeffs());
207 fromField->m_exp[j]->BwdTrans(fromField->m_exp[j]->GetCoeffs(),
208 fromField->m_exp[j]->UpdatePhys());
210 m_f->m_fieldPts->AddField(newPts,
211 fromField->m_fielddef[0]->m_fields[j]);
212 m_f->m_variables.push_back(fromField->m_fielddef[0]->m_fields[j]);
220 clamp_up, def_value);
222 if (!boost::iequals(
m_config[
"cp"].as<string>(),
"NotSet"))
230 boost::ignore_unused(vm);
232 int rank =
m_f->m_comm->GetRank();
233 int nprocs =
m_f->m_comm->GetSize();
236 if (
m_config[
"topts"].as<string>().compare(
"NotSet") != 0)
238 string inFile =
m_config[
"topts"].as<
string>();
240 if (boost::filesystem::path(inFile).extension() ==
".pts")
245 ptsIO->Import(inFile,
m_f->m_fieldPts);
247 else if (boost::filesystem::path(inFile).extension() ==
".csv")
252 csvIO->Import(inFile,
m_f->m_fieldPts);
256 ASSERTL0(
false,
"unknown topts file type");
260 else if (
m_config[
"line"].as<string>().compare(
"NotSet") != 0)
262 vector<NekDouble> values;
264 m_config[
"line"].as<string>(), values),
265 "Failed to interpret line string");
268 "line string should contain 2*Dim+1 values " 269 "N,x0,y0,z0,x1,y1,z1");
272 ASSERTL0(std::modf(values[0], &tmp) == 0.0,
"N is not an integer");
273 ASSERTL0(values[0] > 1,
"N is not a valid number");
275 int dim = (values.size() - 1) / 2;
276 int npts = values[0];
279 int ptsPerProc = npts / nprocs;
280 int extraPts = (rank < nprocs - 1) ? 0: npts % nprocs;
281 int locPts = ptsPerProc + extraPts;
282 int start = rank * ptsPerProc;
283 int end = start + locPts;
287 for (
int i = 0; i < dim; ++i)
290 delta[i] = (values[dim + i + 1] - values[ i + 1]) / (npts - 1);
295 for (
int i = 0, cntLoc = 0; i < npts; ++i)
297 if (i >= start && i < end)
299 for (
int n = 0; n < dim; ++n)
301 pts[n][cntLoc] = values[n+1] + i * delta[n];
313 m_f->m_fieldPts->SetPointsPerEdge(ppe);
315 else if (
m_config[
"plane"].as<string>().compare(
"NotSet") != 0)
317 vector<NekDouble> values;
319 m_config[
"plane"].as<string>(), values),
320 "Failed to interpret plane string");
323 "plane string should contain 4 Dim+2 values " 324 "N1,N2,x0,y0,z0,x1,y1,z1,x2,y2,z2,x3,y3,z3");
327 ASSERTL0(std::modf(values[0], &tmp) == 0.0,
"N1 is not an integer");
328 ASSERTL0(std::modf(values[1], &tmp) == 0.0,
"N2 is not an integer");
330 ASSERTL0(values[0] > 1,
"N1 is not a valid number");
331 ASSERTL0(values[1] > 1,
"N2 is not a valid number");
333 int dim = (values.size() - 2) / 4;
339 int totpts = npts[0] * npts[1];
342 int ptsPerProc = totpts / nprocs;
343 int extraPts = (rank < nprocs - 1) ? 0: totpts % nprocs;
344 int locPts = ptsPerProc + extraPts;
345 int start = rank * ptsPerProc;
346 int end = start + locPts;
351 for (
int i = 0; i < dim; ++i)
354 delta1[i] = (values[2+1*dim + i] - values[2+0*dim + i])/(npts[0]-1);
355 delta2[i] = (values[2+2*dim + i] - values[2+3*dim + i])/(npts[0]-1);
358 for (
int j = 0, cnt = 0, cntLoc = 0; j < npts[1]; ++j)
360 for (
int i = 0; i < npts[0]; ++i, ++cnt)
362 if (cnt >= start && cnt < end)
364 for (
int n = 0; n < dim; ++n)
367 (values[2+n] + i * delta1[n]) *
369 (values[2 + 3*dim + n] + i * delta2[n]) *
378 ppe.push_back(npts[0]);
379 ppe.push_back(npts[1]);
384 m_f->m_fieldPts->SetPointsPerEdge(ppe);
386 else if (
m_config[
"box"].as<string>().compare(
"NotSet") != 0)
388 vector<NekDouble> values;
390 m_config[
"box"].as<string>(), values),
391 "Failed to interpret box string");
394 "box string should contain 9 values " 395 "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
403 int totpts = npts[0]*npts[1]*npts[2];
409 int ptsPerProc = totpts / nprocs;
410 int extraPts = (rank < nprocs - 1) ? 0: totpts % nprocs;
411 int locPts = ptsPerProc + extraPts;
412 int start = rank * ptsPerProc;
413 int end = start + locPts;
415 for (
int i = 0; i < dim; ++i)
418 delta[i] = (values[4 + 2*i] - values[3 + 2*i]) / (npts[i] - 1);
423 for (
int k = 0, cnt = 0, cntLoc = 0; k < npts[2]; ++k)
425 for (
int j = 0; j < npts[1]; ++j)
427 for (
int i = 0; i < npts[0]; ++i, ++cnt)
429 if (cnt >= start && cnt < end)
431 pts[0][cntLoc] = values[3] + i * delta[0];
432 pts[1][cntLoc] = values[5] + j * delta[1];
433 pts[2][cntLoc] = values[7] + k * delta[2];
441 ppe.push_back(npts[0]);
442 ppe.push_back(npts[1]);
443 ppe.push_back(npts[2]);
448 m_f->m_fieldPts->SetPointsPerEdge(ppe);
449 vector<NekDouble> boxdim;
450 boxdim.assign(&values[3], &values[3] + 6);
451 m_f->m_fieldPts->SetBoxSize(boxdim);
456 "Missing target points for ProcessInterpPoints.");
461 vector<MultiRegions::ExpListSharedPtr> &field0,
467 boost::ignore_unused(def_value);
469 ASSERTL0(pts->GetNFields() == field0.size(),
"ptField has too few fields");
471 int nfields = field0.size();
474 if (
m_f->m_comm->GetRank() == 0)
480 if (
m_f->m_comm->GetRank() == 0)
485 for (
int f = 0; f < nfields; ++f)
487 for (
int i = 0; i < pts->GetNpoints(); ++i)
489 if (pts->GetPointVal(f, i) > clamp_up)
491 pts->SetPointVal(f, i, clamp_up);
493 else if (pts->GetPointVal(f, i) < clamp_low)
495 pts->SetPointVal(f, i, clamp_low);
504 int dim = pts->GetDim();
505 int nq1 = pts->GetNpoints();
511 vector<NekDouble> values;
513 "Failed to interpret cp string");
516 "cp string should contain 2 values " 517 "p0 and q (=1/2 rho u^2)");
520 qinv = 1.0/values[1];
522 for(
int i = 0; i < pts->GetNFields(); ++i)
524 if(boost::iequals(pts->GetFieldName(i),
"p"))
529 if(boost::iequals(pts->GetFieldName(i),
"u")||
530 boost::iequals(pts->GetFieldName(i),
"v")||
531 boost::iequals(pts->GetFieldName(i),
"w"))
541 WARNINGL0(
false,
"Did not find velocity components for Cp0");
546 WARNINGL0(
false,
"Failed to find 'p' field to determine cp0");
552 for (f = 0; f < 2; ++f)
557 for (r = 0; r < nq1; r++)
561 data[0][r] = qinv*(pts->GetPointVal(dim + pfield, r) - p0);
566 for(
int i = 0; i < velid.size(); ++i)
568 q += 0.5*pts->GetPointVal(dim + velid[i], r)*
569 pts->GetPointVal(dim + velid[i], r);
571 data[1][r] = qinv*(pts->GetPointVal(dim + pfield, r)+q - p0);
578 pts->AddField(data[0],
"Cp");
579 m_f->m_variables.push_back(
"Cp");
582 pts->AddField(data[1],
"Cp0");
583 m_f->m_variables.push_back(
"Cp0");
589 const int goal)
const
#define ASSERTL0(condition, msg)
virtual ~ProcessInterpPoints()
A class that contains algorithms for interpolation between pts fields, expansions and different meshe...
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
int PrintProgressbar(const int position, const int goal, const std::string message, int lastprogress=-1)
Prints a progressbar.
std::map< std::string, ConfigOption > m_config
List of configuration values.
void SetProgressCallback(FuncPointerT func, ObjectPointerT obj)
sets a callback funtion which gets called every time the interpolation progresses ...
Represents a command-line configuration option.
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
void CreateFieldPts(po::variables_map &vm)
std::shared_ptr< Field > FieldSharedPtr
std::shared_ptr< DomainRange > DomainRangeShPtr
CommFactory & GetCommFactory()
static SessionReaderSharedPtr CreateInstance(int argc, char *argv[])
Creates an instance of the SessionReader class.
std::shared_ptr< PtsIO > PtsIOSharedPtr
virtual void Process(po::variables_map &vm)
Write mesh to output file.
std::shared_ptr< PtsField > PtsFieldSharedPtr
std::pair< ModuleType, std::string > ModuleKey
#define WARNINGL0(condition, msg)
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
FIELD_UTILS_EXPORT void Interpolate(const std::vector< MultiRegions::ExpListSharedPtr > expInField, std::vector< MultiRegions::ExpListSharedPtr > &expOutField, NekDouble def_value=0.0)
Interpolate from an expansion to an expansion.
void PrintProgressbar(const int position, const int goal) const
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Abstract base class for processing modules.
std::shared_ptr< CsvIO > CsvIOSharedPtr
static FieldMetaDataMap NullFieldMetaDataMap
void InterpolateFieldToPts(vector< MultiRegions::ExpListSharedPtr > &field0, LibUtilities::PtsFieldSharedPtr &pts, NekDouble clamp_low, NekDouble clamp_up, NekDouble def_value)
static MeshGraphSharedPtr Read(const LibUtilities::SessionReaderSharedPtr pSession, DomainRangeShPtr rng=NullDomainRangeShPtr, bool fillGraph=true)
std::map< int, ExpansionShPtr > ExpansionMap
ModuleFactory & GetModuleFactory()
FieldSharedPtr m_f
Field object.