50 "Given discrete data at quadrature points project them onto an "
52 "basis and output fld file. Requires frompts and .xml and .fld files.");
58 false,
"NotSet",
"reset any nan value to prescribed value");
61 false,
"NotSet",
"Pts file from which to interpolate field");
73 bool setnantovalue =
false;
76 if (!boost::iequals(
m_config[
"setnantovalue"].as<string>(),
"NotSet"))
86 "ProcessInterpPointDataToFld requires frompts parameter");
87 string inFile =
m_config[
"frompts"].as<
string>().c_str();
92 if (fs::path(inFile).extension() ==
".pts")
96 else if (fs::path(inFile).extension() ==
".csv")
103 "Unsupported file format for the \"frompts\" parameter. "
104 "Supported formats: \".pts\" and \".csv\"");
107 int nFields = fieldPts->GetNFields();
108 ASSERTL0(nFields > 0,
"No field values provided in input");
110 int dim = fieldPts->GetDim();
114 m_f->m_numHomogeneousDir == 0,
115 "ProcessInterpPointDataToFld does not support homogeneous expansion");
117 m_f->m_exp.resize(nFields);
118 for (i = 1; i < nFields; ++i)
120 m_f->m_exp[i] =
m_f->AppendExpList(
m_f->m_numHomogeneousDir);
124 fieldPts->GetPts(pts);
129 for (
int i = 0; i < pts[0].size(); ++i)
131 for (
int j = 0; j < nFields; ++j)
133 if ((boost::math::isnan)(pts[j + dim][i]))
135 pts[j + dim][i] = defvalue;
143 int totcoeffs =
m_f->m_exp[0]->GetNcoeffs();
145 ASSERTL0(pts[0].size() != totcoeffs,
146 "length of points in .pts file is different "
147 "to the number of coefficients in expansion ");
149 for (
int i = 0; i < nFields; ++i)
153 for (
int e = 0; e <
m_f->m_exp[0]->GetNumElmts(); ++e)
155 int ncoeffs =
m_f->m_exp[i]->GetExp(e)->GetNcoeffs();
156 int offset =
m_f->m_exp[i]->GetCoeff_Offset(e);
157 Vmath::Vcopy(ncoeffs, &pts[i + dim][cnt], 1, &coeffs[offset],
160 m_f->m_exp[i]->GetExp(e)->EquiSpacedToCoeffs(
161 coeffs + offset, tmp = coeffs + offset);
164 m_f->m_exp[i]->BwdTrans(
m_f->m_exp[i]->GetCoeffs(),
165 m_f->m_exp[i]->UpdatePhys());
170 int totpoints =
m_f->m_exp[0]->GetTotPoints();
177 m_f->m_exp[0]->GetCoords(coords[0], coords[1], coords[2]);
179 if (pts[0].size() != totpoints)
181 WARNINGL0(
false,
"length of points in .pts file is different to "
182 "the number of quadrature points in xml file");
183 totpoints = min(totpoints, (
int)pts[0].size());
186 for (j = 0; j < dim; ++j)
188 for (i = 0; i < totpoints; ++i)
190 if (fabs(coords[j][i] - pts[j][i]) > 1e-4)
193 "Coordinates do not match within 1e-4: " +
194 boost::lexical_cast<string>(coords[j][i]) +
" versus " +
195 boost::lexical_cast<string>(pts[j][i]) +
" diff " +
196 boost::lexical_cast<string>(
197 fabs(coords[j][i] - pts[j][i]));
204 for (j = 0; j < nFields; ++j)
208 for (i = 0; i < totpoints; ++i)
210 phys[i] = pts[j + dim][i];
214 m_f->m_exp[j]->FwdTransLocalElmt(phys,
215 m_f->m_exp[j]->UpdateCoeffs());
220 for (
int j = 0; j < fieldPts->GetNFields(); ++j)
222 m_f->m_variables.push_back(fieldPts->GetFieldName(j));
#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)
FieldSharedPtr m_f
Field object.
std::map< std::string, ConfigOption > m_config
List of configuration values.
Abstract base class for processing modules.
~ProcessPointDataToFld() override
static ModuleKey className
void v_Process(po::variables_map &vm) override
Write mesh to output file.
static std::shared_ptr< Module > create(FieldSharedPtr f)
Creates an instance of this class.
ProcessPointDataToFld(FieldSharedPtr f)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
void Import(const std::string &inFile, PtsFieldSharedPtr &ptsField, FieldMetaDataMap &fieldmetadatamap=NullFieldMetaDataMap, DomainRangeShPtr &Range=NullDomainRangeShPtr)
Import a pts field from file.
std::shared_ptr< Field > FieldSharedPtr
std::pair< ModuleType, std::string > ModuleKey
ModuleFactory & GetModuleFactory()
std::shared_ptr< PtsField > PtsFieldSharedPtr
CommFactory & GetCommFactory()
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Represents a command-line configuration option.