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 or plane 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");
91 "Specify a plane of N1 x N1 points using "
92 "plane=N1,N2,x0,y0,z0,z1,y1,z1,x2,y2,z2,x3,"
105 cout <<
"Processing point interpolation" << endl;
112 if(
m_config[
"line"].as<string>().compare(
"NotSet") != 0)
114 string help =
m_config[
"line"].as<
string>();
115 vector<NekDouble> values;
117 m_config[
"line"].as<string>().c_str(),values),
118 "Failed to interpret line string");
121 "line string should contain 2Dim+1 values "
122 "N,x0,y0,z0,x1,y1,z1");
124 int dim = (values.size()-1)/2;
125 int npts = values[0];
128 for(
int i = 0; i < dim; ++i)
133 for(
int i = 0; i <
npts; ++i)
135 pts[0][i] = values[1]
136 + i/((
NekDouble)(npts-1))*(values[dim+1] - values[1]);
139 pts[1][i] = values[2]
140 + i/((
NekDouble)(npts-1))*(values[dim+2] - values[2]);
144 pts[2][i] = values[3]
145 + i/((
NekDouble)(npts-1))*(values[dim+3] - values[3]);
153 m_f->m_fieldPts->SetPointsPerEdge(ppe);
157 else if(
m_config[
"plane"].as<string>().compare(
"NotSet") != 0)
159 string help =
m_config[
"plane"].as<
string>();
160 vector<NekDouble> values;
162 m_config[
"plane"].as<string>().c_str(),values),
163 "Failed to interpret plane string");
166 "line string should contain 2Dim+1 values "
167 "N,x0,y0,z0,x1,y1,z1");
170 int dim = (values.size()-1)/4;
172 int npts1 = values[0];
173 int npts2 = values[1];
177 for(
int i = 0; i < dim; ++i)
182 for(
int j = 0; j < npts2; ++j)
184 for(
int i = 0; i < npts1; ++i)
187 (values[2] + i/((
NekDouble)(npts1-1))*(values[dim+2] - values[2]))*(1.0-j/((
NekDouble)(npts2-1))) +
188 (values[3*dim+2] + i/((
NekDouble)(npts1-1))*(values[2*dim+2] - values[3*dim+2]))*(j/((
NekDouble)(npts2-1)));
192 (values[3] + i/((
NekDouble)(npts1-1))*(values[dim+3] - values[3]))*(1.0-j/((
NekDouble)(npts2-1))) +
193 (values[3*dim+3] + i/((
NekDouble)(npts1-1))*(values[2*dim+3] - values[3*dim+3]))*(j/((
NekDouble)(npts2-1)));
198 (values[4] + i/((
NekDouble)(npts1-1))*(values[dim+4] - values[4]))*(1.0-j/((
NekDouble)(npts2-1))) +
199 (values[3*dim+4] + i/((
NekDouble)(npts1-1))*(values[2*dim+4] - values[3*dim+4]))*(j/((
NekDouble)(npts2-1)));
206 ppe.push_back(npts1);
207 ppe.push_back(npts2);
209 m_f->m_fieldPts->SetPointsPerEdge(ppe);
218 std::vector<std::string> files;
220 files.push_back(
m_config[
"fromxml"].as<string>());
227 int coordim =
m_f->m_fieldPts->GetDim();
228 int npts =
m_f->m_fieldPts->GetNpoints();
230 m_f->m_fieldPts->GetPts(pts);
232 rng->m_checkShape =
false;
236 rng->m_doZrange =
true;
240 rng->m_doYrange =
true;
244 rng->m_doXrange =
true;
249 ASSERTL0(
false,
"too many values specfied in range");
262 SpatialDomains::ExpansionMap::const_iterator expIt;
265 for (expIt = expansions.begin(); expIt != expansions.end();
268 ElementGIDs[i++] = expIt->second->m_geomShPtr->GetGlobalID();
271 string fromfld =
m_config[
"fromfld"].as<
string>();
272 fromField->m_fld->Import(fromfld,fromField->m_fielddef,
277 int NumHomogeneousDir = fromField->m_fielddef[0]->m_numHomogeneousDir;
281 fromField->m_graph->SetExpansions(fromField->m_fielddef);
283 int nfields = fromField->m_fielddef[0]->m_fields.size();
285 fromField->m_exp.resize(nfields);
286 fromField->m_exp[0] = fromField->SetUpFirstExpList(NumHomogeneousDir,
true);
288 m_f->m_exp.resize(nfields);
291 for(i = 1; i < nfields; ++i)
293 fromField->m_exp[i] = fromField->AppendExpList(NumHomogeneousDir);
297 for(
int j = 0; j < nfields; ++j)
299 for (i = 0; i < fromField->m_fielddef.size(); i++)
301 fromField->m_exp[j]->ExtractDataToCoeffs(
302 fromField->m_fielddef[i],
303 fromField->m_data[i],
304 fromField->m_fielddef[0]->m_fields[j],
305 fromField->m_exp[j]->UpdateCoeffs());
307 fromField->m_exp[j]->BwdTrans(fromField->m_exp[j]->GetCoeffs(),
308 fromField->m_exp[j]->UpdatePhys());
311 m_f->m_fieldPts->AddField(newPts, fromField->m_fielddef[0]->m_fields[j]);
314 if(fromField->m_session->GetComm()->GetRank() == 0)
316 cout <<
"Interpolating [" << flush;
324 clamp_low, clamp_up, def_value);
326 if(fromField->m_session->GetComm()->GetRank() == 0)
334 vector<MultiRegions::ExpListSharedPtr> &field0,
340 int expdim = field0[0]->GetCoordim(0);
343 int nq1 = pts[0].num_elements();
349 m_f->m_data.resize(field0.size());
351 for (f = 0; f < field0.size(); ++f)
353 m_f->m_data[f].resize(nq1);
356 for (r = 0; r < nq1; r++)
358 coords[0] = pts[0][r];
359 coords[1] = pts[1][r];
362 coords[2] = pts[2][r];
366 elmtid = field0[0]->GetExpIndex(coords, Lcoords, 1e-3);
370 offset = field0[0]->GetPhys_Offset(field0[0]->
371 GetOffset_Elmt_Id(elmtid));
373 for (f = 0; f < field0.size(); ++f)
376 value = field0[f]->GetExp(elmtid)->
377 StdPhysEvaluate(Lcoords, field0[f]->GetPhys() +offset);
379 if ((boost::math::isnan)(value))
381 ASSERTL0(
false,
"new value is not a number");
385 value = (value > clamp_up)? clamp_up :
386 ((value < clamp_low)? clamp_low :
389 m_f->m_data[f][r] = value;
395 for (f = 0; f < field0.size(); ++f)
397 m_f->m_data[f][r] = def_value;
401 if (intpts%1000 == 0)
#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.
map< string, ConfigOption > m_config
List of configuration values.
FieldSharedPtr m_f
Field object.
static SessionReaderSharedPtr CreateInstance(int argc, char *argv[])
Creates an instance of the SessionReader class.
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.
void InterpolateFieldToPts(vector< MultiRegions::ExpListSharedPtr > &field0, Array< OneD, Array< OneD, NekDouble > > &pts, NekDouble clamp_low, NekDouble clamp_up, NekDouble def_value)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.