Nektar++
FilterHilbertFFTPhase.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File FilterHilbertFFTPhase.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// License for the specific language governing rights and limitations under
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: Outputs phase fields during time-stepping.
33//
34///////////////////////////////////////////////////////////////////////////////
35
37
38namespace Nektar
39{
42 "HilbertPhase", FilterHilbertFFTPhase::create);
43
46 const std::shared_ptr<SolverUtils::EquationSystem> &pEquation,
47 const ParamMap &pParams)
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}
94
96{
97}
98
101 [[maybe_unused]] const NekDouble &time)
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();
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}
134
137 [[maybe_unused]] const NekDouble &time)
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());
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}
290
293 &pFields,
294 [[maybe_unused]] const NekDouble &time)
295{
296}
297
299{
300 return true;
301}
302} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
Array< OneD, NekDouble > overlap_phase
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
FilterHilbertFFTPhase(const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_ptr< SolverUtils::EquationSystem > &pEquation, const ParamMap &pParams)
Array< OneD, NekDouble > oldV_zero_mean
void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
Array< OneD, NekDouble > fcoef
LibUtilities::FieldIOSharedPtr m_fld
static SolverUtils::FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_ptr< SolverUtils::EquationSystem > &pEquation, const ParamMap &pParams)
Creates an instance of this class.
static std::string className
Name of the class.
static std::shared_ptr< FieldIO > CreateDefault(const LibUtilities::SessionReaderSharedPtr session)
Returns an object for the default FieldIO method.
Definition: FieldIO.cpp:194
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
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
std::map< std::string, std::string > ParamMap
Definition: Filter.h:66
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:50
std::shared_ptr< SessionReader > SessionReaderSharedPtr
FilterFactory & GetFilterFactory()
double NekDouble