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

Append Reynolds stresses to the average fields. More...

#include <FilterReynoldsStresses.h>

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

Public Member Functions

SOLVER_UTILS_EXPORT FilterReynoldsStresses (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
 
SOLVER_UTILS_EXPORT ~FilterReynoldsStresses () override
 
- Public Member Functions inherited from Nektar::SolverUtils::FilterFieldConvert
SOLVER_UTILS_EXPORT FilterFieldConvert (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
 
SOLVER_UTILS_EXPORT ~FilterFieldConvert () 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 FilterSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
 Creates an instance of this class. More...
 
- Static Public Member Functions inherited from Nektar::SolverUtils::FilterFieldConvert
static FilterSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of the class. More...
 
- Static Public Attributes inherited from Nektar::SolverUtils::FilterFieldConvert
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_FillVariablesName (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields) override
 
void v_ProcessSample (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, const NekDouble &time) override
 
void v_PrepareOutput (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
 
NekDouble v_GetScale () override
 
std::string v_GetFileSuffix () override
 
- Protected Member Functions inherited from Nektar::SolverUtils::FilterFieldConvert
SOLVER_UTILS_EXPORT void v_Initialise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
 
virtual SOLVER_UTILS_EXPORT void v_FillVariablesName (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
 
SOLVER_UTILS_EXPORT void v_Update (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
 
SOLVER_UTILS_EXPORT void v_Finalise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
 
virtual SOLVER_UTILS_EXPORT void v_ProcessSample (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, const NekDouble &time)
 
virtual SOLVER_UTILS_EXPORT void v_PrepareOutput (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
 
virtual SOLVER_UTILS_EXPORT NekDouble v_GetScale ()
 
virtual SOLVER_UTILS_EXPORT std::string v_GetFileSuffix ()
 
void OutputField (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, int dump=-1)
 
SOLVER_UTILS_EXPORT bool v_IsTimeDependent () override
 
void CreateModules (std::vector< std::string > &modcmds)
 
void CreateFields (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
 
void CheckModules (std::vector< ModuleSharedPtr > &modules)
 
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< Array< OneD, NekDouble > > m_fields
 
std::vector< Array< OneD, NekDouble > > m_delta
 
NekDouble m_alpha
 
bool m_movAvg
 
- Protected Attributes inherited from Nektar::SolverUtils::FilterFieldConvert
unsigned int m_numSamples
 
unsigned int m_outputFrequency
 
unsigned int m_sampleFrequency
 
std::string m_outputFile
 
std::string m_restartFile
 
unsigned int m_index
 
unsigned int m_outputIndex
 
bool m_phaseSample
 
NekDouble m_phaseSamplePeriod
 
NekDouble m_phaseSamplePhase
 
NekDouble m_phaseTolerance
 
NekDouble m_dt
 
std::vector< ModuleSharedPtrm_modules
 
LibUtilities::FieldMetaDataMap m_fieldMetaData
 
std::vector< Array< OneD, NekDouble > > m_outFields
 
std::vector< std::string > m_variables
 
FieldSharedPtr m_f
 
po::variables_map m_vm
 
- Protected Attributes inherited from Nektar::SolverUtils::Filter
LibUtilities::SessionReaderSharedPtr m_session
 
const std::weak_ptr< EquationSystemm_equ
 

Friends

class MemoryManager< FilterReynoldsStresses >
 

Additional Inherited Members

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

Detailed Description

Append Reynolds stresses to the average fields.

This class appends the average fields with the Reynolds stresses of the form \( \overline{u' v'} \).

For the default case, this is achieved by calculating \( C_{n} = \Sigma_{i=1}^{n} (u_i - \bar{u}_n)(v_i - \bar{v}_n)\) using the recursive relation:

\[ C_{n} = C_{n-1} + \frac{n}{n-1} (u_n - \bar{u}_n)(v_n - \bar{v}_n) \]

The FilterSampler base class then divides the result by n, leading to the Reynolds stress.

It is also possible to perform the averages using an exponential moving average, in which case either the moving average parameter \( \alpha \) or the time constant \( \tau \) must be prescribed.

Definition at line 42 of file FilterReynoldsStresses.h.

Constructor & Destructor Documentation

◆ FilterReynoldsStresses()

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

Definition at line 64 of file FilterReynoldsStresses.cpp.

68 : FilterFieldConvert(pSession, pEquation, pParams)
69{
70 // Load sampling frequency
71 auto it = pParams.find("SampleFrequency");
72 if (it == pParams.end())
73 {
75 }
76 else
77 {
78 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
79 m_sampleFrequency = round(equ.Evaluate());
80 }
81
82 // Check if should use moving average
83 it = pParams.find("MovingAverage");
84 if (it == pParams.end())
85 {
86 m_movAvg = false;
87 }
88 else
89 {
90 std::string sOption = it->second.c_str();
91 m_movAvg = (boost::iequals(sOption, "true")) ||
92 (boost::iequals(sOption, "yes"));
93 }
94
95 if (!m_movAvg)
96 {
97 return;
98 }
99
100 // Load alpha parameter for moving average
101 it = pParams.find("alpha");
102 if (it == pParams.end())
103 {
104 it = pParams.find("tau");
105 if (it == pParams.end())
106 {
107 ASSERTL0(false, "MovingAverage needs either alpha or tau.");
108 }
109 else
110 {
111 // Load time constant
112 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
113 NekDouble tau = equ.Evaluate();
114 // Load delta T between samples
115 NekDouble dT;
116 m_session->LoadParameter("TimeStep", dT);
117 dT = dT * m_sampleFrequency;
118 // Calculate alpha
119 m_alpha = dT / (tau + dT);
120 }
121 }
122 else
123 {
124 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
125 m_alpha = equ.Evaluate();
126 // Check if tau was also defined
127 it = pParams.find("tau");
128 if (it != pParams.end())
129 {
130 ASSERTL0(false,
131 "Cannot define both alpha and tau in MovingAverage.");
132 }
133 }
134 // Check bounds of m_alpha
135 ASSERTL0(m_alpha > 0 && m_alpha < 1, "Alpha out of bounds.");
136}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
SOLVER_UTILS_EXPORT FilterFieldConvert(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:83
double NekDouble

References ASSERTL0, Nektar::LibUtilities::Equation::Evaluate(), m_alpha, m_movAvg, Nektar::SolverUtils::FilterFieldConvert::m_sampleFrequency, and Nektar::SolverUtils::Filter::m_session.

◆ ~FilterReynoldsStresses()

Nektar::SolverUtils::FilterReynoldsStresses::~FilterReynoldsStresses ( )
override

Definition at line 138 of file FilterReynoldsStresses.cpp.

139{
140}

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 48 of file FilterReynoldsStresses.h.

52 {
55 pSession, pEquation, pParams);
56 return p;
57 }
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_FillVariablesName()

void Nektar::SolverUtils::FilterReynoldsStresses::v_FillVariablesName ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::FilterFieldConvert.

Definition at line 182 of file FilterReynoldsStresses.cpp.

184{
185 size_t dim = pFields.size() - 1;
186 size_t origFields = pFields.size();
187
188 // Fill name of variables
189 for (size_t n = 0; n < origFields; ++n)
190 {
191 m_variables.push_back(pFields[n]->GetSession()->GetVariable(n));
192 }
193 for (int i = 0; i < dim; ++i)
194 {
195 for (int j = i; j < dim; ++j)
196 {
197 std::string var = pFields[i]->GetSession()->GetVariable(i) +
198 pFields[j]->GetSession()->GetVariable(j);
199 m_variables.push_back(var);
200 }
201 }
202}

References Nektar::SolverUtils::FilterFieldConvert::m_variables.

◆ v_GetFileSuffix()

std::string Nektar::SolverUtils::FilterReynoldsStresses::v_GetFileSuffix ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::FilterFieldConvert.

Definition at line 84 of file FilterReynoldsStresses.h.

85 {
86 return "_stress";
87 }

◆ v_GetScale()

NekDouble Nektar::SolverUtils::FilterReynoldsStresses::v_GetScale ( )
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::FilterFieldConvert.

Definition at line 304 of file FilterReynoldsStresses.cpp.

305{
306 if (m_movAvg)
307 {
308 return 1.0;
309 }
310 else
311 {
312 return 1.0 / m_numSamples;
313 }
314}

References m_movAvg, and Nektar::SolverUtils::FilterFieldConvert::m_numSamples.

◆ v_Initialise()

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

Reimplemented from Nektar::SolverUtils::FilterFieldConvert.

Definition at line 142 of file FilterReynoldsStresses.cpp.

145{
146 size_t dim = pFields.size() - 1;
147 size_t nExtraFields = (dim + 1) * dim / 2;
148 size_t origFields = pFields.size();
149 size_t nqtot = pFields[0]->GetTotPoints();
150
151 // Allocate storage
152 m_fields.resize(origFields + nExtraFields);
153 m_delta.resize(dim);
154
155 for (size_t n = 0; n < m_fields.size(); ++n)
156 {
157 m_fields[n] = Array<OneD, NekDouble>(nqtot, 0.0);
158 }
159 for (size_t n = 0; n < m_delta.size(); ++n)
160 {
161 m_delta[n] = Array<OneD, NekDouble>(nqtot, 0.0);
162 }
163
164 // Initialise output arrays
166
167 // Update m_fields if using restart file
168 if (m_numSamples)
169 {
170 for (size_t j = 0; j < m_fields.size(); ++j)
171 {
172 pFields[0]->BwdTrans(m_outFields[j], m_fields[j]);
173 if (pFields[0]->GetWaveSpace())
174 {
175 pFields[0]->HomogeneousBwdTrans(nqtot, m_fields[j],
176 m_fields[j]);
177 }
178 }
179 }
180}
std::vector< Array< OneD, NekDouble > > m_outFields
SOLVER_UTILS_EXPORT void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
std::vector< Array< OneD, NekDouble > > m_delta
std::vector< Array< OneD, NekDouble > > m_fields

References m_delta, m_fields, Nektar::SolverUtils::FilterFieldConvert::m_numSamples, Nektar::SolverUtils::FilterFieldConvert::m_outFields, and Nektar::SolverUtils::FilterFieldConvert::v_Initialise().

◆ v_PrepareOutput()

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

Reimplemented from Nektar::SolverUtils::FilterFieldConvert.

Definition at line 279 of file FilterReynoldsStresses.cpp.

282{
283 size_t dim = pFields.size() - 1;
284
285 m_fieldMetaData["NumberOfFieldDumps"] = std::to_string(m_numSamples);
286
287 // Set wavespace to false, as calculations were performed in physical space
288 bool waveSpace = pFields[0]->GetWaveSpace();
289 pFields[0]->SetWaveSpace(false);
290
291 // Forward transform and put into m_outFields (except pressure)
292 for (size_t i = 0; i < m_fields.size(); ++i)
293 {
294 if (i != dim)
295 {
296 pFields[0]->FwdTransLocalElmt(m_fields[i], m_outFields[i]);
297 }
298 }
299
300 // Restore waveSpace
301 pFields[0]->SetWaveSpace(waveSpace);
302}
LibUtilities::FieldMetaDataMap m_fieldMetaData

References Nektar::SolverUtils::FilterFieldConvert::m_fieldMetaData, m_fields, Nektar::SolverUtils::FilterFieldConvert::m_numSamples, and Nektar::SolverUtils::FilterFieldConvert::m_outFields.

◆ v_ProcessSample()

void Nektar::SolverUtils::FilterReynoldsStresses::v_ProcessSample ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields,
std::vector< Array< OneD, NekDouble > > &  fieldcoeffs,
const NekDouble time 
)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::FilterFieldConvert.

Definition at line 204 of file FilterReynoldsStresses.cpp.

208{
209 size_t i, j, n;
210 size_t nq = pFields[0]->GetTotPoints();
211 size_t dim = pFields.size() - 1;
212 bool waveSpace = pFields[0]->GetWaveSpace();
213 NekDouble nSamples = (NekDouble)m_numSamples;
214
215 // For moving average, take first sample as initial vector
216 NekDouble alpha = m_alpha;
217 if (m_numSamples == 1)
218 {
219 alpha = 1.0;
220 }
221
222 // Define auxiliary constants for averages
223 NekDouble facOld, facAvg, facStress, facDelta;
224 if (m_movAvg)
225 {
226 facOld = 1.0 - alpha;
227 facAvg = alpha;
228 facStress = alpha;
229 facDelta = 1.0;
230 }
231 else
232 {
233 facOld = 1.0;
234 facAvg = 1.0;
235 facStress = nSamples / (nSamples - 1);
236 facDelta = 1.0 / nSamples;
237 }
238
239 Array<OneD, NekDouble> vel(nq);
240 Array<OneD, NekDouble> tmp(nq);
241
242 // Update original velocities in phys space and calculate (\bar{u} - u_n)
243 for (n = 0; n < dim; ++n)
244 {
245 if (waveSpace)
246 {
247 pFields[n]->HomogeneousBwdTrans(nq, pFields[n]->GetPhys(), vel);
248 }
249 else
250 {
251 vel = pFields[n]->GetPhys();
252 }
253 Vmath::Svtsvtp(nq, facAvg, vel, 1, facOld, m_fields[n], 1, m_fields[n],
254 1);
255 Vmath::Svtvm(nq, facDelta, m_fields[n], 1, vel, 1, m_delta[n], 1);
256 }
257 // Update pressure (directly to outFields)
258 Vmath::Svtsvtp(m_outFields[dim].size(), facAvg, pFields[dim]->GetCoeffs(),
259 1, facOld, m_outFields[dim], 1, m_outFields[dim], 1);
260
261 // Ignore Reynolds stress for first sample (its contribution is zero)
262 if (m_numSamples == 1)
263 {
264 return;
265 }
266
267 // Calculate C_{n} = facOld * C_{n-1} + facStress * deltaI * deltaJ
268 for (i = 0, n = dim + 1; i < dim; ++i)
269 {
270 for (j = i; j < dim; ++j, ++n)
271 {
272 Vmath::Vmul(nq, m_delta[i], 1, m_delta[j], 1, tmp, 1);
273 Vmath::Svtsvtp(nq, facStress, tmp, 1, facOld, m_fields[n], 1,
274 m_fields[n], 1);
275 }
276 }
277}
void Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
Svtsvtp (scalar times vector plus scalar times vector):
Definition: Vmath.hpp:473
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Svtvm(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvm (scalar times vector minus vector): z = alpha*x - y.
Definition: Vmath.hpp:424

References m_alpha, m_delta, m_fields, m_movAvg, Nektar::SolverUtils::FilterFieldConvert::m_numSamples, Nektar::SolverUtils::FilterFieldConvert::m_outFields, Vmath::Svtsvtp(), Vmath::Svtvm(), and Vmath::Vmul().

Friends And Related Function Documentation

◆ MemoryManager< FilterReynoldsStresses >

friend class MemoryManager< FilterReynoldsStresses >
friend

Definition at line 1 of file FilterReynoldsStresses.h.

Member Data Documentation

◆ className

std::string Nektar::SolverUtils::FilterReynoldsStresses::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< SolverUtils::EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
Creates an instance of this class.
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:39

Name of the class.

Definition at line 60 of file FilterReynoldsStresses.h.

◆ m_alpha

NekDouble Nektar::SolverUtils::FilterReynoldsStresses::m_alpha
protected

Definition at line 91 of file FilterReynoldsStresses.h.

Referenced by FilterReynoldsStresses(), and v_ProcessSample().

◆ m_delta

std::vector<Array<OneD, NekDouble> > Nektar::SolverUtils::FilterReynoldsStresses::m_delta
protected

Definition at line 90 of file FilterReynoldsStresses.h.

Referenced by v_Initialise(), and v_ProcessSample().

◆ m_fields

std::vector<Array<OneD, NekDouble> > Nektar::SolverUtils::FilterReynoldsStresses::m_fields
protected

Definition at line 89 of file FilterReynoldsStresses.h.

Referenced by v_Initialise(), v_PrepareOutput(), and v_ProcessSample().

◆ m_movAvg

bool Nektar::SolverUtils::FilterReynoldsStresses::m_movAvg
protected

Definition at line 92 of file FilterReynoldsStresses.h.

Referenced by FilterReynoldsStresses(), v_GetScale(), and v_ProcessSample().