46 #include <boost/lexical_cast.hpp>
47 #include <boost/math/special_functions/fpclassify.hpp>
54 ModuleKey ProcessInterpPoints::className =
57 ProcessInterpPoints::create,
58 "Interpolates a set of points to another, requires fromfld and "
59 "fromxml to be defined, a line, plane or block of points can be "
66 false,
"NotSet",
"Xml file from which to interpolate field");
69 "Need to specify fromxml=file.xml");
72 false,
"NotSet",
"Fld file from which to interpolate field");
75 "Need to specify fromfld=file.fld ");
78 ConfigOption(
false,
"-10000000",
"Lower bound for interpolation value");
80 ConfigOption(
false,
"10000000",
"Upper bound for interpolation value");
82 ConfigOption(
false,
"0",
"Default value if point is outside domain");
84 ConfigOption(
false,
"NotSet",
"Specify a line of N points using "
85 "line=N,x0,y0,z0,z1,y1,z1");
87 false,
"NotSet",
"Specify a plane of N1 x N2 points using "
88 "plane=N1,N2,x0,y0,z0,z1,y1,z1,x2,y2,z2,x3,"
91 false,
"NotSet",
"Specify a rectangular box of N1 x N2 x N3 points "
92 "using a box of points limited by box="
93 "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
97 "Parameters p0 and q to determine pressure coefficients "
98 "(box only currently)");
109 if (
m_f->m_comm->TreatAsRankZero())
111 cout <<
"ProcessInterpPoints: interpolating to points..." << endl;
115 int rank =
m_f->m_comm->GetRank();
116 int nprocs =
m_f->m_comm->GetSize();
122 if (
m_config[
"line"].as<string>().compare(
"NotSet") != 0)
124 string help =
m_config[
"line"].as<
string>();
125 vector<NekDouble> values;
127 m_config[
"line"].as<string>().c_str(), values),
128 "Failed to interpret line string");
131 "line string should contain 2 Dim+1 values "
132 "N,x0,y0,z0,x1,y1,z1");
135 ASSERTL0(std::modf(values[0], &tmp) == 0.0,
"N is not an integer");
136 ASSERTL0(values[0] > 1,
"N is not a valid number");
138 int dim = (values.size() - 1) / 2;
139 int npts = values[0];
142 for (
int i = 0; i < dim; ++i)
147 for (
int i = 0; i <
npts; ++i)
151 i / ((
NekDouble)(npts - 1)) * (values[dim + 1] - values[1]);
154 pts[1][i] = values[2] +
156 (values[dim + 2] - values[2]);
160 pts[2][i] = values[3] +
162 (values[dim + 3] - values[3]);
173 m_f->m_fieldPts->SetPointsPerEdge(ppe);
175 else if (
m_config[
"plane"].as<string>().compare(
"NotSet") != 0)
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");
184 "plane string should contain 4 Dim+2 values "
185 "N1,N2,x0,y0,z0,x1,y1,z1,x2,y2,z2,x3,y3,z3");
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");
191 ASSERTL0(values[0] > 1,
"N1 is not a valid number");
192 ASSERTL0(values[1] > 1,
"N2 is not a valid number");
194 int dim = (values.size() - 2) / 4;
196 int npts1 = values[0];
197 int npts2 = values[1];
201 int totpts = npts1 * npts2;
202 int nlocpts = totpts / nprocs;
204 if (rank < nprocs - 1)
206 for (
int i = 0; i < dim; ++i)
214 for (
int j = 0; j < npts2; ++j)
216 for (
int i = 0; i < npts1; ++i)
219 if ((cnt >= rank * nlocpts) &&
220 (cnt < (rank + 1) * nlocpts))
225 (values[dim + 2] - values[2])) *
227 (values[3 * dim + 2] +
229 (values[2 * dim + 2] -
230 values[3 * dim + 2])) *
236 (values[dim + 3] - values[3])) *
238 (values[3 * dim + 3] +
240 (values[2 * dim + 3] -
241 values[3 * dim + 3])) *
249 (values[dim + 4] - values[4])) *
251 (values[3 * dim + 4] +
253 (values[2 * dim + 4] -
254 values[3 * dim + 4])) *
265 totpts = totpts - rank * nlocpts;
267 for (
int i = 0; i < dim; ++i)
275 for (
int j = 0; j < npts2; ++j)
277 for (
int i = 0; i < npts1; ++i)
280 if (cnt >= rank * nlocpts)
285 (values[dim + 2] - values[2])) *
287 (values[3 * dim + 2] +
289 (values[2 * dim + 2] -
290 values[3 * dim + 2])) *
296 (values[dim + 3] - values[3])) *
298 (values[3 * dim + 3] +
300 (values[2 * dim + 3] -
301 values[3 * dim + 3])) *
309 (values[dim + 4] - values[4])) *
311 (values[3 * dim + 4] +
313 (values[2 * dim + 4] -
314 values[3 * dim + 4])) *
325 ppe.push_back(npts1);
326 ppe.push_back(npts2);
331 m_f->m_fieldPts->SetPointsPerEdge(ppe);
333 else if (
m_config[
"box"].as<string>().compare(
"NotSet") != 0)
335 string help =
m_config[
"box"].as<
string>();
336 vector<NekDouble> values;
338 m_config[
"box"].as<string>().c_str(), values),
339 "Failed to interpret box string");
342 "box string should contain 9 values "
343 "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
347 int npts1 = values[0];
348 int npts2 = values[1];
349 int npts3 = values[2];
353 int totpts = npts1 * npts2 * npts3;
354 int nlocpts = totpts / nprocs;
356 if (rank < nprocs - 1)
360 for (
int i = 0; i < dim; ++i)
368 for (
int k = 0; k < npts3; ++k)
370 for (
int j = 0; j < npts2; ++j)
372 for (
int i = 0; i < npts1; ++i)
374 if ((cnt >= rank * nlocpts) &&
375 (cnt < (rank + 1) * nlocpts))
377 pts[0][cntloc] = values[3] +
379 (values[4] - values[3]);
380 pts[1][cntloc] = values[5] +
382 (values[6] - values[5]);
383 pts[2][cntloc] = values[7] +
385 (values[8] - values[7]);
395 totpts = totpts - rank * nlocpts;
397 for (
int i = 0; i < dim; ++i)
405 for (
int k = 0; k < npts3; ++k)
407 for (
int j = 0; j < npts2; ++j)
409 for (
int i = 0; i < npts1; ++i)
411 if (cnt >= rank * nlocpts)
413 pts[0][cntloc] = values[3] +
415 (values[4] - values[3]);
416 pts[1][cntloc] = values[5] +
418 (values[6] - values[5]);
419 pts[2][cntloc] = values[7] +
421 (values[8] - values[7]);
431 ppe.push_back(npts1);
432 ppe.push_back(npts2);
433 ppe.push_back(npts3);
438 m_f->m_fieldPts->SetPointsPerEdge(ppe);
439 vector<NekDouble> boxdim;
440 boxdim.assign(&values[3], &values[3] + 6);
441 m_f->m_fieldPts->SetBoxSize(boxdim);
447 std::vector<std::string> files;
450 fromField->m_session =
457 int coordim =
m_f->m_fieldPts->GetDim();
458 int npts =
m_f->m_fieldPts->GetNpoints();
460 m_f->m_fieldPts->GetPts(pts);
462 rng->m_checkShape =
false;
470 rng->m_doZrange =
true;
473 if (rng->m_zmax == rng->m_zmin)
479 rng->m_doYrange =
true;
483 rng->m_doXrange =
true;
488 ASSERTL0(
false,
"too many values specfied in range");
497 fromField->m_graph->GetExpansions();
500 SpatialDomains::ExpansionMap::const_iterator expIt;
503 for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
505 ElementGIDs[i++] = expIt->second->m_geomShPtr->GetGlobalID();
510 ASSERTL0(i > 0,
"No elements are set. Are the interpolated points "
511 "wihtin the domain given by the xml files?");
513 string fromfld =
m_config[
"fromfld"].as<
string>();
514 m_f->FieldIOForFile(fromfld)->Import(
515 fromfld, fromField->m_fielddef, fromField->m_data,
518 int NumHomogeneousDir = fromField->m_fielddef[0]->m_numHomogeneousDir;
522 fromField->m_graph->SetExpansions(fromField->m_fielddef);
524 int nfields = fromField->m_fielddef[0]->m_fields.size();
526 fromField->m_exp.resize(nfields);
527 fromField->m_exp[0] = fromField->SetUpFirstExpList(NumHomogeneousDir,
true);
529 m_f->m_exp.resize(nfields);
532 for (i = 1; i < nfields; ++i)
534 fromField->m_exp[i] = fromField->AppendExpList(NumHomogeneousDir);
538 for (
int j = 0; j < nfields; ++j)
540 for (i = 0; i < fromField->m_fielddef.size(); i++)
542 fromField->m_exp[j]->ExtractDataToCoeffs(
543 fromField->m_fielddef[i], fromField->m_data[i],
544 fromField->m_fielddef[0]->m_fields[j],
545 fromField->m_exp[j]->UpdateCoeffs());
547 fromField->m_exp[j]->BwdTrans(fromField->m_exp[j]->GetCoeffs(),
548 fromField->m_exp[j]->UpdatePhys());
551 m_f->m_fieldPts->AddField(newPts,
552 fromField->m_fielddef[0]->m_fields[j]);
560 clamp_up, def_value);
562 if (!boost::iequals(
m_config[
"cp"].as<string>(),
"NotSet"))
569 vector<MultiRegions::ExpListSharedPtr> &field0,
575 ASSERTL0(pts->GetNFields() >= field0.size(),
"ptField has too few fields");
577 int nfields = field0.size();
580 if (
m_f->m_comm->GetRank() == 0)
586 if (
m_f->m_comm->GetRank() == 0)
591 for (
int f = 0; f < nfields; ++f)
593 for (
int i = 0; i < pts->GetNpoints(); ++i)
595 if (pts->GetPointVal(f, i) > clamp_up)
597 pts->SetPointVal(f, i, clamp_up);
599 else if (pts->GetPointVal(f, i) < clamp_low)
601 pts->SetPointVal(f, i, clamp_low);
610 int dim = pts->GetDim();
611 int nq1 = pts->GetNpoints();
617 vector<NekDouble> values;
619 m_config[
"cp"].as<string>().c_str(),values),
620 "Failed to interpret cp string");
623 "cp string should contain 2 values "
624 "p0 and q (=1/2 rho u^2)");
627 qinv = 1.0/values[1];
629 for(
int i = 0; i < pts->GetNFields(); ++i)
631 if(boost::iequals(pts->GetFieldName(i),
"p"))
636 if(boost::iequals(pts->GetFieldName(i),
"u")||
637 boost::iequals(pts->GetFieldName(i),
"v")||
638 boost::iequals(pts->GetFieldName(i),
"w"))
648 WARNINGL0(
false,
"Did not find velocity components for Cp0");
653 WARNINGL0(
false,
"Failed to find 'p' field to determine cp0");
659 for (f = 0; f < 2; ++f)
664 for (r = 0; r < nq1; r++)
668 data[0][r] = qinv*(pts->GetPointVal(dim + pfield, r) - p0);
673 for(
int i = 0; i < velid.size(); ++i)
675 q += 0.5*pts->GetPointVal(dim + velid[i], r)*
676 pts->GetPointVal(dim + velid[i], r);
678 data[1][r] = qinv*(pts->GetPointVal(dim + pfield, r)+q - p0);
685 pts->AddField(data[0],
"Cp");
688 pts->AddField(data[1],
"Cp0");
694 const int goal)
const
map< string, ConfigOption > m_config
List of configuration values.
#define ASSERTL0(condition, msg)
virtual ~ProcessInterpPoints()
static bool GenerateOrderedStringVector(const char *const str, std::vector< std::string > &vec)
A class that contains algorithms for interpolation between pts fields, expansions and different meshe...
void SetProgressCallback(FuncPointerT func, ObjectPointerT obj)
sets a callback funtion which gets called every time the interpolation progresses ...
int PrintProgressbar(const int position, const int goal, const string message, int lastprogress=-1)
Prints a progressbar.
static boost::shared_ptr< MeshGraph > Read(const LibUtilities::SessionReaderSharedPtr &pSession, DomainRangeShPtr &rng=NullDomainRangeShPtr)
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
Represents a command-line configuration option.
FIELD_UTILS_EXPORT void Interpolate(const LibUtilities::PtsFieldSharedPtr ptsInField, LibUtilities::PtsFieldSharedPtr &ptsOutField)
Interpolate from a pts field to a pts field.
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.
pair< ModuleType, string > ModuleKey
boost::shared_ptr< PtsField > PtsFieldSharedPtr
static SessionReaderSharedPtr CreateInstance(int argc, char *argv[])
Creates an instance of the SessionReader class.
virtual void Process(po::variables_map &vm)
Write mesh to output file.
boost::shared_ptr< Field > FieldSharedPtr
#define WARNINGL0(condition, msg)
boost::shared_ptr< DomainRange > DomainRangeShPtr
static PtsFieldSharedPtr NullPtsField
void PrintProgressbar(const int position, const int goal) const
static bool GenerateUnOrderedVector(const char *const str, std::vector< NekDouble > &vec)
Abstract base class for processing modules.
static FieldMetaDataMap NullFieldMetaDataMap
void InterpolateFieldToPts(vector< MultiRegions::ExpListSharedPtr > &field0, LibUtilities::PtsFieldSharedPtr &pts, NekDouble clamp_low, NekDouble clamp_up, NekDouble def_value)
std::map< int, ExpansionShPtr > ExpansionMap
ModuleFactory & GetModuleFactory()
FieldSharedPtr m_f
Field object.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.