39 #include <boost/geometry.hpp>
46 #include <boost/math/special_functions/fpclassify.hpp>
48 namespace bg = boost::geometry;
49 namespace bgi = boost::geometry::index;
59 ProcessInterpField::create,
60 "Interpolates one field to another, requires fromxml, "
61 "fromfld to be defined");
67 false,
"NotSet",
"Xml file form which to interpolate field");
69 false,
"NotSet",
"Fld file form which to interpolate field");
72 ConfigOption(
false,
"-10000000",
"Lower bound for interpolation value");
74 ConfigOption(
false,
"10000000",
"Upper bound for interpolation value");
76 ConfigOption(
false,
"0",
"Default value if point is outside domain");
87 if (
m_f->m_comm->TreatAsRankZero())
89 cout <<
"ProcessInterpField: Interpolating field..." << endl;
95 std::vector<std::string> files;
106 int coordim =
m_f->m_exp[0]->GetCoordim(0);
107 int npts =
m_f->m_exp[0]->GetTotPoints();
110 for (
int i = 0; i < coordim; ++i)
115 for (
int i = coordim; i < 3; ++i)
120 m_f->m_exp[0]->GetCoords(coords[0], coords[1], coords[2]);
122 rng->m_checkShape =
false;
126 rng->m_doZrange =
true;
130 rng->m_doYrange =
true;
134 rng->m_doXrange =
true;
139 ASSERTL0(
false,
"too many values specfied in range");
152 if (!expansions.size())
158 SpatialDomains::ExpansionMap::const_iterator expIt;
161 for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
163 ElementGIDs[i++] = expIt->second->m_geomShPtr->GetGlobalID();
166 string fromfld =
m_config[
"fromfld"].as<
string>();
167 m_f->FieldIOForFile(fromfld)->Import(
171 int NumHomogeneousDir =
m_fromField->m_fielddef[0]->m_numHomogeneousDir;
177 int nfields =
m_fromField->m_fielddef[0]->m_fields.size();
181 m_fromField->SetUpFirstExpList(NumHomogeneousDir,
true);
183 m_f->m_exp.resize(nfields);
186 for (i = 1; i < nfields; ++i)
188 m_f->m_exp[i] =
m_f->AppendExpList(NumHomogeneousDir);
193 for (
int j = 0; j < nfields; ++j)
195 for (i = 0; i <
m_fromField->m_fielddef.size(); i++)
206 int nq1 =
m_f->m_exp[0]->GetTotPoints();
214 m_f->m_exp[0]->GetCoords(x1, y1);
216 else if (coordim == 3)
218 m_f->m_exp[0]->GetCoords(x1, y1, z1);
225 for (
int i = 0; i < nfields; i++)
227 for (
int j = 0; j < nq1; ++j)
229 m_f->m_exp[i]->UpdatePhys()[j] = def_value;
234 if (
m_f->m_comm->GetRank() == 0)
239 if (
m_f->m_comm->GetRank() == 0)
244 for (
int i = 0; i < nfields; ++i)
246 for (
int j = 0; j < nq1; ++j)
248 if (
m_f->m_exp[i]->GetPhys()[j] > clamp_up)
250 m_f->m_exp[i]->UpdatePhys()[j] = clamp_up;
252 else if (
m_f->m_exp[i]->GetPhys()[j] < clamp_low)
254 m_f->m_exp[i]->UpdatePhys()[j] = clamp_low;
260 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef =
261 m_f->m_exp[0]->GetFieldDefinitions();
262 std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
264 for (
int j = 0; j < nfields; ++j)
266 m_f->m_exp[j]->FwdTrans(
m_f->m_exp[j]->GetPhys(),
267 m_f->m_exp[j]->UpdateCoeffs());
268 for (i = 0; i < FieldDef.size(); ++i)
270 FieldDef[i]->m_fields.push_back(
272 m_f->m_exp[j]->AppendFieldData(FieldDef[i], FieldData[i]);
276 m_f->m_fielddef = FieldDef;
277 m_f->m_data = FieldData;
281 const int goal)
const
map< string, ConfigOption > m_config
List of configuration values.
#define ASSERTL0(condition, msg)
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 ...
static Array< OneD, NekDouble > NullNekDouble1DArray
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
void PrintProgressbar(const int position, const int goal) const
static SessionReaderSharedPtr CreateInstance(int argc, char *argv[])
Creates an instance of the SessionReader class.
boost::shared_ptr< Field > FieldSharedPtr
boost::shared_ptr< DomainRange > DomainRangeShPtr
FieldSharedPtr m_fromField
virtual ~ProcessInterpField()
virtual void Process(po::variables_map &vm)
Write mesh to output file.
Abstract base class for processing modules.
static FieldMetaDataMap NullFieldMetaDataMap
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.