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

#include <FilterHilbertFFTPhase.h>

Inheritance diagram for Nektar::FilterHilbertFFTPhase:
[legend]

Public Member Functions

 FilterHilbertFFTPhase (const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_ptr< SolverUtils::EquationSystem > &pEquation, const ParamMap &pParams)
 
 ~FilterHilbertFFTPhase () override
 
- Public Member Functions inherited from Nektar::SolverUtils::Filter
SOLVER_UTILS_EXPORT Filter (const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_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 ()
 
SOLVER_UTILS_EXPORT std::string SetupOutput (const std::string ext, const ParamMap &pParams)
 
SOLVER_UTILS_EXPORT std::string SetupOutput (const std::string ext, const std::string inname)
 
SOLVER_UTILS_EXPORT void SetUpdateOnInitialise (bool flag)
 

Static Public Member Functions

static SolverUtils::FilterSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_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
 
- Protected Member Functions inherited from Nektar::SolverUtils::Filter
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
 
virtual SOLVER_UTILS_EXPORT std::string v_SetupOutput (const std::string ext, const ParamMap &pParams)
 
virtual SOLVER_UTILS_EXPORT std::string v_SetupOutput (const std::string ext, const std::string inname)
 

Private Attributes

unsigned int m_index
 
unsigned int m_outputIndex
 
unsigned int vCounter
 
unsigned int m_outputFrequency
 
int m_window
 
unsigned int m_overlap
 
bool linear_trans = true
 
bool provided_mean = false
 
std::string m_outputFile
 
LibUtilities::FieldIOSharedPtr m_fld
 
unsigned int npoints
 
NekDouble m_TimeStep
 
double m_mean
 
Array< OneD, NekDoubleoldV
 
Array< OneD, NekDoubleoldV_zero_mean
 
Array< OneD, NekDoublefcoef
 
Array< OneD, NekDoubleh
 
Array< OneD, NekDoubleoverlap_phase
 
fftw_plan plan_forward
 
fftw_plan plan_backward
 

Friends

class MemoryManager< FilterHilbertFFTPhase >
 

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
 
bool m_updateOnInitialise = true
 

Detailed Description

Definition at line 45 of file FilterHilbertFFTPhase.h.

Constructor & Destructor Documentation

◆ FilterHilbertFFTPhase()

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

Definition at line 44 of file FilterHilbertFFTPhase.cpp.

48 : Filter(pSession, pEquation)
49{
50 std::string ext = "";
51 m_outputFile = Filter::SetupOutput(ext, pParams);
52
53 ASSERTL0(pParams.find("OutputFrequency") != pParams.end(),
54 "Missing parameter 'OutputFrequency'."); // Sampling frequency =
55 // per number of timesteps
56 LibUtilities::Equation equ(m_session->GetInterpreter(),
57 pParams.find("OutputFrequency")->second);
58 m_outputFrequency = floor(equ.Evaluate());
59
60 ASSERTL0(pParams.find("WindowSize") != pParams.end(),
61 "Missing parameter 'WindowSize'.");
62 LibUtilities::Equation equ1(m_session->GetInterpreter(),
63 pParams.find("WindowSize")->second);
64 m_window = floor(equ1.Evaluate()); // No. of output steps taken in each FFT
65 ASSERTL0(m_window % 2 == 0, "'WindowSize' must be even.");
66
67 ASSERTL0(pParams.find("OverlapSize") != pParams.end(),
68 "Missing parameter 'OverlapSize'.");
69 LibUtilities::Equation equ2(m_session->GetInterpreter(),
70 pParams.find("OverlapSize")->second);
71 m_overlap = floor(equ2.Evaluate());
73 "'OverlapSize' must be smaller than 'WindowSize'.");
74
75 if (pParams.find("MeanV") != pParams.end())
76 {
77 provided_mean = true;
78 LibUtilities::Equation equ3(m_session->GetInterpreter(),
79 pParams.find("MeanV")->second);
80 m_mean = equ3.Evaluate();
81 }
82
83 if (pParams.find("LinearTransitionOverlap") != pParams.end())
84 {
86 (pParams.find("LinearTransitionOverlap")->second == "true");
87 }
88
89 m_outputIndex = 0;
90 m_index = 0;
91
93}
#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:194
SOLVER_UTILS_EXPORT std::string SetupOutput(const std::string ext, const ParamMap &pParams)
Definition: Filter.h:139
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:93
SOLVER_UTILS_EXPORT Filter(const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_ptr< EquationSystem > &pEquation)

References ASSERTL0, Nektar::LibUtilities::FieldIO::CreateDefault(), Nektar::LibUtilities::Equation::Evaluate(), linear_trans, m_fld, m_index, m_mean, m_outputFile, m_outputFrequency, m_outputIndex, m_overlap, Nektar::SolverUtils::Filter::m_session, m_window, provided_mean, and Nektar::SolverUtils::Filter::SetupOutput().

◆ ~FilterHilbertFFTPhase()

Nektar::FilterHilbertFFTPhase::~FilterHilbertFFTPhase ( )
override

Definition at line 95 of file FilterHilbertFFTPhase.cpp.

96{
97}

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 51 of file FilterHilbertFFTPhase.h.

55 {
58 pSession, pEquation, pParams);
59 return p;
60 }
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:52

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

◆ v_Finalise()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 291 of file FilterHilbertFFTPhase.cpp.

295{
296}

◆ v_Initialise()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 99 of file FilterHilbertFFTPhase.cpp.

102{
103 m_index = 0;
104 m_outputIndex = 0;
105 vCounter = 0;
106
107 m_TimeStep = pFields[0]->GetSession()->GetParameter("TimeStep");
108
109 npoints = pFields[0]->GetPhys().size();
110 oldV = Array<OneD, NekDouble>(npoints * m_window);
111 oldV_zero_mean = Array<OneD, NekDouble>(npoints * m_window);
112 fcoef = Array<OneD, NekDouble>(npoints * m_window);
113 h = Array<OneD, NekDouble>(npoints * m_window);
114 overlap_phase = Array<OneD, NekDouble>(npoints * m_overlap);
115 std::fill_n(&overlap_phase[0], npoints * m_overlap, 0.0);
116
117 // Initialise FFTW
118 fftw_r2r_kind *fwd_kind = new fftw_r2r_kind[npoints];
119 std::fill_n(fwd_kind, npoints, FFTW_R2HC);
120 plan_forward = fftw_plan_many_r2r(
122 &fcoef[0], &m_window, npoints, 1, &fwd_kind[0], FFTW_ESTIMATE);
123
124 fftw_r2r_kind *bwd_kind = new fftw_r2r_kind[npoints];
125 std::fill_n(bwd_kind, npoints, FFTW_HC2R);
126 plan_backward = fftw_plan_many_r2r(1, &m_window, npoints, &fcoef[0],
127 &m_window, npoints, 1, &h[0], &m_window,
128 npoints, 1, &bwd_kind[0], FFTW_ESTIMATE);
129
130 delete[] fwd_kind;
131 delete[] bwd_kind;
132 v_Update(pFields, 0.0);
133}
Array< OneD, NekDouble > overlap_phase
void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
Array< OneD, NekDouble > oldV_zero_mean
Array< OneD, NekDouble > fcoef

References fcoef, h, m_index, m_outputIndex, m_overlap, m_TimeStep, m_window, npoints, oldV, oldV_zero_mean, overlap_phase, plan_backward, plan_forward, v_Update(), and vCounter.

◆ v_IsTimeDependent()

bool Nektar::FilterHilbertFFTPhase::v_IsTimeDependent ( )
overrideprotectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 298 of file FilterHilbertFFTPhase.cpp.

299{
300 return true;
301}

◆ v_Update()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 135 of file FilterHilbertFFTPhase.cpp.

138{
139 if (m_index++ % m_outputFrequency > 0)
140 {
141 return;
142 }
143
144 // Store old data
145 Array<OneD, NekDouble> currentV = pFields[0]->GetPhys();
146 memcpy(&oldV[vCounter * npoints], &currentV[0],
147 npoints * sizeof(NekDouble));
148 vCounter++;
149
150 // Do Hilbert Transform and output phase
151 if (vCounter == m_window)
152 {
153 // Zero mean
154 if (!provided_mean)
155 {
156 for (int j = 0; j < npoints; ++j)
157 {
158 m_mean = 0.0;
159 for (int i = 0; i < m_window; ++i)
160 {
161 m_mean += oldV[i * npoints + j];
162 }
163 m_mean = m_mean / (double)m_window; // time-average
164
165 for (int i = 0; i < m_window; ++i)
166 {
167 oldV_zero_mean[i * npoints + j] =
168 oldV[i * npoints + j] - m_mean;
169 }
170 }
171 }
172 else // minus mean directly
173 {
174 for (int i = 0; i < npoints * m_window; ++i)
175 {
176 oldV_zero_mean[i] = oldV[i] - m_mean;
177 }
178 }
179
180 // Do Hilbert Transform:
181 fftw_execute(plan_forward);
182 // Rearrange out: imag->real, -real->imag
183 double tmp_fcoef;
184 for (int i = 1; i < m_window / 2; ++i)
185 {
186 for (int j = 0; j < npoints; ++j)
187 {
188 tmp_fcoef = -fcoef[i * npoints + j] / m_window;
189 fcoef[i * npoints + j] =
190 fcoef[(m_window - i) * npoints + j] / m_window;
191 fcoef[(m_window - i) * npoints + j] = tmp_fcoef;
192 }
193 }
194 std::fill_n(&fcoef[0], npoints, 0.0);
195 std::fill_n(&fcoef[(m_window / 2) * npoints], npoints, 0.0);
196 fftw_execute(plan_backward);
197
198 // Output batch
199 Array<OneD, NekDouble> phase_coeff(pFields[0]->GetNcoeffs());
200 Array<OneD, NekDouble> phase(npoints);
201
202 for (int t = 0; t < m_window - m_overlap; ++t)
203 {
204 std::stringstream vTmpFilename;
205 std::string vOutputFilename;
206 std::string ext = ".chk";
207
208 vTmpFilename
209 << fs::path(m_outputFile).replace_extension("").string() << "_"
210 << m_outputIndex << ext;
211 vOutputFilename = Filter::SetupOutput(ext, vTmpFilename.str());
212
213 // Get FieldDefinitions and allocate Fielddata
214 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef =
215 pFields[0]->GetFieldDefinitions();
216 std::vector<std::vector<NekDouble>> FieldData(FieldDef.size());
217
218 // Calculate phase
219 if (t < m_overlap)
220 {
221 if (linear_trans)
222 {
223 for (int i = 0; i < npoints; ++i)
224 {
225 phase[i] = ((double)(m_overlap - t) / m_overlap) *
226 overlap_phase[t * npoints + i] +
227 ((double)t / m_overlap) *
228 atan2(h[t * npoints + i],
229 oldV_zero_mean[t * npoints + i]);
230 }
231 }
232 else // simple average
233 {
234 for (int i = 0; i < npoints; ++i)
235 {
236 phase[i] = 0.5 * overlap_phase[t * npoints + i] +
237 0.5 * atan2(h[t * npoints + i],
238 oldV_zero_mean[t * npoints + i]);
239 }
240 }
241 }
242 else
243 {
244 for (int i = 0; i < npoints; ++i)
245 {
246 phase[i] = atan2(h[t * npoints + i],
247 oldV_zero_mean[t * npoints + i]);
248 }
249 }
250
251 // Convert to coefficient space
252 pFields[0]->FwdTransLocalElmt(phase, phase_coeff);
253
254 // for each item in FieldDef
255 // push "phase" to FieldDef[i]->m_fields
256 // call AppendFieldData(FieldDef[i], FieldData[i], data)
257 for (int i = 0; i < FieldDef.size(); ++i)
258 {
259 FieldDef[i]->m_fields.push_back("phase");
260 pFields[0]->AppendFieldData(FieldDef[i], FieldData[i],
261 phase_coeff);
262 }
263
264 // Update time in field info if required
265 LibUtilities::FieldMetaDataMap fieldMetaDataMap;
266 fieldMetaDataMap["Time"] = boost::lexical_cast<std::string>(
268
269 m_fld->Write(vOutputFilename, FieldDef, FieldData,
270 fieldMetaDataMap);
272 }
273
274 // Save phase for next iter
275 for (int t = m_window - m_overlap; t < m_window; ++t)
276 {
277 for (int i = 0; i < npoints; ++i)
278 {
279 overlap_phase[(t - m_window + m_overlap) * npoints + i] =
280 atan2(h[t * npoints + i], oldV_zero_mean[t * npoints + i]);
281 }
282 }
283
284 // Move overlap part to front
285 memcpy(&oldV[0], &oldV[(m_window - m_overlap) * npoints],
286 npoints * m_overlap * sizeof(NekDouble));
288 }
289}
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:50
double NekDouble

References fcoef, h, linear_trans, m_fld, m_index, m_mean, m_outputFile, m_outputFrequency, m_outputIndex, m_overlap, m_TimeStep, m_window, npoints, oldV, oldV_zero_mean, overlap_phase, plan_backward, plan_forward, provided_mean, Nektar::SolverUtils::Filter::SetupOutput(), and vCounter.

Referenced by v_Initialise().

Friends And Related Function Documentation

◆ MemoryManager< FilterHilbertFFTPhase >

friend class MemoryManager< FilterHilbertFFTPhase >
friend

Definition at line 1 of file FilterHilbertFFTPhase.h.

Member Data Documentation

◆ className

std::string Nektar::FilterHilbertFFTPhase::className
static
Initial value:
=
static SolverUtils::FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_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.
FilterFactory & GetFilterFactory()

Name of the class.

Definition at line 63 of file FilterHilbertFFTPhase.h.

◆ fcoef

Array<OneD, NekDouble> Nektar::FilterHilbertFFTPhase::fcoef
private

Definition at line 106 of file FilterHilbertFFTPhase.h.

Referenced by v_Initialise(), and v_Update().

◆ h

Array<OneD, NekDouble> Nektar::FilterHilbertFFTPhase::h
private

Definition at line 107 of file FilterHilbertFFTPhase.h.

Referenced by v_Initialise(), and v_Update().

◆ linear_trans

bool Nektar::FilterHilbertFFTPhase::linear_trans = true
private

Definition at line 93 of file FilterHilbertFFTPhase.h.

Referenced by FilterHilbertFFTPhase(), and v_Update().

◆ m_fld

LibUtilities::FieldIOSharedPtr Nektar::FilterHilbertFFTPhase::m_fld
private

Definition at line 96 of file FilterHilbertFFTPhase.h.

Referenced by FilterHilbertFFTPhase(), and v_Update().

◆ m_index

unsigned int Nektar::FilterHilbertFFTPhase::m_index
private

Definition at line 85 of file FilterHilbertFFTPhase.h.

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

◆ m_mean

double Nektar::FilterHilbertFFTPhase::m_mean
private

Definition at line 103 of file FilterHilbertFFTPhase.h.

Referenced by FilterHilbertFFTPhase(), and v_Update().

◆ m_outputFile

std::string Nektar::FilterHilbertFFTPhase::m_outputFile
private

Definition at line 95 of file FilterHilbertFFTPhase.h.

Referenced by FilterHilbertFFTPhase(), and v_Update().

◆ m_outputFrequency

unsigned int Nektar::FilterHilbertFFTPhase::m_outputFrequency
private

Definition at line 90 of file FilterHilbertFFTPhase.h.

Referenced by FilterHilbertFFTPhase(), and v_Update().

◆ m_outputIndex

unsigned int Nektar::FilterHilbertFFTPhase::m_outputIndex
private

Definition at line 86 of file FilterHilbertFFTPhase.h.

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

◆ m_overlap

unsigned int Nektar::FilterHilbertFFTPhase::m_overlap
private

Definition at line 92 of file FilterHilbertFFTPhase.h.

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

◆ m_TimeStep

NekDouble Nektar::FilterHilbertFFTPhase::m_TimeStep
private

Definition at line 100 of file FilterHilbertFFTPhase.h.

Referenced by v_Initialise(), and v_Update().

◆ m_window

int Nektar::FilterHilbertFFTPhase::m_window
private

Definition at line 91 of file FilterHilbertFFTPhase.h.

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

◆ npoints

unsigned int Nektar::FilterHilbertFFTPhase::npoints
private

Definition at line 99 of file FilterHilbertFFTPhase.h.

Referenced by v_Initialise(), and v_Update().

◆ oldV

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

Definition at line 104 of file FilterHilbertFFTPhase.h.

Referenced by v_Initialise(), and v_Update().

◆ oldV_zero_mean

Array<OneD, NekDouble> Nektar::FilterHilbertFFTPhase::oldV_zero_mean
private

Definition at line 105 of file FilterHilbertFFTPhase.h.

Referenced by v_Initialise(), and v_Update().

◆ overlap_phase

Array<OneD, NekDouble> Nektar::FilterHilbertFFTPhase::overlap_phase
private

Definition at line 108 of file FilterHilbertFFTPhase.h.

Referenced by v_Initialise(), and v_Update().

◆ plan_backward

fftw_plan Nektar::FilterHilbertFFTPhase::plan_backward
private

Definition at line 111 of file FilterHilbertFFTPhase.h.

Referenced by v_Initialise(), and v_Update().

◆ plan_forward

fftw_plan Nektar::FilterHilbertFFTPhase::plan_forward
private

Definition at line 110 of file FilterHilbertFFTPhase.h.

Referenced by v_Initialise(), and v_Update().

◆ provided_mean

bool Nektar::FilterHilbertFFTPhase::provided_mean = false
private

Definition at line 94 of file FilterHilbertFFTPhase.h.

Referenced by FilterHilbertFFTPhase(), and v_Update().

◆ vCounter

unsigned int Nektar::FilterHilbertFFTPhase::vCounter
private

Definition at line 87 of file FilterHilbertFFTPhase.h.

Referenced by v_Initialise(), and v_Update().