Nektar++
SessionFunction.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: SessionFunction.cpp
4 //
5 // For more information, please see: http://www.nektar.info/
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2017 Kilian Lackhove
10 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
11 // Department of Aeronautics, Imperial College London (UK), and Scientific
12 // Computing and Imaging Institute, University of Utah (USA).
13 //
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Session Function
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 
38 
40 
41 #include <boost/format.hpp>
42 
43 using namespace std;
44 
45 namespace Nektar
46 {
47 namespace SolverUtils
48 {
49 
50 /**
51  * Representation of a FUNCTION defined in the session xml file.
52  *
53  * @param session The session where the function was defined.
54  * @param field The field the function is defined on.
55  * @param functionName The name of the function.
56  * @param toCache Store the evaluated function for later use.
57  */
58 SessionFunction::SessionFunction(
60  const MultiRegions::ExpListSharedPtr &field, std::string functionName,
61  bool toCache)
62  : m_session(session), m_field(field), m_name(functionName),
63  m_toCache(toCache)
64 {
65  ASSERTL0(m_session->DefinesFunction(m_name),
66  "Function '" + m_name + "' does not exist.");
67 }
68 
69 /**
70  * Evaluates a function defined in the xml session file at each quadrature
71  * point.
72  *
73  * @param pArray The array into which to write the values.
74  * @param pTime The time at which to evaluate the function.
75  * @param domain The domain to evaluate the function in.
76  */
78  const NekDouble pTime, const int domain)
79 {
80  std::vector<std::string> vFieldNames = m_session->GetVariables();
81 
82  for (int i = 0; i < vFieldNames.size(); i++)
83  {
84  Evaluate(vFieldNames[i], pArray[i], pTime, domain);
85  }
86 }
87 
88 /**
89  * Evaluates a function defined in the xml session file at each quadrature
90  * point.
91  *
92  * @param pFieldNames The names of the fields to evaluate the function for.
93  * @param pArray The array into which to write the values.
94  * @param pTime The time at which to evaluate the function.
95  * @param domain The domain to evaluate the function in.
96  */
97 void SessionFunction::Evaluate(std::vector<std::string> pFieldNames,
99  const NekDouble &pTime, const int domain)
100 {
101  ASSERTL1(pFieldNames.size() == pArray.size(),
102  "Function '" + m_name +
103  "' variable list size mismatch with array storage.");
104 
105  for (int i = 0; i < pFieldNames.size(); i++)
106  {
107  Evaluate(pFieldNames[i], pArray[i], pTime, domain);
108  }
109 }
110 
111 /**
112  * Evaluates a function defined in the xml session file at each quadrature
113  * point.
114  *
115  * @param pFieldNames The names of the fields to evaluate the function for.
116  * @param pFields The fields into which to write the values.
117  * @param pTime The time at which to evaluate the function.
118  * @param domain The domain to evaluate the function in.
119  */
121  std::vector<std::string> pFieldNames,
123  const NekDouble &pTime, const int domain)
124 {
125  ASSERTL0(pFieldNames.size() == pFields.size(),
126  "Field list / name list size mismatch.");
127 
128  for (int i = 0; i < pFieldNames.size(); i++)
129  {
130  Evaluate(pFieldNames[i], pFields[i]->UpdatePhys(), pTime, domain);
131  if (pFields[i]->GetWaveSpace())
132  {
133  pFields[i]->HomogeneousFwdTrans(pFields[i]->GetTotPoints(),
134  pFields[i]->GetPhys(),
135  pFields[i]->UpdatePhys());
136  }
137  pFields[i]->FwdTransLocalElmt(pFields[i]->GetPhys(),
138  pFields[i]->UpdateCoeffs());
139  }
140 }
141 
142 /**
143  * Evaluates a function defined in the xml session file at each quadrature
144  * point.
145  *
146  * @param pFieldName The name of the field to evaluate the function for.
147  * @param pArray The array into which to write the values.
148  * @param pTime The time at which to evaluate the function.
149  * @param domain The domain to evaluate the function in.
150  */
151 void SessionFunction::Evaluate(std::string pFieldName,
152  Array<OneD, NekDouble> &pArray,
153  const NekDouble &pTime, const int domain)
154 {
156  m_session->GetFunctionType(m_name, pFieldName, domain);
157 
158  unsigned int nq = m_field->GetNpoints();
159 
160  std::pair<std::string, int> key(pFieldName, domain);
161  // sorry
162  if ((m_arrays.find(key) != m_arrays.end()) &&
164  ((m_lastCached.find(key) != m_lastCached.end()) &&
165  (pTime - m_lastCached[key] < NekConstants::kNekZeroTol))))
166  {
167  // found cached field
168  if (pArray.size() < nq)
169  {
170  pArray = Array<OneD, NekDouble>(nq);
171  }
172  Vmath::Vcopy(nq, m_arrays[key], 1, pArray, 1);
173 
174  return;
175  }
176 
178  {
179  EvaluateExp(pFieldName, pArray, pTime, domain);
180  }
181  else if (vType == LibUtilities::eFunctionTypeFile ||
183  {
184  std::string filename =
185  m_session->GetFunctionFilename(m_name, pFieldName, domain);
186 
187  if (boost::filesystem::path(filename).extension() == ".pts" ||
188  boost::filesystem::path(filename).extension() == ".csv")
189  {
190  EvaluatePts(pFieldName, pArray, pTime, domain);
191  }
192  else
193  {
194  EvaluateFld(pFieldName, pArray, pTime, domain);
195  }
196  }
197  else
198  {
199  ASSERTL0(false, "unknown eFunctionType");
200  }
201 
202  if (m_toCache)
203  {
204  m_arrays[key] = Array<OneD, NekDouble>(nq);
205  Vmath::Vcopy(nq, pArray, 1, m_arrays[key], 1);
206 
207  m_lastCached[key] = pTime;
208  }
209 }
210 
211 /**
212  * @brief Provide a description of a function for a given field name.
213  *
214  * @param pFieldName Field name.
215  * @param domain The domain to evaluate the function in.
216  */
217 std::string SessionFunction::Describe(std::string pFieldName, const int domain)
218 {
219  std::string retVal;
220 
222  m_session->GetFunctionType(m_name, pFieldName, domain);
224  {
226  m_session->GetFunction(m_name, pFieldName, domain);
227  retVal = ffunc->GetExpression();
228  }
229  else if (vType == LibUtilities::eFunctionTypeFile ||
231  {
232  std::string filename =
233  m_session->GetFunctionFilename(m_name, pFieldName, domain);
234  retVal = "from file " + filename;
235  }
236  else
237  {
238  ASSERTL0(false, "unknown eFunctionType");
239  }
240 
241  return retVal;
242 }
243 
244 /**
245  * Evaluates a function from expression
246  *
247  * @param pFieldName The name of the field to evaluate the function for.
248  * @param pArray The array into which to write the values.
249  * @param pTime The time at which to evaluate the function.
250  * @param domain The domain to evaluate the function in.
251  */
252 void SessionFunction::EvaluateExp(string pFieldName,
253  Array<OneD, NekDouble> &pArray,
254  const NekDouble &pTime, const int domain)
255 {
256  unsigned int nq = m_field->GetNpoints();
257  if (pArray.size() < nq)
258  {
259  pArray = Array<OneD, NekDouble>(nq);
260  }
261 
262  Array<OneD, NekDouble> x0(nq);
263  Array<OneD, NekDouble> x1(nq);
264  Array<OneD, NekDouble> x2(nq);
265 
266  // Get the coordinates (assuming all fields have the same
267  // discretisation)
268  m_field->GetCoords(x0, x1, x2);
270  m_session->GetFunction(m_name, pFieldName, domain);
271 
272  ffunc->Evaluate(x0, x1, x2, pTime, pArray);
273 }
274 
275 /**
276  * Evaluates a function from fld file
277  *
278  * @param pFieldName The name of the field to evaluate the function for.
279  * @param pArray The array into which to write the values.
280  * @param pTime The time at which to evaluate the function.
281  * @param domain The domain to evaluate the function in.
282  */
283 void SessionFunction::EvaluateFld(string pFieldName,
284  Array<OneD, NekDouble> &pArray,
285  const NekDouble &pTime, const int domain)
286 {
287  unsigned int nq = m_field->GetNpoints();
288  if (pArray.size() < nq)
289  {
290  pArray = Array<OneD, NekDouble>(nq);
291  }
292 
293  std::string filename =
294  m_session->GetFunctionFilename(m_name, pFieldName, domain);
295  std::string fileVar =
296  m_session->GetFunctionFilenameVariable(m_name, pFieldName, domain);
297 
298  if (fileVar.length() == 0)
299  {
300  fileVar = pFieldName;
301  }
302 
303  // In case of eFunctionTypeTransientFile, generate filename from
304  // format string
306  m_session->GetFunctionType(m_name, pFieldName, domain);
308  {
309  try
310  {
311 #if (defined _WIN32 && _MSC_VER < 1900)
312  // We need this to make sure boost::format has always
313  // two digits in the exponents of Scientific notation.
314  unsigned int old_exponent_format;
315  old_exponent_format = _set_output_format(_TWO_DIGIT_EXPONENT);
316  filename = boost::str(boost::format(filename) % pTime);
317  _set_output_format(old_exponent_format);
318 #else
319  filename = boost::str(boost::format(filename) % pTime);
320 #endif
321  }
322  catch (...)
323  {
324  ASSERTL0(false, "Invalid Filename in function \"" + m_name +
325  "\", variable \"" + fileVar + "\"")
326  }
327  }
328 
329  // Define list of global element ids
330  int numexp = m_field->GetExpSize();
331  Array<OneD, int> ElementGIDs(numexp);
332  for (int i = 0; i < numexp; ++i)
333  {
334  ElementGIDs[i] = m_field->GetExp(i)->GetGeom()->GetGlobalID();
335  }
336 
337  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
338  std::vector<std::vector<NekDouble>> FieldData;
341  fldIO->Import(filename, FieldDef, FieldData,
343 
344  int idx = -1;
345  Array<OneD, NekDouble> vCoeffs(m_field->GetNcoeffs(), 0.0);
346  // Loop over all the expansions
347  for (int i = 0; i < FieldDef.size(); ++i)
348  {
349  // Find the index of the required field in the
350  // expansion segment
351  for (int j = 0; j < FieldDef[i]->m_fields.size(); ++j)
352  {
353  if (FieldDef[i]->m_fields[j] == fileVar)
354  {
355  idx = j;
356  }
357  }
358 
359  if (idx >= 0)
360  {
361  m_field->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
362  FieldDef[i]->m_fields[idx], vCoeffs);
363  }
364  else
365  {
366  cout << "Field " + fileVar + " not found." << endl;
367  }
368  }
369 
370  bool wavespace = m_field->GetWaveSpace();
371  m_field->SetWaveSpace(false);
372  m_field->BwdTrans(vCoeffs, pArray);
373  m_field->SetWaveSpace(wavespace);
374 }
375 
376 /**
377  * Evaluates a function from pts file
378  *
379  * @param pFieldName The name of the field to evaluate the function for.
380  * @param pArray The array into which to write the values.
381  * @param pTime The time at which to evaluate the function.
382  * @param domain The domain to evaluate the function in.
383  */
384 void SessionFunction::EvaluatePts(string pFieldName,
385  Array<OneD, NekDouble> &pArray,
386  const NekDouble &pTime, const int domain)
387 {
388  unsigned int nq = m_field->GetNpoints();
389  if (pArray.size() < nq)
390  {
391  pArray = Array<OneD, NekDouble>(nq);
392  }
393 
394  std::string filename =
395  m_session->GetFunctionFilename(m_name, pFieldName, domain);
396  std::string fileVar =
397  m_session->GetFunctionFilenameVariable(m_name, pFieldName, domain);
398 
399  if (fileVar.length() == 0)
400  {
401  fileVar = pFieldName;
402  }
403 
404  // In case of eFunctionTypeTransientFile, generate filename from
405  // format string
407  m_session->GetFunctionType(m_name, pFieldName, domain);
409  {
410  try
411  {
412 #if (defined _WIN32 && _MSC_VER < 1900)
413  // We need this to make sure boost::format has always
414  // two digits in the exponents of Scientific notation.
415  unsigned int old_exponent_format;
416  old_exponent_format = _set_output_format(_TWO_DIGIT_EXPONENT);
417  filename = boost::str(boost::format(filename) % pTime);
418  _set_output_format(old_exponent_format);
419 #else
420  filename = boost::str(boost::format(filename) % pTime);
421 #endif
422  }
423  catch (...)
424  {
425  ASSERTL0(false, "Invalid Filename in function \"" + m_name +
426  "\", variable \"" + fileVar + "\"")
427  }
428  }
429 
431  // check if we already loaded this file. For transient files,
432  // funcFilename != filename so we can make sure we only keep the
433  // latest pts field per funcFilename.
434  std::string funcFilename =
435  m_session->GetFunctionFilename(m_name, pFieldName, domain);
436 
438  if (boost::filesystem::path(filename).extension() == ".pts")
439  {
440  LibUtilities::PtsIO ptsIO(m_session->GetComm());
441  ptsIO.Import(filename, inPts);
442  }
443  else if (boost::filesystem::path(filename).extension() == ".csv")
444  {
445  LibUtilities::CsvIO csvIO(m_session->GetComm());
446  csvIO.Import(filename, inPts);
447  }
448  else
449  {
450  ASSERTL1(false, "Unsupported file type");
451  }
452 
453  Array<OneD, Array<OneD, NekDouble>> pts(inPts->GetDim() +
454  inPts->GetNFields());
455  for (int i = 0; i < inPts->GetDim() + inPts->GetNFields(); ++i)
456  {
457  pts[i] = Array<OneD, NekDouble>(nq);
458  }
459  if (inPts->GetDim() == 1)
460  {
461  m_field->GetCoords(pts[0]);
462  }
463  else if (inPts->GetDim() == 2)
464  {
465  m_field->GetCoords(pts[0], pts[1]);
466  }
467  else if (inPts->GetDim() == 3)
468  {
469  m_field->GetCoords(pts[0], pts[1], pts[2]);
470  }
472  inPts->GetDim(), inPts->GetFieldNames(), pts);
473 
475  interp;
476  if (m_interpolators.find(funcFilename) != m_interpolators.end())
477  {
478  interp = m_interpolators[funcFilename];
479  }
480  else
481  {
482  interp = FieldUtils::Interpolator<
483  std::vector<MultiRegions::ExpListSharedPtr>>(
485  if (m_session->GetComm()->GetRank() == 0)
486  {
488  this);
489  }
490  interp.CalcWeights(inPts, outPts);
491  if (m_session->GetComm()->GetRank() == 0)
492  {
493  cout << endl;
494  if (m_session->DefinesCmdLineArgument("verbose"))
495  {
496  interp.PrintStatistics();
497  }
498  }
499  }
500 
501  if (m_toCache)
502  {
503  m_interpolators[funcFilename] = interp;
504  }
505 
506  // TODO: only interpolate the field we actually want
507  interp.Interpolate(inPts, outPts);
508 
509  int fieldInd;
510  vector<string> fieldNames = outPts->GetFieldNames();
511  for (fieldInd = 0; fieldInd < fieldNames.size(); ++fieldInd)
512  {
513  if (outPts->GetFieldName(fieldInd) == fileVar)
514  {
515  break;
516  }
517  }
518  ASSERTL0(fieldInd != fieldNames.size(), "field not found");
519 
520  pArray = outPts->GetPts(fieldInd + outPts->GetDim());
521 }
522 
523 // end of namespaces
524 } // namespace SolverUtils
525 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
A class that contains algorithms for interpolation between pts fields, expansions and different meshe...
FIELD_UTILS_EXPORT void Interpolate(const T expInField, T &expOutField, NekDouble def_value=0.0)
Interpolate from an expansion to an expansion.
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
Definition: FieldIO.cpp:226
void PrintStatistics()
Returns if the weights have already been computed.
void SetProgressCallback(FuncPointerT func, ObjectPointerT obj)
sets a callback funtion which gets called every time the interpolation progresses
void CalcWeights(const LibUtilities::PtsFieldSharedPtr ptsInField, LibUtilities::PtsFieldSharedPtr &ptsOutField, bool reuseTree=false)
Compute interpolation weights without doing any interpolation.
void Import(const std::string &inFile, PtsFieldSharedPtr &ptsField, FieldMetaDataMap &fieldmetadatamap=NullFieldMetaDataMap, DomainRangeShPtr &Range=NullDomainRangeShPtr)
Import a pts field from file.
Definition: PtsIO.cpp:68
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::map< std::string, FieldUtils::Interpolator< std::vector< MultiRegions::ExpListSharedPtr > > > m_interpolators
Interpolator for pts file input for a variable & domain combination.
MultiRegions::ExpListSharedPtr m_field
The expansion we want to evaluate this function for.
bool m_toCache
Store resulting arrays (and interpolators)
std::map< std::pair< std::string, int >, Array< OneD, NekDouble > > m_arrays
Cached result arrays.
SOLVER_UTILS_EXPORT void EvaluateFld(std::string pFieldName, Array< OneD, NekDouble > &pArray, const NekDouble &pTime=0.0, const int domain=0)
SOLVER_UTILS_EXPORT void EvaluatePts(std::string pFieldName, Array< OneD, NekDouble > &pArray, const NekDouble &pTime=0.0, const int domain=0)
Evaluates a function from pts file.
SOLVER_UTILS_EXPORT void Evaluate(Array< OneD, Array< OneD, NekDouble >> &pArray, const NekDouble pTime=0.0, const int domain=0)
Evaluates a function defined in the xml session file at each quadrature point.
SOLVER_UTILS_EXPORT void PrintProgressbar(const int position, const int goal) const
std::map< std::pair< std::string, int >, NekDouble > m_lastCached
Last time the cache for this variable & domain combo was updated.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT std::string Describe(std::string pFieldName, const int domain=0)
Provide a description of a function for a given field name.
SOLVER_UTILS_EXPORT void EvaluateExp(std::string pFieldName, Array< OneD, NekDouble > &pArray, const NekDouble &pTime=0.0, const int domain=0)
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:327
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:129
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:53
std::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:190
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
static const NekDouble kNekZeroTol
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255