Nektar++
Core/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
41
42#include <boost/format.hpp>
43
44using namespace std;
45
46namespace Nektar::SolverUtils
47{
48
49/**
50 * Representation of a FUNCTION defined in the session xml file.
51 *
52 * @param session The session where the function was defined.
53 * @param field The field the function is defined on.
54 * @param functionName The name of the function.
55 * @param toCache Store the evaluated function for later use.
56 */
59 const MultiRegions::ExpListSharedPtr &field, std::string functionName,
60 bool toCache)
61 : m_session(session), m_field(field), m_name(functionName),
62 m_toCache(toCache)
63{
64 ASSERTL0(m_session->DefinesFunction(m_name),
65 "Function '" + m_name + "' does not exist.");
66}
67
68/**
69 * Evaluates a function defined in the xml session file at each quadrature
70 * point.
71 *
72 * @param pArray The array into which to write the values.
73 * @param pTime The time at which to evaluate the function.
74 * @param domain The domain to evaluate the function in.
75 */
77 const NekDouble pTime, const int domain)
78{
79 std::vector<std::string> vFieldNames = m_session->GetVariables();
80
81 for (int i = 0; i < vFieldNames.size(); i++)
82 {
83 Evaluate(vFieldNames[i], pArray[i], pTime, domain);
84 }
85}
86
87/**
88 * Evaluates a function defined in the xml session file at each quadrature
89 * point.
90 *
91 * @param pFieldNames The names of the fields to evaluate the function for.
92 * @param pArray The array into which to write the values.
93 * @param pTime The time at which to evaluate the function.
94 * @param domain The domain to evaluate the function in.
95 */
96void SessionFunction::Evaluate(std::vector<std::string> pFieldNames,
98 const NekDouble &pTime, const int domain)
99{
100 ASSERTL1(pFieldNames.size() == pArray.size(),
101 "Function '" + m_name +
102 "' variable list size mismatch with array storage.");
103
104 for (int i = 0; i < pFieldNames.size(); i++)
105 {
106 Evaluate(pFieldNames[i], pArray[i], pTime, domain);
107 }
108}
109
110/**
111 * Evaluates a function defined in the xml session file at each quadrature
112 * point.
113 *
114 * @param pFieldNames The names of the fields to evaluate the function for.
115 * @param pFields The fields into which to write the values.
116 * @param pTime The time at which to evaluate the function.
117 * @param domain The domain to evaluate the function in.
118 */
120 std::vector<std::string> pFieldNames,
122 const NekDouble &pTime, const int domain)
123{
124 ASSERTL0(pFieldNames.size() == pFields.size(),
125 "Field list / name list size mismatch.");
126
127 for (int i = 0; i < pFieldNames.size(); i++)
128 {
129 Evaluate(pFieldNames[i], pFields[i]->UpdatePhys(), pTime, domain);
130 if (pFields[i]->GetWaveSpace())
131 {
132 pFields[i]->HomogeneousFwdTrans(pFields[i]->GetTotPoints(),
133 pFields[i]->GetPhys(),
134 pFields[i]->UpdatePhys());
135 }
136 pFields[i]->FwdTransLocalElmt(pFields[i]->GetPhys(),
137 pFields[i]->UpdateCoeffs());
138 }
139}
140
141/**
142 * Evaluates a function defined in the xml session file at each quadrature
143 * point.
144 *
145 * @param pFieldName The name of the field to evaluate the function for.
146 * @param pArray The array into which to write the values.
147 * @param pTime The time at which to evaluate the function.
148 * @param domain The domain to evaluate the function in.
149 */
150void SessionFunction::Evaluate(std::string pFieldName,
152 const NekDouble &pTime, const int domain)
153{
155 m_session->GetFunctionType(m_name, pFieldName, domain);
156
157 unsigned int nq = m_field->GetNpoints();
158
159 std::pair<std::string, int> key(pFieldName, domain);
160 // sorry
161 if ((m_arrays.find(key) != m_arrays.end()) &&
163 ((m_lastCached.find(key) != m_lastCached.end()) &&
164 (pTime - m_lastCached[key] < NekConstants::kNekZeroTol))))
165 {
166 // found cached field
167 if (pArray.size() < nq)
168 {
169 pArray = Array<OneD, NekDouble>(nq);
170 }
171 Vmath::Vcopy(nq, m_arrays[key], 1, pArray, 1);
172
173 return;
174 }
175
177 {
178 EvaluateExp(pFieldName, pArray, pTime, domain);
179 }
180 else if (vType == LibUtilities::eFunctionTypeFile ||
182 {
183 std::string filename =
184 m_session->GetFunctionFilename(m_name, pFieldName, domain);
185
186 if (fs::path(filename).extension() == ".pts" ||
187 fs::path(filename).extension() == ".csv")
188 {
189 EvaluatePts(pFieldName, pArray, pTime, domain);
190 }
191 else
192 {
193 EvaluateFld(pFieldName, pArray, pTime, domain);
194 }
195 }
196 else
197 {
198 ASSERTL0(false, "unknown eFunctionType");
199 }
200
201 if (m_toCache)
202 {
204 Vmath::Vcopy(nq, pArray, 1, m_arrays[key], 1);
205
206 m_lastCached[key] = pTime;
207 }
208}
209
210/**
211 * @brief Provide a description of a function for a given field name.
212 *
213 * @param pFieldName Field name.
214 * @param domain The domain to evaluate the function in.
215 */
216std::string SessionFunction::Describe(std::string pFieldName, const int domain)
217{
218 std::string retVal;
219
221 m_session->GetFunctionType(m_name, pFieldName, domain);
223 {
225 m_session->GetFunction(m_name, pFieldName, domain);
226 retVal = ffunc->GetExpression();
227 }
228 else if (vType == LibUtilities::eFunctionTypeFile ||
230 {
231 std::string filename =
232 m_session->GetFunctionFilename(m_name, pFieldName, domain);
233 retVal = "from file " + filename;
234 }
235 else
236 {
237 ASSERTL0(false, "unknown eFunctionType");
238 }
239
240 return retVal;
241}
242
243/**
244 * Evaluates a function from expression
245 *
246 * @param pFieldName The name of the field to evaluate the function for.
247 * @param pArray The array into which to write the values.
248 * @param pTime The time at which to evaluate the function.
249 * @param domain The domain to evaluate the function in.
250 */
251void SessionFunction::EvaluateExp(string pFieldName,
253 const NekDouble &pTime, const int domain)
254{
255 unsigned int nq = m_field->GetNpoints();
256 if (pArray.size() < nq)
257 {
258 pArray = Array<OneD, NekDouble>(nq);
259 }
260
264
265 // Get the coordinates (assuming all fields have the same
266 // discretisation)
267 m_field->GetCoords(x0, x1, x2);
269 m_session->GetFunction(m_name, pFieldName, domain);
270
271 ffunc->Evaluate(x0, x1, x2, pTime, pArray);
272}
273
274/**
275 * Evaluates a function from fld file
276 *
277 * @param pFieldName The name of the field to evaluate the function for.
278 * @param pArray The array into which to write the values.
279 * @param pTime The time at which to evaluate the function.
280 * @param domain The domain to evaluate the function in.
281 */
282void SessionFunction::EvaluateFld(string pFieldName,
284 const NekDouble &pTime, const int domain)
285{
286 unsigned int nq = m_field->GetNpoints();
287 if (pArray.size() < nq)
288 {
289 pArray = Array<OneD, NekDouble>(nq);
290 }
291
292 std::string filename =
293 m_session->GetFunctionFilename(m_name, pFieldName, domain);
294 std::string fileVar =
295 m_session->GetFunctionFilenameVariable(m_name, pFieldName, domain);
296
297 if (fileVar.length() == 0)
298 {
299 fileVar = pFieldName;
300 }
301
302 // In case of eFunctionTypeTransientFile, generate filename from
303 // format string
305 m_session->GetFunctionType(m_name, pFieldName, domain);
307 {
308 try
309 {
310#if (defined _WIN32 && _MSC_VER < 1900)
311 // We need this to make sure boost::format has always
312 // two digits in the exponents of Scientific notation.
313 unsigned int old_exponent_format;
314 old_exponent_format = _set_output_format(_TWO_DIGIT_EXPONENT);
315 filename = boost::str(boost::format(filename) % pTime);
316 _set_output_format(old_exponent_format);
317#else
318 filename = boost::str(boost::format(filename) % pTime);
319#endif
320 }
321 catch (...)
322 {
323 ASSERTL0(false, "Invalid Filename in function \"" + m_name +
324 "\", variable \"" + fileVar + "\"")
325 }
326 }
327
328 // Define list of global element ids
329 int numexp = m_field->GetExpSize();
330 Array<OneD, int> ElementGIDs(numexp);
331 for (int i = 0; i < numexp; ++i)
332 {
333 ElementGIDs[i] = m_field->GetExp(i)->GetGeom()->GetGlobalID();
334 }
335
336 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
337 std::vector<std::vector<NekDouble>> FieldData;
340 fldIO->Import(filename, FieldDef, FieldData,
342
343 int idx = -1;
344 Array<OneD, NekDouble> vCoeffs(m_field->GetNcoeffs(), 0.0);
345 // Loop over all the expansions
346 for (int i = 0; i < FieldDef.size(); ++i)
347 {
348 // Find the index of the required field in the
349 // expansion segment
350 for (int j = 0; j < FieldDef[i]->m_fields.size(); ++j)
351 {
352 if (FieldDef[i]->m_fields[j] == fileVar)
353 {
354 idx = j;
355 }
356 }
357
358 if (idx >= 0)
359 {
360 m_field->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
361 FieldDef[i]->m_fields[idx], vCoeffs);
362 }
363 else
364 {
365 cout << "Field " + fileVar + " not found." << endl;
366 }
367 }
368
369 bool wavespace = m_field->GetWaveSpace();
370 m_field->SetWaveSpace(false);
371 m_field->BwdTrans(vCoeffs, pArray);
372 m_field->SetWaveSpace(wavespace);
373}
374
375/**
376 * Evaluates a function from pts file
377 *
378 * @param pFieldName The name of the field to evaluate the function for.
379 * @param pArray The array into which to write the values.
380 * @param pTime The time at which to evaluate the function.
381 * @param domain The domain to evaluate the function in.
382 */
383void SessionFunction::EvaluatePts(string pFieldName,
385 const NekDouble &pTime, const int domain)
386{
387 unsigned int nq = m_field->GetNpoints();
388 if (pArray.size() < nq)
389 {
390 pArray = Array<OneD, NekDouble>(nq);
391 }
392
393 std::string filename =
394 m_session->GetFunctionFilename(m_name, pFieldName, domain);
395 std::string fileVar =
396 m_session->GetFunctionFilenameVariable(m_name, pFieldName, domain);
397
398 if (fileVar.length() == 0)
399 {
400 fileVar = pFieldName;
401 }
402
403 // In case of eFunctionTypeTransientFile, generate filename from
404 // format string
406 m_session->GetFunctionType(m_name, pFieldName, domain);
408 {
409 try
410 {
411#if (defined _WIN32 && _MSC_VER < 1900)
412 // We need this to make sure boost::format has always
413 // two digits in the exponents of Scientific notation.
414 unsigned int old_exponent_format;
415 old_exponent_format = _set_output_format(_TWO_DIGIT_EXPONENT);
416 filename = boost::str(boost::format(filename) % pTime);
417 _set_output_format(old_exponent_format);
418#else
419 filename = boost::str(boost::format(filename) % pTime);
420#endif
421 }
422 catch (...)
423 {
424 ASSERTL0(false, "Invalid Filename in function \"" + m_name +
425 "\", variable \"" + fileVar + "\"")
426 }
427 }
428
430 // check if we already loaded this file. For transient files,
431 // funcFilename != filename so we can make sure we only keep the
432 // latest pts field per funcFilename.
433 std::string funcFilename =
434 m_session->GetFunctionFilename(m_name, pFieldName, domain);
435
437 if (fs::path(filename).extension() == ".pts")
438 {
439 LibUtilities::PtsIO ptsIO(m_session->GetComm());
440 ptsIO.Import(filename, inPts);
441 }
442 else if (fs::path(filename).extension() == ".csv")
443 {
444 LibUtilities::CsvIO csvIO(m_session->GetComm());
445 csvIO.Import(filename, inPts);
446 }
447 else
448 {
449 ASSERTL1(false, "Unsupported file type");
450 }
451
452 Array<OneD, Array<OneD, NekDouble>> pts(inPts->GetDim() +
453 inPts->GetNFields());
454 for (int i = 0; i < inPts->GetDim() + inPts->GetNFields(); ++i)
455 {
456 pts[i] = Array<OneD, NekDouble>(nq);
457 }
458 if (inPts->GetDim() == 1)
459 {
460 m_field->GetCoords(pts[0]);
461 }
462 else if (inPts->GetDim() == 2)
463 {
464 m_field->GetCoords(pts[0], pts[1]);
465 }
466 else if (inPts->GetDim() == 3)
467 {
468 m_field->GetCoords(pts[0], pts[1], pts[2]);
469 }
471 inPts->GetDim(), inPts->GetFieldNames(), pts);
472
474 interp;
475 if (m_interpolators.find(funcFilename) != m_interpolators.end())
476 {
477 interp = m_interpolators[funcFilename];
478 }
479 else
480 {
482 std::vector<MultiRegions::ExpListSharedPtr>>(
484 if (m_session->GetComm()->GetRank() == 0)
485 {
487 this);
488 }
489 interp.CalcWeights(inPts, outPts);
490 if (m_session->GetComm()->GetRank() == 0)
491 {
492 cout << endl;
493 if (m_session->DefinesCmdLineArgument("verbose"))
494 {
495 interp.PrintStatistics();
496 }
497 }
498 }
499
500 if (m_toCache)
501 {
502 m_interpolators[funcFilename] = interp;
503 }
504
505 // TODO: only interpolate the field we actually want
506 interp.Interpolate(inPts, outPts);
507
508 int fieldInd;
509 vector<string> fieldNames = outPts->GetFieldNames();
510 for (fieldInd = 0; fieldInd < fieldNames.size(); ++fieldInd)
511 {
512 if (outPts->GetFieldName(fieldInd) == fileVar)
513 {
514 break;
515 }
516 }
517 ASSERTL0(fieldInd != fieldNames.size(), "field not found");
518
519 pArray = outPts->GetPts(fieldInd + outPts->GetDim());
520}
521
522// end of namespaces
523} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
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., NekDouble tolerance=NekConstants::kFindDistanceMin)
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:224
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:66
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 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 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 SessionFunction(const LibUtilities::SessionReaderSharedPtr &session, const MultiRegions::ExpListSharedPtr &field, std::string functionName, bool toCache=false)
Representation of a FUNCTION defined in the session xml file.
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:322
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:125
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:51
std::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:184
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
static const NekDouble kNekZeroTol
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
STL namespace.