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

#include <FilterLagrangianPoints.h>

Inheritance diagram for Nektar::SolverUtils::FilterLagrangianPoints:
[legend]

Classes

struct  MovingFrame
 

Public Member Functions

SOLVER_UTILS_EXPORT FilterLagrangianPoints (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
 
SOLVER_UTILS_EXPORT ~FilterLagrangianPoints () override
 
- Public Member Functions inherited from Nektar::SolverUtils::EvaluatePoints
SOLVER_UTILS_EXPORT EvaluatePoints ()
 
virtual ~EvaluatePoints ()
 
SOLVER_UTILS_EXPORT void PartitionMobilePts (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
 
SOLVER_UTILS_EXPORT void EvaluateMobilePhys (const MultiRegions::ExpListSharedPtr &pField, std::vector< Array< OneD, NekDouble > > &PhysicsData, NekDouble time)
 
SOLVER_UTILS_EXPORT void PassMobilePhysToStatic (std::map< int, std::set< int > > &callbackUpdateMobCoords)
 
SOLVER_UTILS_EXPORT void PassStaticCoordsToMobile (std::map< int, std::set< int > > &callbackUpdateMobCoords)
 
SOLVER_UTILS_EXPORT void CopyStaticPtsToMobile ()
 
SOLVER_UTILS_EXPORT void CopyMobilePtsToStatic ()
 
SOLVER_UTILS_EXPORT void Initialise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time, std::vector< std::string > &defaultValues)
 
SOLVER_UTILS_EXPORT void SetUpCommInfo ()
 
- Public Member Functions inherited from Nektar::SolverUtils::Filter
SOLVER_UTILS_EXPORT Filter (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation)
 
virtual SOLVER_UTILS_EXPORT ~Filter ()
 
SOLVER_UTILS_EXPORT void Initialise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
 
SOLVER_UTILS_EXPORT void Update (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
 
SOLVER_UTILS_EXPORT void Finalise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
 
SOLVER_UTILS_EXPORT bool IsTimeDependent ()
 

Static Public Member Functions

static FilterSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
 

Static Public Attributes

static std::string className
 

Protected Member Functions

void ExtraPhysicsVars (std::vector< std::string > &extraVars)
 
SOLVER_UTILS_EXPORT void v_Initialise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
 
void v_Update (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
 
void v_Finalise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
 
bool v_IsTimeDependent () override
 
void v_ModifyVelocity (Array< OneD, NekDouble > gcoords, NekDouble time, Array< OneD, NekDouble > vel) override
 
void GetPhysicsData (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::vector< Array< OneD, NekDouble > > &PhysicsData)
 
void OutputSamplePoints (NekDouble time)
 
void OutputStatPoints (NekDouble time)
 
- Protected Member Functions inherited from Nektar::SolverUtils::EvaluatePoints
void UpdateMobileCoords (std::map< int, Array< OneD, NekDouble > > &globalCoords)
 
void GatherMobilePhysics (std::map< int, Array< OneD, NekDouble > > &revData)
 
void PartitionLocalPoints (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::set< int > &notfound)
 
void PartitionExchangeNonlocalPoints (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::set< int > &notfound)
 
void SyncColumnComm ()
 
NekDouble GetDefaultValue (const int i, const Array< OneD, const NekDouble > x, const NekDouble time)
 
void Pack2Int (const int &a, const int &b, double &d)
 
void unPack2Int (int &a, int &b, const double &d)
 
void Pack3Int (const int &a, const int &b, const int &c, double &d)
 
void unPack3Int (int &a, int &b, int &c, const double &d)
 
void Pack2Short (const int &a, const int &b, int &c)
 
void unPack2Short (int &a, int &b, const int &c)
 
virtual void v_ModifyVelocity (Array< OneD, NekDouble > gcoords, NekDouble time, Array< OneD, NekDouble > vel)
 
virtual void v_Initialise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)=0
 
virtual void v_Update (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)=0
 
virtual void v_Finalise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)=0
 
virtual bool v_IsTimeDependent ()=0
 

Protected Attributes

std::vector< std::string > m_defaultValues
 
struct Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame m_frame
 
std::vector< NekDoublem_box
 
std::map< std::string, std::string > v_params
 
unsigned int m_index = 0
 
unsigned int m_outputFrequency
 
std::string m_outputFile
 
bool m_outputL2Norm
 
bool m_isMovablePoints
 
unsigned int m_outputSampleFrequency
 
std::set< int > m_samplePointIDs
 
std::ofstream m_ofstreamSamplePoints
 
- Protected Attributes inherited from Nektar::SolverUtils::EvaluatePoints
int m_spacedim
 
int m_rank
 
int m_rowRank
 
int m_columnRank
 
bool m_isHomogeneous1D
 
int m_nPhysics
 
LibUtilities::CommSharedPtr m_comm
 
StationaryPointsSharedPtr m_staticPts
 
std::map< int, MobilePointSharedPtrm_mobilePts
 
std::map< int, LibUtilities::EquationSharedPtrm_defaultValueFunction
 
std::map< int, int > m_recvMobInfo
 
std::map< int, int > m_recvStatInfo
 
- Protected Attributes inherited from Nektar::SolverUtils::Filter
LibUtilities::SessionReaderSharedPtr m_session
 
const std::weak_ptr< EquationSystemm_equ
 

Friends

class MemoryManager< FilterLagrangianPoints >
 

Additional Inherited Members

- Public Types inherited from Nektar::SolverUtils::Filter
typedef std::map< std::string, std::string > ParamMap
 

Detailed Description

Definition at line 97 of file FilterLagrangianPoints.h.

Constructor & Destructor Documentation

◆ FilterLagrangianPoints()

Nektar::SolverUtils::FilterLagrangianPoints::FilterLagrangianPoints ( const LibUtilities::SessionReaderSharedPtr pSession,
const std::weak_ptr< EquationSystem > &  pEquation,
const std::map< std::string, std::string > &  pParams 
)

Definition at line 447 of file FilterLagrangianPoints.cpp.

451 : Filter(pSession, pEquation), v_params(pParams)
452{
453}
SOLVER_UTILS_EXPORT Filter(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation)
Definition: Filter.cpp:45
std::map< std::string, std::string > v_params

◆ ~FilterLagrangianPoints()

Nektar::SolverUtils::FilterLagrangianPoints::~FilterLagrangianPoints ( )
override

Definition at line 458 of file FilterLagrangianPoints.cpp.

References m_ofstreamSamplePoints.

Member Function Documentation

◆ create()

static FilterSharedPtr Nektar::SolverUtils::FilterLagrangianPoints::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const std::weak_ptr< EquationSystem > &  pEquation,
const std::map< std::string, std::string > &  pParams 
)
inlinestatic

Definition at line 102 of file FilterLagrangianPoints.h.

106 {
109 pSession, pEquation, pParams);
110 return p;
111 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< Filter > FilterSharedPtr
A shared pointer to a Driver object.
Definition: Filter.h:51

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ ExtraPhysicsVars()

void Nektar::SolverUtils::FilterLagrangianPoints::ExtraPhysicsVars ( std::vector< std::string > &  extraVars)
protected

Definition at line 678 of file FilterLagrangianPoints.cpp.

680{
681 if (m_spacedim == 2)
682 {
683 extraVars.push_back("u_x");
684 extraVars.push_back("u_y");
685 extraVars.push_back("v_x");
686 extraVars.push_back("v_y");
687 extraVars.push_back("LW_z");
688 }
689 else if (m_spacedim == 3)
690 {
691 extraVars.push_back("u_x");
692 extraVars.push_back("u_y");
693 extraVars.push_back("u_z");
694 extraVars.push_back("v_x");
695 extraVars.push_back("v_y");
696 extraVars.push_back("v_z");
697 extraVars.push_back("w_x");
698 extraVars.push_back("w_y");
699 extraVars.push_back("w_z");
700 extraVars.push_back("LW_x");
701 extraVars.push_back("LW_y");
702 extraVars.push_back("LW_z");
703 }
704}

References Nektar::SolverUtils::EvaluatePoints::m_spacedim.

Referenced by v_Initialise().

◆ GetPhysicsData()

void Nektar::SolverUtils::FilterLagrangianPoints::GetPhysicsData ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
std::vector< Array< OneD, NekDouble > > &  PhysicsData 
)
protected

Definition at line 706 of file FilterLagrangianPoints.cpp.

709{
710 int npts = pFields[0]->GetTotPoints();
711 for (int i = 0; i < m_spacedim; ++i)
712 {
713 PhysicsData.push_back(pFields[i]->UpdatePhys());
714 }
715 int offset = PhysicsData.size();
716 for (int d = 0; d < m_spacedim; ++d)
717 {
718 Array<OneD, Array<OneD, NekDouble>> data(m_spacedim);
719 for (int i = 0; i < m_spacedim; ++i)
720 {
721 data[i] = Array<OneD, NekDouble>(npts, 0.);
722 }
723 if (m_spacedim == 2)
724 {
725 pFields[d]->PhysDeriv(pFields[d]->UpdatePhys(), data[0], data[1]);
726 }
727 else if (m_spacedim == 3)
728 {
729 pFields[d]->PhysDeriv(pFields[d]->UpdatePhys(), data[0], data[1],
730 data[2]);
731 }
732 for (int i = 0; i < m_spacedim; ++i)
733 {
734 PhysicsData.push_back(data[i]);
735 }
736 }
737 int vortexdim = m_spacedim == 3 ? 3 : 1;
738 Array<OneD, Array<OneD, NekDouble>> vorticity(vortexdim);
739 Array<OneD, Array<OneD, NekDouble>> Ldata(vortexdim);
740 Array<OneD, Array<OneD, NekDouble>> temp(m_spacedim);
741 Array<OneD, NekDouble> temps(npts, 0.);
742 for (int d = 0; d < vortexdim; ++d)
743 {
744 Ldata[d] = Array<OneD, NekDouble>(npts, 0.);
745 vorticity[d] = Array<OneD, NekDouble>(npts, 0.);
746 }
747 for (int d = 0; d < m_spacedim; ++d)
748 {
749 temp[d] = Array<OneD, NekDouble>(npts, 0.);
750 }
751 if (m_spacedim == 2)
752 {
753 // 0 ux, 1 uy, 2 vx, 3 vy
754 Vmath::Vsub(npts, PhysicsData[offset + 2], 1, PhysicsData[offset + 1],
755 1, vorticity[0], 1); // vx-uy
756 }
757 else if (m_spacedim == 3)
758 {
759 // 0 ux, 1 uy, 2 uz, 3 vx, 4 vy, 5 vz, 6 wx, 7 wy, 8 wz
760 Vmath::Vsub(npts, PhysicsData[offset + 7], 1, PhysicsData[offset + 5],
761 1, vorticity[0], 1); // wy-vz
762 Vmath::Vsub(npts, PhysicsData[offset + 2], 1, PhysicsData[offset + 6],
763 1, vorticity[1], 1); // uz-wx
764 Vmath::Vsub(npts, PhysicsData[offset + 3], 1, PhysicsData[offset + 1],
765 1, vorticity[2], 1); // vx-uy
766 }
767 for (int d = 0; d < vortexdim; ++d)
768 {
769 if (m_spacedim == 2)
770 {
771 pFields[0]->PhysDeriv(vorticity[d], temp[0], temp[1]);
772 }
773 else if (m_spacedim == 3)
774 {
775 pFields[0]->PhysDeriv(vorticity[d], temp[0], temp[1], temp[2]);
776 }
777 for (int j = 0; j < m_spacedim; ++j)
778 {
779 pFields[0]->PhysDeriv(j, temp[j], temps);
780 Vmath::Vadd(npts, temps, 1, Ldata[d], 1, Ldata[d], 1);
781 }
782 PhysicsData.push_back(Ldata[d]);
783 }
784}
std::vector< double > d(NPUPPER *NPUPPER)
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220

References Nektar::UnitTests::d(), Nektar::SolverUtils::EvaluatePoints::m_spacedim, Vmath::Vadd(), and Vmath::Vsub().

Referenced by v_Update().

◆ OutputSamplePoints()

void Nektar::SolverUtils::FilterLagrangianPoints::OutputSamplePoints ( NekDouble  time)
protected

Definition at line 823 of file FilterLagrangianPoints.cpp.

824{
825 m_ofstreamSamplePoints << "ZONE T=\"" << time
826 << "\", I=" << m_samplePointIDs.size() << endl;
827 for (auto &pid : m_samplePointIDs)
828 {
829 m_ofstreamSamplePoints << pid << " " << time << " ";
830 Array<OneD, NekDouble> data(m_spacedim, 0.);
831 m_staticPts->GetCoordsByPID(pid, data);
832 for (int d = 0; d < m_spacedim; ++d)
833 {
834 m_ofstreamSamplePoints << data[d] + m_frame.m_frameDisp[d] << " ";
835 }
836 m_staticPts->GetPhysicsByPID(pid, data);
837 for (int d = 0; d < m_spacedim; ++d)
838 {
839 m_ofstreamSamplePoints << data[d] + m_frame.m_frameVel[d] << " ";
840 }
842 }
843}
StationaryPointsSharedPtr m_staticPts
struct Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame m_frame

References Nektar::UnitTests::d(), m_frame, Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame::m_frameDisp, Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame::m_frameVel, m_ofstreamSamplePoints, m_samplePointIDs, Nektar::SolverUtils::EvaluatePoints::m_spacedim, and Nektar::SolverUtils::EvaluatePoints::m_staticPts.

Referenced by v_Update().

◆ OutputStatPoints()

void Nektar::SolverUtils::FilterLagrangianPoints::OutputStatPoints ( NekDouble  time)
protected

Definition at line 845 of file FilterLagrangianPoints.cpp.

846{
847 boost::format filename(m_outputFile);
848 filename % m_rank % time;
849 std::map<std::string, NekDouble> params;
850 if (m_frame.m_frameVelFunction.size())
851 {
852 std::vector<std::string> COORDS = {"x", "y", "z"};
853 std::vector<std::string> VELOCI = {"u", "v", "w"};
854 for (auto it : m_frame.m_frameVelFunction)
855 {
856 params[VELOCI[it.first]] = m_frame.m_frameVel[it.first];
857 }
858 for (auto it : m_frame.m_frameDispFunction)
859 {
860 params[COORDS[it.first]] = m_frame.m_frameDisp[it.first];
861 }
862 }
863 m_staticPts->OutputData(filename.str(), m_outputL2Norm, params);
864}
std::map< int, LibUtilities::EquationSharedPtr > m_frameVelFunction
std::map< int, LibUtilities::EquationSharedPtr > m_frameDispFunction

References CellMLToNektar.pycml::format, m_frame, Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame::m_frameDisp, Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame::m_frameDispFunction, Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame::m_frameVel, Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame::m_frameVelFunction, m_outputFile, m_outputL2Norm, Nektar::SolverUtils::EvaluatePoints::m_rank, and Nektar::SolverUtils::EvaluatePoints::m_staticPts.

Referenced by v_Update().

◆ v_Finalise()

void Nektar::SolverUtils::FilterLagrangianPoints::v_Finalise ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 866 of file FilterLagrangianPoints.cpp.

870{
871}

◆ v_Initialise()

void Nektar::SolverUtils::FilterLagrangianPoints::v_Initialise ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 466 of file FilterLagrangianPoints.cpp.

469{
470 auto it = v_params.find("DefaultValues");
471 if (it != v_params.end())
472 {
474 }
476 NekDouble dt;
477 m_session->LoadParameter("TimeStep", dt);
478
479 it = v_params.find("FrameVelocity");
480 if (it != v_params.end())
481 {
482 std::vector<std::string> strs;
483 ParseUtils::GenerateVector(it->second, strs);
484 for (int i = strs.size(); i < m_spacedim; ++i)
485 {
486 strs.push_back("0");
487 }
488 for (int i = 0; i < m_spacedim; ++i)
489 {
492 pFields[0]->GetSession()->GetInterpreter(), strs[i]);
493 }
494 }
495
496 it = v_params.find("FrameDisplacement");
497 if (it != v_params.end())
498 {
499 std::vector<std::string> strs;
500 ParseUtils::GenerateVector(it->second, strs);
501 for (int i = strs.size(); i < m_spacedim; ++i)
502 {
503 strs.push_back("0");
504 }
505 for (int i = 0; i < m_spacedim; ++i)
506 {
509 pFields[0]->GetSession()->GetInterpreter(), strs[i]);
510 }
511 }
512
513 int intOrder = 1;
514 it = v_params.find("IntOrder");
515 if (it != v_params.end())
516 {
517 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
518 intOrder = round(equ.Evaluate());
519 }
520 if (intOrder == 0)
521 {
522 intOrder = 1;
523 m_isMovablePoints = false;
524 }
525 else
526 {
527 m_isMovablePoints = true;
528 }
529
530 it = v_params.find("RootOutputL2Norm");
531 if (it != v_params.end())
532 {
533 std::string value = it->second;
534 m_outputL2Norm = m_rank == 0 && (value == "1" || value == "Yes");
535 }
536 else
537 {
538 m_outputL2Norm = false;
539 }
540
541 it = v_params.find("OutputFile");
542 std::string outputFilename;
543 if (it == v_params.end())
544 {
545 outputFilename = m_session->GetSessionName();
546 }
547 else
548 {
549 ASSERTL0(it->second.length() > 0, "Missing parameter 'OutputFile'.");
550 outputFilename = it->second;
551 }
552 m_outputFile = outputFilename + "_R%05d_T%010.6lf.plt";
553
554 // OutputFrequency
555 it = v_params.find("OutputFrequency");
556 if (it == v_params.end())
557 {
559 }
560 else
561 {
562 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
563 m_outputFrequency = round(equ.Evaluate());
564 }
565
566 it = v_params.find("OutputSampleFrequency");
567 if (it == v_params.end())
568 {
570 }
571 else
572 {
573 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
574 m_outputSampleFrequency = round(equ.Evaluate());
575 }
576
577 it = v_params.find("Box_" + std::to_string(m_rank));
578 if (it != v_params.end())
579 {
580 std::vector<NekDouble> values;
581 ParseUtils::GenerateVector(it->second, values);
582 if (values.size() < m_spacedim * 3 + 1)
583 {
585 "Box formate error "
586 "(offset,Nx,Ny,Nx,xmin,xmax,ymin,ymax,zmin,zmax): " +
587 it->second);
588 }
589 m_frame.Update(time);
590 std::vector<int> Np(m_spacedim, 1);
591 m_box.resize(2 * m_spacedim, 0.);
592 int idOffset = std::round(values[0]);
593 for (int d = 0; d < m_spacedim; ++d)
594 {
595 Np[d] = std::round(values[d + 1]);
596 }
597 for (int d = 0; d < 2 * m_spacedim; ++d)
598 {
599 m_box[d] = values[1 + d + m_spacedim] - m_frame.m_frameDisp[d / 2];
600 }
601 std::vector<std::string> extraVars;
602 ExtraPhysicsVars(extraVars);
604 m_rank, m_spacedim, intOrder, idOffset, dt, Np, m_box, extraVars);
606
607 it = v_params.find("OutputSamplePoints");
608 std::vector<unsigned int> samplepts;
609 if (it != v_params.end())
610 {
611 ParseUtils::GenerateSeqVector(it->second, samplepts);
612 }
613 std::set<int> sampleIds(samplepts.begin(), samplepts.end());
614 it = v_params.find("OutputSamplePointsCondition");
615 LibUtilities::EquationSharedPtr sampleCondition;
616 if (it != v_params.end())
617 {
618 sampleCondition =
620 pFields[0]->GetSession()->GetInterpreter(), it->second);
621 }
622
623 for (int p = 0; p < m_staticPts->GetTotPoints(); ++p)
624 {
625 Array<OneD, NekDouble> data(m_spacedim, 0.);
626 m_staticPts->GetCoords(p, data);
627 int globalId = m_staticPts->LocalToGlobal(p);
628 NekDouble x = data[0] + m_frame.m_frameDisp[0],
629 y = data[1] + m_frame.m_frameDisp[1];
630 NekDouble z =
631 m_spacedim > 2 ? data[2] : 0. + m_frame.m_frameDisp[2];
632 bool issample =
633 (sampleCondition == nullptr)
634 ? false
635 : (sampleCondition->Evaluate(x, y, z, time) > 0);
636 if (sampleIds.find(globalId) != sampleIds.end() || issample)
637 {
638 m_samplePointIDs.insert(globalId);
639 }
640 }
641 }
642 if (!m_samplePointIDs.empty())
643 {
644 std::string sampleFilename = outputFilename + "_Sample_R%05d.dat";
645 boost::format filename(sampleFilename);
646 filename % m_rank;
647 m_ofstreamSamplePoints.open(filename.str(), std::ofstream::out);
648 std::vector<std::string> COORDS = {"x", "y", "z"};
649 std::vector<std::string> VELOCI = {"u", "v", "w"};
650 m_ofstreamSamplePoints << "variables = pid, time ";
651 for (int d = 0; d < m_spacedim; ++d)
652 {
653 m_ofstreamSamplePoints << COORDS[d] << " ";
654 }
655 for (int d = 0; d < m_spacedim; ++d)
656 {
657 m_ofstreamSamplePoints << VELOCI[d] << " ";
658 }
660 }
662 v_Update(pFields, time);
663}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
Nektar::ErrorUtil::NekError NekError
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:130
static bool GenerateSeqVector(const std::string &str, std::vector< unsigned int > &out)
Takes a comma-separated compressed string and converts it to entries in a vector.
Definition: ParseUtils.cpp:104
SOLVER_UTILS_EXPORT void Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time, std::vector< std::string > &defaultValues)
SOLVER_UTILS_EXPORT void SetUpCommInfo()
SOLVER_UTILS_EXPORT void CopyStaticPtsToMobile()
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:83
void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
void ExtraPhysicsVars(std::vector< std::string > &extraVars)
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:125
std::vector< double > z(NPUPPER)
double NekDouble

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::SolverUtils::EvaluatePoints::CopyStaticPtsToMobile(), Nektar::UnitTests::d(), Nektar::LibUtilities::Equation::Evaluate(), ExtraPhysicsVars(), CellMLToNektar.pycml::format, Nektar::ParseUtils::GenerateSeqVector(), Nektar::ParseUtils::GenerateVector(), Nektar::SolverUtils::EvaluatePoints::Initialise(), m_box, m_defaultValues, m_frame, Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame::m_frameDisp, Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame::m_frameDispFunction, Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame::m_frameVelFunction, m_isMovablePoints, m_ofstreamSamplePoints, m_outputFile, m_outputFrequency, m_outputL2Norm, m_outputSampleFrequency, Nektar::SolverUtils::EvaluatePoints::m_rank, m_samplePointIDs, Nektar::SolverUtils::Filter::m_session, Nektar::SolverUtils::EvaluatePoints::m_spacedim, Nektar::SolverUtils::EvaluatePoints::m_staticPts, CellMLToNektar.cellml_metadata::p, Nektar::SolverUtils::EvaluatePoints::SetUpCommInfo(), Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame::Update(), v_params, v_Update(), and Nektar::UnitTests::z().

◆ v_IsTimeDependent()

bool Nektar::SolverUtils::FilterLagrangianPoints::v_IsTimeDependent ( )
overrideprotectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 873 of file FilterLagrangianPoints.cpp.

874{
875 return true;
876}

◆ v_ModifyVelocity()

void Nektar::SolverUtils::FilterLagrangianPoints::v_ModifyVelocity ( Array< OneD, NekDouble gcoords,
NekDouble  time,
Array< OneD, NekDouble vel 
)
overrideprotectedvirtual

◆ v_Update()

void Nektar::SolverUtils::FilterLagrangianPoints::v_Update ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 786 of file FilterLagrangianPoints.cpp.

789{
790 if (m_isMovablePoints || m_index == 0)
791 {
792 PartitionMobilePts(pFields);
793 }
794 std::vector<Array<OneD, NekDouble>> PhysicsData;
795 GetPhysicsData(pFields, PhysicsData);
796 m_frame.Update(time);
797 EvaluateMobilePhys(pFields[0], PhysicsData, time);
798 std::map<int, std::set<int>> callbackUpdateMobCoords;
799 PassMobilePhysToStatic(callbackUpdateMobCoords);
800 // output Lagrangian coordinates
801 if (m_index % m_outputFrequency == 0 && m_staticPts != nullptr)
802 {
803 OutputStatPoints(time);
804 }
805 // output sample points
806 if (m_index % m_outputSampleFrequency == 0 && !m_samplePointIDs.empty())
807 {
808 OutputSamplePoints(time);
809 }
810 // update Lagrangian coordinates
812 {
813 if (m_staticPts != nullptr)
814 {
815 m_staticPts->TimeAdvance(m_index + 1);
816 }
817 // return the updated global coordinates to mobile points
818 PassStaticCoordsToMobile(callbackUpdateMobCoords);
819 }
820 ++m_index;
821}
SOLVER_UTILS_EXPORT void EvaluateMobilePhys(const MultiRegions::ExpListSharedPtr &pField, std::vector< Array< OneD, NekDouble > > &PhysicsData, NekDouble time)
SOLVER_UTILS_EXPORT void PassMobilePhysToStatic(std::map< int, std::set< int > > &callbackUpdateMobCoords)
SOLVER_UTILS_EXPORT void PartitionMobilePts(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
SOLVER_UTILS_EXPORT void PassStaticCoordsToMobile(std::map< int, std::set< int > > &callbackUpdateMobCoords)
void GetPhysicsData(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::vector< Array< OneD, NekDouble > > &PhysicsData)

References Nektar::SolverUtils::EvaluatePoints::EvaluateMobilePhys(), GetPhysicsData(), m_frame, m_index, m_isMovablePoints, m_outputFrequency, m_outputSampleFrequency, m_samplePointIDs, Nektar::SolverUtils::EvaluatePoints::m_staticPts, OutputSamplePoints(), OutputStatPoints(), Nektar::SolverUtils::EvaluatePoints::PartitionMobilePts(), Nektar::SolverUtils::EvaluatePoints::PassMobilePhysToStatic(), Nektar::SolverUtils::EvaluatePoints::PassStaticCoordsToMobile(), and Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame::Update().

Referenced by v_Initialise().

Friends And Related Function Documentation

◆ MemoryManager< FilterLagrangianPoints >

friend class MemoryManager< FilterLagrangianPoints >
friend

Definition at line 94 of file FilterLagrangianPoints.h.

Member Data Documentation

◆ className

std::string Nektar::SolverUtils::FilterLagrangianPoints::className
static
Initial value:
=
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:39

Definition at line 113 of file FilterLagrangianPoints.h.

◆ m_box

std::vector<NekDouble> Nektar::SolverUtils::FilterLagrangianPoints::m_box
protected

Definition at line 161 of file FilterLagrangianPoints.h.

Referenced by v_Initialise().

◆ m_defaultValues

std::vector<std::string> Nektar::SolverUtils::FilterLagrangianPoints::m_defaultValues
protected

Definition at line 142 of file FilterLagrangianPoints.h.

Referenced by v_Initialise().

◆ m_frame

struct Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame Nektar::SolverUtils::FilterLagrangianPoints::m_frame
protected

◆ m_index

unsigned int Nektar::SolverUtils::FilterLagrangianPoints::m_index = 0
protected

Definition at line 163 of file FilterLagrangianPoints.h.

Referenced by v_Update().

◆ m_isMovablePoints

bool Nektar::SolverUtils::FilterLagrangianPoints::m_isMovablePoints
protected

Definition at line 167 of file FilterLagrangianPoints.h.

Referenced by v_Initialise(), and v_Update().

◆ m_ofstreamSamplePoints

std::ofstream Nektar::SolverUtils::FilterLagrangianPoints::m_ofstreamSamplePoints
protected

◆ m_outputFile

std::string Nektar::SolverUtils::FilterLagrangianPoints::m_outputFile
protected

Definition at line 165 of file FilterLagrangianPoints.h.

Referenced by OutputStatPoints(), and v_Initialise().

◆ m_outputFrequency

unsigned int Nektar::SolverUtils::FilterLagrangianPoints::m_outputFrequency
protected

Definition at line 164 of file FilterLagrangianPoints.h.

Referenced by v_Initialise(), and v_Update().

◆ m_outputL2Norm

bool Nektar::SolverUtils::FilterLagrangianPoints::m_outputL2Norm
protected

Definition at line 166 of file FilterLagrangianPoints.h.

Referenced by OutputStatPoints(), and v_Initialise().

◆ m_outputSampleFrequency

unsigned int Nektar::SolverUtils::FilterLagrangianPoints::m_outputSampleFrequency
protected

Definition at line 168 of file FilterLagrangianPoints.h.

Referenced by v_Initialise(), and v_Update().

◆ m_samplePointIDs

std::set<int> Nektar::SolverUtils::FilterLagrangianPoints::m_samplePointIDs
protected

Definition at line 169 of file FilterLagrangianPoints.h.

Referenced by OutputSamplePoints(), v_Initialise(), and v_Update().

◆ v_params

std::map<std::string, std::string> Nektar::SolverUtils::FilterLagrangianPoints::v_params
protected

Definition at line 162 of file FilterLagrangianPoints.h.

Referenced by v_Initialise().