43 #include <boost/math/special_functions/fpclassify.hpp>
49 ModuleKey ProcessInterpPoints::className =
52 ProcessInterpPoints::create,
53 "Interpolates a set of points to another, requires fromfld and "
54 "fromxml to be defined, a line, plane or block of points can be "
63 "Xml file from which to interpolate field");
66 "Need to specify fromxml=file.xml");
70 "Fld file from which to interpolate field");
73 "Need to specify fromfld=file.fld ");
77 "Lower bound for interpolation value");
80 "Upper bound for interpolation value");
83 "Default value if point is outside domain");
86 "Specify a line of N points using "
87 "line=N,x0,y0,z0,z1,y1,z1");
90 "Specify a plane of N1 x N2 points using "
91 "plane=N1,N2,x0,y0,z0,z1,y1,z1,x2,y2,z2,x3,"
95 "Specify a rectangular box of N1 x N2 x N3 points "
96 "using a box of points limited by box="
97 "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
101 "Parameters p0 and q to determine pressure coefficients "
102 "(box only currently)");
114 if(
m_f->m_comm->TreatAsRankZero())
116 cout <<
"ProcessInterpPoints: interpolating to points..." << endl;
120 int rank =
m_f->m_comm->GetRank();
121 int nprocs =
m_f->m_comm->GetSize();
127 if(
m_config[
"line"].as<string>().compare(
"NotSet") != 0)
129 string help =
m_config[
"line"].as<
string>();
130 vector<NekDouble> values;
132 m_config[
"line"].as<string>().c_str(),values),
133 "Failed to interpret line string");
136 "line string should contain 2 Dim+1 values "
137 "N,x0,y0,z0,x1,y1,z1");
140 ASSERTL0(std::modf(values[0], &tmp) == 0.0,
"N is not an integer");
141 ASSERTL0(values[0] > 1,
"N is not a valid number");
143 int dim = (values.size()-1)/2;
144 int npts = values[0];
147 for(
int i = 0; i < dim; ++i)
152 for(
int i = 0; i <
npts; ++i)
154 pts[0][i] = values[1]
155 + i/((
NekDouble)(npts-1))*(values[dim+1] - values[1]);
158 pts[1][i] = values[2]
159 + i/((
NekDouble)(npts-1))*(values[dim+2] - values[2]);
163 pts[2][i] = values[3]
164 + i/((
NekDouble)(npts-1))*(values[dim+3] - values[3]);
173 m_f->m_fieldPts->SetPointsPerEdge(ppe);
176 else if(
m_config[
"plane"].as<string>().compare(
"NotSet") != 0)
178 string help =
m_config[
"plane"].as<
string>();
179 vector<NekDouble> values;
181 m_config[
"plane"].as<string>().c_str(),values),
182 "Failed to interpret plane string");
185 "plane string should contain 4 Dim+2 values "
186 "N1,N2,x0,y0,z0,x1,y1,z1,x2,y2,z2,x3,y3,z3");
189 ASSERTL0(std::modf(values[0], &tmp) == 0.0,
"N1 is not an integer");
190 ASSERTL0(std::modf(values[1], &tmp) == 0.0,
"N2 is not an integer");
192 ASSERTL0(values[0] > 1,
"N1 is not a valid number");
193 ASSERTL0(values[1] > 1,
"N2 is not a valid number");
195 int dim = (values.size()-2)/4;
197 int npts1 = values[0];
198 int npts2 = values[1];
202 int totpts = npts1*npts2;
203 int nlocpts = totpts/nprocs;
207 for(
int i = 0; i < dim; ++i)
215 for(
int j = 0; j < npts2; ++j)
217 for(
int i = 0; i < npts1; ++i)
220 if((cnt >= rank*nlocpts)&&(cnt < (rank+1)*nlocpts))
223 (values[2] + i/((
NekDouble)(npts1-1))*(values[dim+2] - values[2]))*(1.0-j/((
NekDouble)(npts2-1))) +
224 (values[3*dim+2] + i/((
NekDouble)(npts1-1))*(values[2*dim+2] - values[3*dim+2]))*(j/((
NekDouble)(npts2-1)));
227 (values[3] + i/((
NekDouble)(npts1-1))*(values[dim+3] - values[3]))*(1.0-j/((
NekDouble)(npts2-1))) +
228 (values[3*dim+3] + i/((
NekDouble)(npts1-1))*(values[2*dim+3] - values[3*dim+3]))*(j/((
NekDouble)(npts2-1)));
233 (values[4] + i/((
NekDouble)(npts1-1))*(values[dim+4] - values[4]))*(1.0-j/((
NekDouble)(npts2-1))) +
234 (values[3*dim+4] + i/((
NekDouble)(npts1-1))*(values[2*dim+4] - values[3*dim+4]))*(j/((
NekDouble)(npts2-1)));
244 totpts = totpts - rank*nlocpts;
246 for(
int i = 0; i < dim; ++i)
254 for(
int j = 0; j < npts2; ++j)
256 for(
int i = 0; i < npts1; ++i)
259 if(cnt >= rank*nlocpts)
262 (values[2] + i/((
NekDouble)(npts1-1))*(values[dim+2] - values[2]))*(1.0-j/((
NekDouble)(npts2-1))) +
263 (values[3*dim+2] + i/((
NekDouble)(npts1-1))*(values[2*dim+2] - values[3*dim+2]))*(j/((
NekDouble)(npts2-1)));
266 (values[3] + i/((
NekDouble)(npts1-1))*(values[dim+3] - values[3]))*(1.0-j/((
NekDouble)(npts2-1))) +
267 (values[3*dim+3] + i/((
NekDouble)(npts1-1))*(values[2*dim+3] - values[3*dim+3]))*(j/((
NekDouble)(npts2-1)));
272 (values[4] + i/((
NekDouble)(npts1-1))*(values[dim+4] - values[4]))*(1.0-j/((
NekDouble)(npts2-1))) +
273 (values[3*dim+4] + i/((
NekDouble)(npts1-1))*(values[2*dim+4] - values[3*dim+4]))*(j/((
NekDouble)(npts2-1)));
283 ppe.push_back(npts1);
284 ppe.push_back(npts2);
287 m_f->m_fieldPts->SetPointsPerEdge(ppe);
290 else if(
m_config[
"box"].as<string>().compare(
"NotSet") != 0)
292 string help =
m_config[
"box"].as<
string>();
293 vector<NekDouble> values;
295 m_config[
"box"].as<string>().c_str(),values),
296 "Failed to interpret box string");
299 "box string should contain 9 values "
300 "N1,N2,N3,xmin,xmax,ymin,ymax,zmin,zmax");
304 int npts1 = values[0];
305 int npts2 = values[1];
306 int npts3 = values[2];
310 int totpts = npts1*npts2*npts3;
311 int nlocpts = totpts/nprocs;
317 for(
int i = 0; i < dim; ++i)
325 for(
int k = 0; k < npts3; ++k)
327 for(
int j = 0; j < npts2; ++j)
329 for(
int i = 0; i < npts1; ++i)
331 if((cnt >= rank*nlocpts)&&(cnt < (rank+1)*nlocpts))
333 pts[0][cntloc] = values[3] +
334 i/((
NekDouble)(npts1-1))*(values[4]-values[3]);
335 pts[1][cntloc] = values[5] +
336 j/((
NekDouble)(npts2-1))*(values[6]-values[5]);
337 pts[2][cntloc] = values[7] +
338 k/((
NekDouble)(npts3-1))*(values[8]-values[7]);
348 totpts = totpts - rank*nlocpts;
350 for(
int i = 0; i < dim; ++i)
358 for(
int k = 0; k < npts3; ++k)
360 for(
int j = 0; j < npts2; ++j)
362 for(
int i = 0; i < npts1; ++i)
364 if(cnt >= rank*nlocpts)
366 pts[0][cntloc] = values[3] +
367 i/((
NekDouble)(npts1-1))*(values[4]-values[3]);
368 pts[1][cntloc] = values[5] +
369 j/((
NekDouble)(npts2-1))*(values[6]-values[5]);
370 pts[2][cntloc] = values[7] +
371 k/((
NekDouble)(npts3-1))*(values[8]-values[7]);
381 ppe.push_back(npts1);
382 ppe.push_back(npts2);
383 ppe.push_back(npts3);
386 m_f->m_fieldPts->SetPointsPerEdge(ppe);
387 vector<NekDouble> boxdim;
388 boxdim.assign(&values[3],&values[3]+6);
389 m_f->m_fieldPts->SetBoxSize(boxdim);
396 std::vector<std::string> files;
398 files.push_back(
m_config[
"fromxml"].as<string>());
405 int coordim =
m_f->m_fieldPts->GetDim();
406 int npts =
m_f->m_fieldPts->GetNpoints();
408 m_f->m_fieldPts->GetPts(pts);
410 rng->m_checkShape =
false;
418 rng->m_doZrange =
true;
421 if(rng->m_zmax == rng->m_zmin)
427 rng->m_doYrange =
true;
431 rng->m_doXrange =
true;
436 ASSERTL0(
false,
"too many values specfied in range");
450 SpatialDomains::ExpansionMap::const_iterator expIt;
453 for (expIt = expansions.begin(); expIt != expansions.end();
456 ElementGIDs[i++] = expIt->second->m_geomShPtr->GetGlobalID();
461 ASSERTL0(i > 0,
"No elements are set. Are the interpolated points "
462 "wihtin the domain given by the xml files?");
464 string fromfld =
m_config[
"fromfld"].as<
string>();
465 fromField->m_fld->Import(fromfld,fromField->m_fielddef,
470 int NumHomogeneousDir = fromField->m_fielddef[0]->m_numHomogeneousDir;
474 fromField->m_graph->SetExpansions(fromField->m_fielddef);
476 int nfields = fromField->m_fielddef[0]->m_fields.size();
478 fromField->m_exp.resize(nfields);
479 fromField->m_exp[0] = fromField->SetUpFirstExpList(NumHomogeneousDir,
true);
481 m_f->m_exp.resize(nfields);
484 for(i = 1; i < nfields; ++i)
486 fromField->m_exp[i] = fromField->AppendExpList(NumHomogeneousDir);
490 for(
int j = 0; j < nfields; ++j)
492 for (i = 0; i < fromField->m_fielddef.size(); i++)
494 fromField->m_exp[j]->ExtractDataToCoeffs(
495 fromField->m_fielddef[i],
496 fromField->m_data[i],
497 fromField->m_fielddef[0]->m_fields[j],
498 fromField->m_exp[j]->UpdateCoeffs());
500 fromField->m_exp[j]->BwdTrans(fromField->m_exp[j]->GetCoeffs(),
501 fromField->m_exp[j]->UpdatePhys());
504 m_f->m_fieldPts->AddField(newPts, fromField->m_fielddef[0]->m_fields[j]);
509 cout <<
"Interpolating on proc 0 [" << flush;
517 clamp_low, clamp_up, def_value, !rank);
526 vector<MultiRegions::ExpListSharedPtr> &field0,
533 int dim = pts.num_elements();
536 int nq1 = pts[0].num_elements();
540 int nfields = field0.size();
545 if(!boost::iequals(
m_config[
"cp"].as<string>(),
"NotSet"))
548 vector<NekDouble> values;
550 m_config[
"cp"].as<string>().c_str(),values),
551 "Failed to interpret cp string");
554 "cp string should contain 2 values "
555 "p0 and q (=1/2 rho u^2)");
558 qinv = 1.0/values[1];
562 for(
int i = 0; i < fPts->GetNFields(); ++i)
564 if(boost::iequals(fPts->GetFieldName(i),
"p"))
569 if(boost::iequals(fPts->GetFieldName(i),
"u")||
570 boost::iequals(fPts->GetFieldName(i),
"v")||
571 boost::iequals(fPts->GetFieldName(i),
"w"))
580 m_f->m_fieldPts->AddField(newPts,
"Cp");
586 m_f->m_fieldPts->AddField(newPts,
"Cp0");
591 WARNINGL0(
false,
"Did not find velocity components for Cp0");
596 WARNINGL0(
false,
"Failed to find 'p' field to determine cp0");
602 m_f->m_data.resize(nfields);
604 for (f = 0; f < nfields; ++f)
606 m_f->m_data[f].resize(nq1);
609 for (r = 0; r < nq1; r++)
611 coords[0] = pts[0][r];
612 coords[1] = pts[1][r];
615 coords[2] = pts[2][r];
619 elmtid = field0[0]->GetExpIndex(coords, Lcoords, 1e-3);
623 offset = field0[0]->GetPhys_Offset(field0[0]->
624 GetOffset_Elmt_Id(elmtid));
626 for (f = 0; f < field0.size(); ++f)
629 value = field0[f]->GetExp(elmtid)->
630 StdPhysEvaluate(Lcoords, field0[f]->GetPhys() +offset);
632 if ((boost::math::isnan)(value))
634 ASSERTL0(
false,
"new value is not a number");
638 value = (value > clamp_up)? clamp_up :
639 ((value < clamp_low)? clamp_low :
642 m_f->m_data[f][r] = value;
648 for (f = 0; f < field0.size(); ++f)
650 m_f->m_data[f][r] = def_value;
654 if (intpts%1000 == 0&&isRoot)
662 m_f->m_data[nfields-2][r] = qinv*(
m_f->m_data[pfield][r] - p0);
667 for(
int i = 0; i < velid.size(); ++i)
669 q += 0.5*
m_f->m_data[velid[i]][r]*
m_f->m_data[velid[i]][r];
672 m_f->m_data[nfields-1][r] = qinv*(
m_f->m_data[pfield][r]+q - p0);
#define ASSERTL0(condition, msg)
pair< ModuleType, string > ModuleKey
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.
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 InterpolateFieldToPts(vector< MultiRegions::ExpListSharedPtr > &field0, Array< OneD, Array< OneD, NekDouble > > &pts, NekDouble clamp_low, NekDouble clamp_up, NekDouble def_value, bool isRoot)
map< string, ConfigOption > m_config
List of configuration values.
FieldSharedPtr m_f
Field object.
boost::shared_ptr< PtsField > PtsFieldSharedPtr
static SessionReaderSharedPtr CreateInstance(int argc, char *argv[])
Creates an instance of the SessionReader class.
#define WARNINGL0(condition, msg)
boost::shared_ptr< DomainRange > DomainRangeShPtr
boost::shared_ptr< Field > FieldSharedPtr
Represents a command-line configuration option.
static PtsFieldSharedPtr NullPtsField
virtual ~ProcessInterpPoints()
static bool GenerateUnOrderedVector(const char *const str, std::vector< NekDouble > &vec)
static FieldMetaDataMap NullFieldMetaDataMap
std::map< int, ExpansionShPtr > ExpansionMap
ModuleFactory & GetModuleFactory()
Abstract base class for processing modules.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.