Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Attributes | Friends | List of all members
Nektar::FilterOffsetPhase Class Reference

#include <FilterOffsetPhase.h>

Inheritance diagram for Nektar::FilterOffsetPhase:
[legend]

Public Member Functions

 FilterOffsetPhase (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const ParamMap &pParams)
 
 ~FilterOffsetPhase () override
 
- 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 SolverUtils::FilterSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const ParamMap &pParams)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of the class. More...
 

Protected Member Functions

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
 
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
 

Private Attributes

unsigned int m_index
 
unsigned int m_outputIndex
 
unsigned int m_outputFrequency
 
std::string m_outputFile
 
CellModelSharedPtr m_cell
 
LibUtilities::FieldIOSharedPtr m_fld
 
double meanV = -57.2986
 
Array< OneD, NekDoubleoldV
 

Friends

class MemoryManager< FilterOffsetPhase >
 

Additional Inherited Members

- Public Types inherited from Nektar::SolverUtils::Filter
typedef std::map< std::string, std::string > ParamMap
 
- Protected Attributes inherited from Nektar::SolverUtils::Filter
LibUtilities::SessionReaderSharedPtr m_session
 
const std::weak_ptr< EquationSystemm_equ
 

Detailed Description

Definition at line 47 of file FilterOffsetPhase.h.

Constructor & Destructor Documentation

◆ FilterOffsetPhase()

Nektar::FilterOffsetPhase::FilterOffsetPhase ( const LibUtilities::SessionReaderSharedPtr pSession,
const std::weak_ptr< SolverUtils::EquationSystem > &  pEquation,
const ParamMap pParams 
)

Definition at line 44 of file FilterOffsetPhase.cpp.

48 : Filter(pSession, pEquation)
49{
50 if (pParams.find("OutputFile") == pParams.end())
51 {
52 m_outputFile = m_session->GetSessionName();
53 }
54 else
55 {
56 ASSERTL0(!(pParams.find("OutputFile")->second.empty()),
57 "Missing parameter 'OutputFile'.");
58 m_outputFile = pParams.find("OutputFile")->second;
59 }
60 ASSERTL0(pParams.find("OutputFrequency") != pParams.end(),
61 "Missing parameter 'OutputFrequency'.");
62 LibUtilities::Equation equ(m_session->GetInterpreter(),
63 pParams.find("OutputFrequency")->second);
64 m_outputFrequency = floor(equ.Evaluate());
65
66 if (pParams.find("MeanV") != pParams.end())
67 {
68 LibUtilities::Equation equ(m_session->GetInterpreter(),
69 pParams.find("MeanV")->second);
70 meanV = equ.Evaluate();
71 }
72
73 m_outputIndex = 0;
74 m_index = 0;
76}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
LibUtilities::FieldIOSharedPtr m_fld
static std::shared_ptr< FieldIO > CreateDefault(const LibUtilities::SessionReaderSharedPtr session)
Returns an object for the default FieldIO method.
Definition: FieldIO.cpp:195
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:83
SOLVER_UTILS_EXPORT Filter(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation)
Definition: Filter.cpp:45

References ASSERTL0, Nektar::LibUtilities::FieldIO::CreateDefault(), Nektar::LibUtilities::Equation::Evaluate(), m_fld, m_index, m_outputFile, m_outputFrequency, m_outputIndex, Nektar::SolverUtils::Filter::m_session, and meanV.

◆ ~FilterOffsetPhase()

Nektar::FilterOffsetPhase::~FilterOffsetPhase ( )
override

Definition at line 78 of file FilterOffsetPhase.cpp.

79{
80}

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 53 of file FilterOffsetPhase.h.

57 {
60 pSession, pEquation, pParams);
61 return p;
62 }
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.

◆ v_Finalise()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 145 of file FilterOffsetPhase.cpp.

149{
150}

◆ v_Initialise()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 82 of file FilterOffsetPhase.cpp.

85{
86 m_index = 0;
87 m_outputIndex = 0;
88
89 oldV = pFields[0]->GetPhys(); // Allocate storage for past data
90
91 v_Update(pFields, 0.0);
92}
Array< OneD, NekDouble > oldV
void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override

References m_index, m_outputIndex, oldV, and v_Update().

◆ v_IsTimeDependent()

bool Nektar::FilterOffsetPhase::v_IsTimeDependent ( )
overrideprotectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 152 of file FilterOffsetPhase.cpp.

153{
154 return true;
155}

◆ v_Update()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 94 of file FilterOffsetPhase.cpp.

97{
98 if (m_index++ % m_outputFrequency > 0)
99 {
100 return;
101 }
102
103 std::stringstream vOutputFilename;
104 vOutputFilename << m_outputFile << "_" << m_outputIndex << ".chk";
105
106 // Get FieldDefinitions
107 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef =
108 pFields[0]->GetFieldDefinitions();
109 std::vector<std::vector<NekDouble>> FieldData(FieldDef.size());
110
111 Array<OneD, NekDouble> currentV = pFields[0]->GetPhys();
112 Array<OneD, NekDouble> phase(currentV.size());
113 Array<OneD, NekDouble> phase_coeff(pFields[0]->GetNcoeffs());
114
115 // Calculate phase angle in physical space:
116 for (int i = 0; i < currentV.size(); ++i)
117 {
118 phase[i] = atan2(currentV[i] - meanV,
119 oldV[i] - meanV); // subtract mean -57.2986
120 }
121
122 // Convert to coefficient space:
123 pFields[0]->FwdTransLocalElmt(phase, phase_coeff);
124
125 // for each item in FieldDef
126 // push "phase" to FieldDef[i]->m_fields
127 // call AppendFieldData(FieldDef[i], FieldData[i], data)
128 for (int i = 0; i < FieldDef.size(); ++i)
129 {
130 FieldDef[i]->m_fields.push_back("phase");
131 pFields[0]->AppendFieldData(FieldDef[i], FieldData[i], phase_coeff);
132 }
133
134 // Save currentV
135 oldV = currentV;
136
137 // Update time in field info if required
138 LibUtilities::FieldMetaDataMap fieldMetaDataMap;
139 fieldMetaDataMap["Time"] = boost::lexical_cast<std::string>(time);
140
141 m_fld->Write(vOutputFilename.str(), FieldDef, FieldData, fieldMetaDataMap);
143}
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:50

References m_fld, m_index, m_outputFile, m_outputFrequency, m_outputIndex, meanV, and oldV.

Referenced by v_Initialise().

Friends And Related Function Documentation

◆ MemoryManager< FilterOffsetPhase >

friend class MemoryManager< FilterOffsetPhase >
friend

Definition at line 1 of file FilterOffsetPhase.h.

Member Data Documentation

◆ className

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

Name of the class.

Definition at line 65 of file FilterOffsetPhase.h.

◆ m_cell

CellModelSharedPtr Nektar::FilterOffsetPhase::m_cell
private

Definition at line 90 of file FilterOffsetPhase.h.

◆ m_fld

LibUtilities::FieldIOSharedPtr Nektar::FilterOffsetPhase::m_fld
private

Definition at line 91 of file FilterOffsetPhase.h.

Referenced by FilterOffsetPhase(), and v_Update().

◆ m_index

unsigned int Nektar::FilterOffsetPhase::m_index
private

Definition at line 86 of file FilterOffsetPhase.h.

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

◆ m_outputFile

std::string Nektar::FilterOffsetPhase::m_outputFile
private

Definition at line 89 of file FilterOffsetPhase.h.

Referenced by FilterOffsetPhase(), and v_Update().

◆ m_outputFrequency

unsigned int Nektar::FilterOffsetPhase::m_outputFrequency
private

Definition at line 88 of file FilterOffsetPhase.h.

Referenced by FilterOffsetPhase(), and v_Update().

◆ m_outputIndex

unsigned int Nektar::FilterOffsetPhase::m_outputIndex
private

Definition at line 87 of file FilterOffsetPhase.h.

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

◆ meanV

double Nektar::FilterOffsetPhase::meanV = -57.2986
private

Definition at line 92 of file FilterOffsetPhase.h.

Referenced by FilterOffsetPhase(), and v_Update().

◆ oldV

Array<OneD, NekDouble> Nektar::FilterOffsetPhase::oldV
private

Definition at line 94 of file FilterOffsetPhase.h.

Referenced by v_Initialise(), and v_Update().