Nektar++
Public Member Functions | Protected Member Functions | Protected Attributes | Friends | List of all members
Nektar::SolverUtils::FileFieldInterpolator Class Reference

#include <FileSolution.h>

Public Member Functions

void InitObject (const std::string functionName, LibUtilities::SessionReaderSharedPtr pSession, const Array< OneD, const MultiRegions::ExpListSharedPtr > pFields, std::set< std::string > &variables, std::map< std::string, int > &series, std::map< std::string, NekDouble > &time)
 
void InterpolateField (const std::string variable, Array< OneD, NekDouble > &outarray, NekDouble time)
 
void InterpolateField (const int v, Array< OneD, NekDouble > &outarray, const NekDouble time)
 
NekDouble GetStartTime ()
 

Protected Member Functions

DNekBlkMatSharedPtr GetFloquetBlockMatrix (int nexp)
 
 FileFieldInterpolator ()
 
 ~FileFieldInterpolator ()
 
void DFT (const std::string file, const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const bool timefromfile)
 
void ImportFldBase (std::string pInfile, const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, int slice, std::map< std::string, NekDouble > &params)
 Import Base flow. More...
 

Protected Attributes

LibUtilities::SessionReaderSharedPtr m_session
 
int m_start
 number of slices More...
 
int m_skip
 
int m_slices
 
NekDouble m_timeStart
 period length More...
 
NekDouble m_period
 
int m_isperiodic
 
int m_interporder
 
std::map< int, Array< OneD, NekDouble > > m_interp
 interpolation vector More...
 
std::map< std::string, int > m_variableMap
 variables More...
 

Friends

class MemoryManager< FileFieldInterpolator >
 

Detailed Description

Definition at line 52 of file FileSolution.h.

Constructor & Destructor Documentation

◆ FileFieldInterpolator()

Nektar::SolverUtils::FileFieldInterpolator::FileFieldInterpolator ( )
protected

Definition at line 451 of file FileSolution.cpp.

452{
453}

◆ ~FileFieldInterpolator()

Nektar::SolverUtils::FileFieldInterpolator::~FileFieldInterpolator ( )
protected

Definition at line 447 of file FileSolution.cpp.

448{
449}

Member Function Documentation

◆ DFT()

void Nektar::SolverUtils::FileFieldInterpolator::DFT ( const std::string  file,
const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
const bool  timefromfile 
)
protected

Definition at line 656 of file FileSolution.cpp.

660{
661 for (auto it : m_variableMap)
662 {
663 int ncoeffs = pFields[it.second]->GetNcoeffs();
664 m_interp[it.second] = Array<OneD, NekDouble>(ncoeffs * m_slices, 0.0);
665 }
666
667 // Import the slides into auxiliary vector
668 // The base flow should be stored in the form "filename_%d.ext"
669 // A subdirectory can also be included, such as "dir/filename_%d.ext"
670 size_t found = file.find("%d");
671 std::map<std::string, NekDouble> params;
672 if (found != string::npos)
673 {
674 ASSERTL0(file.find("%d", found + 1) == string::npos,
675 "There are more than one '%d'.");
676 int nstart = m_start;
677 std::vector<NekDouble> times;
678 for (int i = 0; i < m_slices; ++i)
679 {
680 int filen = nstart + i * m_skip;
681 auto fmt = boost::format(file) % filen;
682 ImportFldBase(fmt.str(), pFields, i, params);
683 if (m_session->GetComm()->GetRank() == 0)
684 {
685 cout << "read base flow file " << fmt.str() << endl;
686 }
687 if (timefromfile && params.count("time"))
688 {
689 times.push_back(params["time"]);
690 }
691 }
692 if (timefromfile && times.size() == m_slices && m_slices > 1)
693 {
694 m_timeStart = times[0];
695 if (m_interporder < 1)
696 {
697 m_period = m_slices * (times[m_slices - 1] - times[0]) /
698 (m_slices - 1.);
699 }
700 else
701 {
702 m_period = times[m_slices - 1] - times[0];
703 }
704 }
705 }
706 else if (m_slices == 1)
707 {
708 ImportFldBase(file.c_str(), pFields, 0, params);
709 }
710 else
711 {
712 ASSERTL0(
713 false,
714 "Since N_slices is specified, the filename provided for function "
715 "'BaseFlow' must include exactly one instance of the format "
716 "specifier '%d', to index the time-slices.");
717 }
718
719 if (!m_isperiodic || m_slices == 1)
720 {
721 return;
722 }
723
724 // Discrete Fourier Transform of the fields
725 for (auto it : m_interp)
726 {
727 int npoints = pFields[it.first]->GetNcoeffs();
728#ifdef NEKTAR_USING_FFTW
729
730 // Discrete Fourier Transform using FFTW
731 Array<OneD, NekDouble> fft_in(npoints * m_slices);
732 Array<OneD, NekDouble> fft_out(npoints * m_slices);
733
734 Array<OneD, NekDouble> m_tmpIN(m_slices);
735 Array<OneD, NekDouble> m_tmpOUT(m_slices);
736
737 // Shuffle the data
738 for (int j = 0; j < m_slices; ++j)
739 {
740 Vmath::Vcopy(npoints, &(it.second)[j * npoints], 1, &(fft_in[j]),
741 m_slices);
742 }
743
745 m_slices);
746
747 // FFT Transform
748 for (int i = 0; i < npoints; i++)
749 {
750 m_FFT->FFTFwdTrans(m_tmpIN = fft_in + i * m_slices,
751 m_tmpOUT = fft_out + i * m_slices);
752 }
753
754 // Reshuffle data
755 for (int s = 0; s < m_slices; ++s)
756 {
757 Vmath::Vcopy(npoints, &fft_out[s], m_slices,
758 &(it.second)[s * npoints], 1);
759 }
760
761 Vmath::Zero(fft_in.size(), &fft_in[0], 1);
762 Vmath::Zero(fft_out.size(), &fft_out[0], 1);
763#else
764 // Discrete Fourier Transform using MVM
765 DNekBlkMatSharedPtr blkmat;
766 blkmat = GetFloquetBlockMatrix(npoints);
767
768 int nrows = blkmat->GetRows();
769 int ncols = blkmat->GetColumns();
770
771 Array<OneD, NekDouble> sortedinarray(ncols);
772 Array<OneD, NekDouble> sortedoutarray(nrows);
773
774 // Shuffle the data
775 for (int j = 0; j < m_slices; ++j)
776 {
777 Vmath::Vcopy(npoints, &(it.second)[j * npoints], 1,
778 &(sortedinarray[j]), m_slices);
779 }
780
781 // Create NekVectors from the given data arrays
782 NekVector<NekDouble> in(ncols, sortedinarray, eWrapper);
783 NekVector<NekDouble> out(nrows, sortedoutarray, eWrapper);
784
785 // Perform matrix-vector multiply.
786 out = (*blkmat) * in;
787
788 // Reshuffle data
789 for (int s = 0; s < m_slices; ++s)
790 {
791 Vmath::Vcopy(npoints, &sortedoutarray[s], m_slices,
792 &(it.second)[s * npoints], 1);
793 }
794
795 for (int r = 0; r < sortedinarray.size(); ++r)
796 {
797 sortedinarray[0] = 0;
798 sortedoutarray[0] = 0;
799 }
800
801#endif
802 }
803}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
LibUtilities::SessionReaderSharedPtr m_session
Definition: FileSolution.h:72
std::map< std::string, int > m_variableMap
variables
Definition: FileSolution.h:85
void ImportFldBase(std::string pInfile, const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, int slice, std::map< std::string, NekDouble > &params)
Import Base flow.
std::map< int, Array< OneD, NekDouble > > m_interp
interpolation vector
Definition: FileSolution.h:83
DNekBlkMatSharedPtr GetFloquetBlockMatrix(int nexp)
NektarFFTFactory & GetNektarFFTFactory()
Definition: NektarFFT.cpp:65
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:77
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::eWrapper, CellMLToNektar.pycml::format, GetFloquetBlockMatrix(), Nektar::LibUtilities::GetNektarFFTFactory(), ImportFldBase(), m_interp, m_interporder, m_isperiodic, m_period, m_session, m_skip, m_slices, m_start, m_timeStart, m_variableMap, Vmath::Vcopy(), and Vmath::Zero().

Referenced by InitObject().

◆ GetFloquetBlockMatrix()

DNekBlkMatSharedPtr Nektar::SolverUtils::FileFieldInterpolator::GetFloquetBlockMatrix ( int  nexp)
protected

Definition at line 621 of file FileSolution.cpp.

622{
623 DNekMatSharedPtr loc_mat;
624 DNekBlkMatSharedPtr BlkMatrix;
625
626 Array<OneD, unsigned int> nrows(nexp);
627 Array<OneD, unsigned int> ncols(nexp);
628
629 nrows = Array<OneD, unsigned int>(nexp, m_slices);
630 ncols = Array<OneD, unsigned int>(nexp, m_slices);
631
632 MatrixStorage blkmatStorage = eDIAGONAL;
633 BlkMatrix = MemoryManager<DNekBlkMat>::AllocateSharedPtr(nrows, ncols,
634 blkmatStorage);
635
636 const LibUtilities::PointsKey Pkey(m_slices,
638 const LibUtilities::BasisKey BK(LibUtilities::eFourier, m_slices, Pkey);
639 StdRegions::StdSegExp StdSeg(BK);
640
641 StdRegions::StdMatrixKey matkey(StdRegions::eFwdTrans,
642 StdSeg.DetShapeType(), StdSeg);
643
644 loc_mat = StdSeg.GetStdMatrix(matkey);
645
646 // set up array of block matrices.
647 for (int i = 0; i < nexp; ++i)
648 {
649 BlkMatrix->SetBlock(i, i, loc_mat);
650 }
651
652 return BlkMatrix;
653}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:74
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::eDIAGONAL, Nektar::LibUtilities::eFourier, Nektar::LibUtilities::eFourierEvenlySpaced, Nektar::StdRegions::eFwdTrans, Nektar::StdRegions::StdExpansion::GetStdMatrix(), and m_slices.

Referenced by DFT().

◆ GetStartTime()

NekDouble Nektar::SolverUtils::FileFieldInterpolator::GetStartTime ( )

Definition at line 805 of file FileSolution.cpp.

806{
807 return m_timeStart;
808}

References m_timeStart.

◆ ImportFldBase()

void Nektar::SolverUtils::FileFieldInterpolator::ImportFldBase ( std::string  pInfile,
const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
int  pSlice,
std::map< std::string, NekDouble > &  params 
)
protected

Import Base flow.

Import field from infile and load into m_fields. This routine will also perform a BwdTrans to ensure data is in both the physical and coefficient storage.

Parameters
pInFileFilename to read.
pFieldsArray of expansion lists

Definition at line 462 of file FileSolution.cpp.

466{
467 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
468 std::vector<std::vector<NekDouble>> FieldData;
469
470 int numexp = pFields[0]->GetExpSize();
471 Array<OneD, int> ElementGIDs(numexp);
472
473 // Define list of global element ids
474 for (int i = 0; i < numexp; ++i)
475 {
476 ElementGIDs[i] = pFields[0]->GetExp(i)->GetGeom()->GetGlobalID();
477 }
478
479 // Get Homogeneous
482 fld->Import(pInfile, FieldDef, FieldData,
484
485 int nFileVar = FieldDef[0]->m_fields.size();
486 for (int j = 0; j < nFileVar; ++j)
487 {
488 std::string var = FieldDef[0]->m_fields[j];
489 if (!m_variableMap.count(var))
490 {
491 continue;
492 }
493 int ncoeffs = pFields[m_variableMap[var]]->GetNcoeffs();
494 Array<OneD, NekDouble> tmp_coeff(ncoeffs, 0.);
495 for (int i = 0; i < FieldDef.size(); ++i)
496 {
497 pFields[m_variableMap[var]]->ExtractDataToCoeffs(
498 FieldDef[i], FieldData[i], FieldDef[i]->m_fields[j], tmp_coeff);
499 }
500 Vmath::Vcopy(ncoeffs, &tmp_coeff[0], 1,
501 &m_interp[m_variableMap[var]][pSlice * ncoeffs], 1);
502 }
503
504 LibUtilities::FieldMetaDataMap fieldMetaDataMap;
505 fld->ImportFieldMetaData(pInfile, fieldMetaDataMap);
506 // check to see if time defined
507 if (fieldMetaDataMap != LibUtilities::NullFieldMetaDataMap)
508 {
509 auto iter = fieldMetaDataMap.find("Time");
510 if (iter != fieldMetaDataMap.end())
511 {
512 params["time"] = std::stod(iter->second);
513 }
514 }
515}
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
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:322
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:50
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:51

References Nektar::LibUtilities::FieldIO::CreateForFile(), m_interp, m_session, m_variableMap, Nektar::LibUtilities::NullFieldMetaDataMap, and Vmath::Vcopy().

Referenced by DFT().

◆ InitObject()

void Nektar::SolverUtils::FileFieldInterpolator::InitObject ( const std::string  functionName,
LibUtilities::SessionReaderSharedPtr  pSession,
const Array< OneD, const MultiRegions::ExpListSharedPtr pFields,
std::set< std::string > &  variables,
std::map< std::string, int > &  series,
std::map< std::string, NekDouble > &  time 
)

Definition at line 341 of file FileSolution.cpp.

347{
348 m_session = pSession;
349 for (size_t i = 0; i < m_session->GetVariables().size(); ++i)
350 {
351 std::string var = m_session->GetVariable(i);
352 if (variables.count(var))
353 {
354 ASSERTL0(m_session->DefinesFunction(functionName, var) &&
355 (m_session->GetFunctionType(functionName, var) ==
357 functionName + "(" + var + ") is not defined as a file.");
358 m_variableMap[var] = i;
359 }
360 }
361
362 if (series.count("start"))
363 {
364 m_start = series["start"];
365 }
366 else
367 {
368 m_start = 0;
369 }
370
371 if (series.count("skip"))
372 {
373 m_skip = series["skip"];
374 }
375 else
376 {
377 m_skip = 1;
378 }
379
380 if (series.count("slices"))
381 {
382 m_slices = series["slices"];
383 }
384 else
385 {
386 m_slices = 1;
387 }
388
389 if (series.count("order"))
390 {
391 m_interporder = series["order"];
392 }
393 else
394 {
395 m_interporder = 0;
396 }
397
398 if (series.count("isperiodic"))
399 {
400 m_isperiodic = series["isperiodic"];
401 }
402 else
403 {
404 m_isperiodic = 1;
405 }
406
407 bool timefromfile = false;
408 if (times.count("period"))
409 {
410 m_period = times["period"];
411 }
412 else if (m_slices > 1)
413 {
414 timefromfile = true;
415 }
416 if (times.count("start"))
417 {
418 m_timeStart = times["start"];
419 }
420 else
421 {
422 m_timeStart = 0.;
423 }
424
425 if (m_session->GetComm()->GetRank() == 0)
426 {
427 cout << "baseflow info : interpolation order " << m_interporder
428 << ", period " << m_period << ", periodicity ";
429 if (m_isperiodic)
430 {
431 cout << "yes\n";
432 }
433 else
434 {
435 cout << "no\n";
436 }
437 cout << "baseflow info : files from " << m_start << " to "
438 << (m_start + (m_slices - 1) * m_skip) << " (skip " << m_skip
439 << ") with " << (m_slices - (m_interporder > 1))
440 << " time intervals" << endl;
441 }
442
443 string file = m_session->GetFunctionFilename("Solution", 0);
444 DFT(file, pFields, timefromfile);
445}
void DFT(const std::string file, const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const bool timefromfile)

References ASSERTL0, DFT(), Nektar::LibUtilities::eFunctionTypeFile, m_interporder, m_isperiodic, m_period, m_session, m_skip, m_slices, m_start, m_timeStart, and m_variableMap.

◆ InterpolateField() [1/2]

void Nektar::SolverUtils::FileFieldInterpolator::InterpolateField ( const int  v,
Array< OneD, NekDouble > &  outarray,
const NekDouble  time 
)

Definition at line 528 of file FileSolution.cpp.

531{
532 // doesnot have this variable
533 if (!m_interp.count(v))
534 {
535 return;
536 }
537 // one slice, steady solution
538 int npoints = m_interp[v].size() / m_slices;
539 if (m_slices == 1)
540 {
541 Vmath::Vcopy(npoints, &m_interp[v][0], 1, &outarray[0], 1);
542 return;
543 }
544 // unsteady solution
545 time -= m_timeStart;
546 if (m_isperiodic && time > m_period)
547 {
548 time = fmod(time, m_period);
549 if (time < 0.)
550 {
551 time += m_period;
552 }
553 }
554 if (m_interporder < 1)
555 {
556 NekDouble BetaT = 2 * M_PI * fmod(time, m_period) / m_period;
557 NekDouble phase;
558 Array<OneD, NekDouble> auxiliary(npoints);
559
560 Vmath::Vcopy(npoints, &m_interp[v][0], 1, &outarray[0], 1);
561 Vmath::Svtvp(npoints, cos(0.5 * m_slices * BetaT),
562 &m_interp[v][npoints], 1, &outarray[0], 1, &outarray[0],
563 1);
564
565 for (int i = 2; i < m_slices; i += 2)
566 {
567 phase = (i >> 1) * BetaT;
568
569 Vmath::Svtvp(npoints, cos(phase), &m_interp[v][i * npoints], 1,
570 &outarray[0], 1, &outarray[0], 1);
571 Vmath::Svtvp(npoints, -sin(phase), &m_interp[v][(i + 1) * npoints],
572 1, &outarray[0], 1, &outarray[0], 1);
573 }
574 }
575 else
576 {
577 NekDouble x = time;
578 x = x / m_period * (m_slices - 1);
579 int ix = x;
580 if (ix < 0)
581 {
582 ix = 0;
583 }
584 if (ix > m_slices - 2)
585 {
586 ix = m_slices - 2;
587 }
588 int padleft = max(0, m_interporder / 2 - 1);
589 if (padleft > ix)
590 {
591 padleft = ix;
592 }
593 int padright = m_interporder - 1 - padleft;
594 if (padright > m_slices - 1 - ix)
595 {
596 padright = m_slices - 1 - ix;
597 }
598 padleft = m_interporder - 1 - padright;
599 Array<OneD, NekDouble> coeff(m_interporder, 1.);
600 for (int i = 0; i < m_interporder; ++i)
601 {
602 for (int j = 0; j < m_interporder; ++j)
603 {
604 if (i != j)
605 {
606 coeff[i] *= (x - ix + padleft - (NekDouble)j) /
607 ((NekDouble)i - (NekDouble)j);
608 }
609 }
610 }
611 Vmath::Zero(npoints, &outarray[0], 1);
612 for (int i = ix - padleft; i < ix + padright + 1; ++i)
613 {
614 Vmath::Svtvp(npoints, coeff[i - ix + padleft],
615 &m_interp[v][i * npoints], 1, &outarray[0], 1,
616 &outarray[0], 1);
617 }
618 }
619}
double NekDouble
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
Definition: Vmath.hpp:396

References m_interp, m_interporder, m_isperiodic, m_period, m_slices, m_timeStart, Vmath::Svtvp(), Vmath::Vcopy(), and Vmath::Zero().

◆ InterpolateField() [2/2]

void Nektar::SolverUtils::FileFieldInterpolator::InterpolateField ( const std::string  variable,
Array< OneD, NekDouble > &  outarray,
NekDouble  time 
)

Definition at line 517 of file FileSolution.cpp.

520{
521 if (!m_variableMap.count(variable))
522 {
523 return;
524 }
525 InterpolateField(m_variableMap[variable], outarray, time);
526}
void InterpolateField(const std::string variable, Array< OneD, NekDouble > &outarray, NekDouble time)

References InterpolateField(), and m_variableMap.

Referenced by InterpolateField().

Friends And Related Function Documentation

◆ MemoryManager< FileFieldInterpolator >

friend class MemoryManager< FileFieldInterpolator >
friend

Definition at line 51 of file FileSolution.h.

Member Data Documentation

◆ m_interp

std::map<int, Array<OneD, NekDouble> > Nektar::SolverUtils::FileFieldInterpolator::m_interp
protected

interpolation vector

Definition at line 83 of file FileSolution.h.

Referenced by DFT(), ImportFldBase(), and InterpolateField().

◆ m_interporder

int Nektar::SolverUtils::FileFieldInterpolator::m_interporder
protected

Definition at line 81 of file FileSolution.h.

Referenced by DFT(), InitObject(), and InterpolateField().

◆ m_isperiodic

int Nektar::SolverUtils::FileFieldInterpolator::m_isperiodic
protected

Definition at line 80 of file FileSolution.h.

Referenced by DFT(), InitObject(), and InterpolateField().

◆ m_period

NekDouble Nektar::SolverUtils::FileFieldInterpolator::m_period
protected

Definition at line 79 of file FileSolution.h.

Referenced by DFT(), InitObject(), and InterpolateField().

◆ m_session

LibUtilities::SessionReaderSharedPtr Nektar::SolverUtils::FileFieldInterpolator::m_session
protected

Definition at line 72 of file FileSolution.h.

Referenced by DFT(), ImportFldBase(), and InitObject().

◆ m_skip

int Nektar::SolverUtils::FileFieldInterpolator::m_skip
protected

Definition at line 75 of file FileSolution.h.

Referenced by DFT(), and InitObject().

◆ m_slices

int Nektar::SolverUtils::FileFieldInterpolator::m_slices
protected

Definition at line 76 of file FileSolution.h.

Referenced by DFT(), GetFloquetBlockMatrix(), InitObject(), and InterpolateField().

◆ m_start

int Nektar::SolverUtils::FileFieldInterpolator::m_start
protected

number of slices

Definition at line 74 of file FileSolution.h.

Referenced by DFT(), and InitObject().

◆ m_timeStart

NekDouble Nektar::SolverUtils::FileFieldInterpolator::m_timeStart
protected

period length

Definition at line 78 of file FileSolution.h.

Referenced by DFT(), GetStartTime(), InitObject(), and InterpolateField().

◆ m_variableMap

std::map<std::string, int> Nektar::SolverUtils::FileFieldInterpolator::m_variableMap
protected

variables

Definition at line 85 of file FileSolution.h.

Referenced by DFT(), ImportFldBase(), InitObject(), and InterpolateField().