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

#include <FilterBodyFittedVelocity.h>

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

Public Member Functions

SOLVER_UTILS_EXPORT FilterBodyFittedVelocity (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
 
SOLVER_UTILS_EXPORT ~FilterBodyFittedVelocity () 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< 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 Types

enum  FilterType { eOriginal , eMax , eMin }
 
enum  ProblemType { eCompressible , eIncompressible , eOthers }
 

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::string m_bodyFittedCooriateFile
 
FilterType m_filterType
 
ProblemType m_problemType
 
unsigned int m_spaceDim
 
unsigned int m_nFields
 
unsigned int m_nAddFields
 
unsigned int m_nVars
 
std::vector< std::string > m_bfsVars
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_bfcsDir
 
std::vector< Array< OneD, NekDouble > > m_curFieldsVels_Car
 
std::vector< Array< OneD, NekDouble > > m_curFieldsVels
 
std::vector< Array< OneD, NekDouble > > m_outFieldsVels
 
std::vector< Array< OneD, NekDouble > > m_curFieldsThermalVars
 
std::vector< Array< OneD, NekDouble > > m_outFieldsThermalVars
 
- 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
 

Private Attributes

bool m_initialized
 

Friends

class MemoryManager< FilterBodyFittedVelocity >
 

Additional Inherited Members

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

Detailed Description

Definition at line 44 of file FilterBodyFittedVelocity.h.

Member Enumeration Documentation

◆ FilterType

◆ ProblemType

Constructor & Destructor Documentation

◆ FilterBodyFittedVelocity()

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

Definition at line 54 of file FilterBodyFittedVelocity.cpp.

57 : FilterFieldConvert(pSession, pEquation, pParams)
58{
59 // Load sampling frequency
60 auto it = pParams.find("SampleFrequency");
61 if (it == pParams.end())
62 {
64 }
65 else
66 {
67 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
68 m_sampleFrequency = round(equ.Evaluate());
69 }
70
71 // Load flag for filter type
72 it = pParams.find("OriginalOrMaxOrMin");
73 if (it == pParams.end())
74 {
76 }
77 else
78 {
79 std::string sOptionType = it->second.c_str();
80 if (boost::iequals(sOptionType, "maximum") ||
81 boost::iequals(sOptionType, "max"))
82 {
84 }
85 else if (boost::iequals(sOptionType, "minumun") ||
86 boost::iequals(sOptionType, "min"))
87 {
89 }
90 else if (boost::iequals(sOptionType, "original"))
91 {
93 }
94 else
95 {
96 WARNINGL0(false, "Detailed filter type is not found, use original");
98 }
99 }
100
101 // Load bodyFittedCoordinate from a file generated by FieldConvert -m
102 // bodyFittedVelocity
103 it = pParams.find("BodyFittedCoordinateFile");
104 ASSERTL0(it->second.length() > 0,
105 "Missing parameter 'BodyFittedCoordinateFile'.");
106
107 if (it->second.find_last_of('.') != std::string::npos)
108 {
109
110 m_bodyFittedCooriateFile = it->second;
111 }
112 else
113 {
114 std::stringstream outname;
115 outname << it->second << ".fld";
116 m_bodyFittedCooriateFile = outname.str();
117 }
118
120}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:215
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

References ASSERTL0, eMax, eMin, eOriginal, Nektar::LibUtilities::Equation::Evaluate(), m_bodyFittedCooriateFile, m_filterType, m_initialized, Nektar::SolverUtils::FilterFieldConvert::m_restartFile, Nektar::SolverUtils::FilterFieldConvert::m_sampleFrequency, Nektar::SolverUtils::Filter::m_session, and WARNINGL0.

◆ ~FilterBodyFittedVelocity()

Nektar::SolverUtils::FilterBodyFittedVelocity::~FilterBodyFittedVelocity ( )
override

Definition at line 122 of file FilterBodyFittedVelocity.cpp.

123{
124}

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 50 of file FilterBodyFittedVelocity.h.

54 {
57 pSession, pEquation, pParams);
58 return p;
59 }
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::FilterBodyFittedVelocity::v_FillVariablesName ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pFields)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::FilterFieldConvert.

Definition at line 275 of file FilterBodyFittedVelocity.cpp.

277{
278 // Fill in the var names using the routine in base class
280
281 // Check type of problem
282 std::vector<std::string> varNames =
283 pFields[0]->GetSession()->GetVariables();
284 if (boost::iequals(varNames[0], "u") && boost::iequals(varNames[1], "v"))
285 {
287 m_spaceDim = pFields.size() - 1;
288 }
289 else if (boost::iequals(varNames[0], "rho") &&
290 boost::iequals(varNames[1], "rhou"))
291 {
293 m_spaceDim = pFields.size() - 2;
294 }
295 else
296 {
297 WARNINGL0(false, "Problem type is not claear. Please check");
299 m_spaceDim = 2;
300 }
301
302 m_nVars = m_variables.size(); // Include extra vars
303 m_nFields = pFields.size(); // Conservative vars for cfs
304 m_nAddFields = 3;
305
306 // Fill in the body-fitted velocity names
307 // bfc represents body-fitted coordinate
308 m_variables.push_back("distanceToWall");
309 m_variables.push_back("u_bfc");
310 m_variables.push_back("v_bfc");
311 if (m_spaceDim == 3)
312 {
313 m_variables.push_back("w_bfc");
314 ++m_nAddFields;
315 }
316 else if (m_spaceDim != 2)
317 {
318 ASSERTL0(false, "Unsupported dimension");
319 }
320}
virtual SOLVER_UTILS_EXPORT void v_FillVariablesName(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)

References ASSERTL0, eCompressible, eIncompressible, eOthers, m_nAddFields, m_nFields, m_nVars, m_problemType, m_spaceDim, Nektar::SolverUtils::FilterFieldConvert::m_variables, Nektar::SolverUtils::FilterFieldConvert::v_FillVariablesName(), and WARNINGL0.

◆ v_GetFileSuffix()

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

Reimplemented from Nektar::SolverUtils::FilterFieldConvert.

Definition at line 118 of file FilterBodyFittedVelocity.h.

119 {
120 if (m_filterType == eMax)
121 {
122 return "_max";
123 }
124 else if (m_filterType == eMin)
125 {
126 return "_min";
127 }
128 else
129 {
130 return "_original";
131 }
132 }

References eMax, eMin, and m_filterType.

◆ v_GetScale()

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

Reimplemented from Nektar::SolverUtils::FilterFieldConvert.

Definition at line 494 of file FilterBodyFittedVelocity.cpp.

495{
496 return 1.0;
497}

◆ v_Initialise()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 126 of file FilterBodyFittedVelocity.cpp.

129{
130 // Initialise output arrays
132
133 // Generate the body-fitted coordinate system and distance field
134 // bfcsDir[i][j][k]
135 // i - dir vec: 0-main tangential, 1-normal, (2-minor tangential)
136 // j - component: 0-x, 1-y, (2-z)
137 // k - point id
138 const int npoints = pFields[0]->GetTotPoints();
139 m_bfcsDir = Array<OneD, Array<OneD, Array<OneD, NekDouble>>>(m_spaceDim);
140 for (int i = 0; i < m_spaceDim; ++i)
141 {
142 m_bfcsDir[i] = Array<OneD, Array<OneD, NekDouble>>(m_spaceDim);
143 for (int j = 0; j < m_spaceDim; ++j)
144 {
145 m_bfcsDir[i][j] = Array<OneD, NekDouble>(npoints);
146 }
147 }
148
149 // Initialize the body fitted coordinate in the file
150 // m_bfsVars holds the var names for the coordinate system and will be
151 // used to find corresponding field to load from the input file.
152 // Ref: FilterFieldConvert.cpp (v_Initialise) & ProcessInterpField.cpp
153 // (Process)
154 m_bfsVars.resize(m_spaceDim * m_spaceDim); // 4 or 9
155 if (m_spaceDim == 2)
156 {
157 m_bfsVars[0] = "bfc_dir0_x";
158 m_bfsVars[1] = "bfc_dir0_y";
159 m_bfsVars[2] = "bfc_dir1_x";
160 m_bfsVars[3] = "bfc_dir1_y";
161 }
162 else
163 {
164 m_bfsVars[0] = "bfc_dir0_x";
165 m_bfsVars[1] = "bfc_dir0_y";
166 m_bfsVars[2] = "bfc_dir0_z";
167 m_bfsVars[3] = "bfc_dir1_x";
168 m_bfsVars[4] = "bfc_dir1_y";
169 m_bfsVars[5] = "bfc_dir1_z";
170 m_bfsVars[6] = "bfc_dir2_x";
171 m_bfsVars[7] = "bfc_dir2_y";
172 m_bfsVars[8] = "bfc_dir2_z";
173 }
174
175 // Load file
176 std::vector<LibUtilities::FieldDefinitionsSharedPtr> fieldDef;
177 std::vector<std::vector<NekDouble>> fieldData;
178 LibUtilities::FieldMetaDataMap fieldMetaData;
181 fld->Import(m_bodyFittedCooriateFile, fieldDef, fieldData, fieldMetaData);
182
183 // Get a tmp array to hold the coeffs from the input file
184 Array<OneD, NekDouble> tmp_coeffs;
185 tmp_coeffs = Array<OneD, NekDouble>(pFields[0]->GetNcoeffs(), 0.0);
186
187 // Extract data and convert to physical space
188 for (int k = 0; k < m_bfsVars.size(); ++k)
189 {
190 // In the fieldDef, for example fieldDef[0] is the quad and
191 // fieldDef[1] is the tri region
192 // Each of fieldDef[i] contains all the fields we expected in an fld
193 for (int i = 0; i < fieldData.size(); ++i)
194 {
195 pFields[0]->ExtractDataToCoeffs(fieldDef[i], fieldData[i],
196 m_bfsVars[k], tmp_coeffs);
197 }
198
199 // Bwd transorm the data and put them into m_bfcsDir[][]
200 int ii, jj;
201 jj = k % m_spaceDim;
202 ii = k / m_spaceDim;
203
204 pFields[0]->BwdTrans(tmp_coeffs, m_bfcsDir[ii][jj]);
205 }
206
207 // Allocate storage
209 m_curFieldsVels.resize(m_spaceDim); // m_spaceDim = m_nAddFields-1
211
212 for (int n = 0; n < m_spaceDim; ++n)
213 {
214 m_curFieldsVels_Car[n] = Array<OneD, NekDouble>(npoints, 0.0);
215 m_curFieldsVels[n] = Array<OneD, NekDouble>(npoints, 0.0);
216 m_outFieldsVels[n] = Array<OneD, NekDouble>(npoints, 0.0);
217 }
218
219 // Allocate storage for rho,p,T
221 {
222 m_curFieldsThermalVars.resize(3);
223 m_outFieldsThermalVars.resize(3);
224
225 for (int n = 0; n < 3; ++n)
226 {
227 m_curFieldsThermalVars[n] = Array<OneD, NekDouble>(npoints, 0.0);
228 m_outFieldsThermalVars[n] = Array<OneD, NekDouble>(npoints, 0.0);
229 }
230 }
231
232 // m_outFields contains the initialization values
233 // 2d: rho,rhou,rhov,E,u,v,p,T,s,a,Mach,Sensor,distanceToWall,u_bfc,v_bfc
234 if (m_initialized)
235 {
236 for (int n = 0; n < m_spaceDim; ++n)
237 {
238 pFields[0]->BwdTrans(m_outFields[m_nVars + 1 + n],
239 m_outFieldsVels[n]);
240 if (pFields[0]->GetWaveSpace())
241 {
242 pFields[0]->HomogeneousBwdTrans(npoints, m_outFieldsVels[n],
243 m_outFieldsVels[n]);
244 }
245 }
246
248 {
249 // For cfs, initialize thermal vars
250 // shift_thermal for [rho, p, T] in m_outFields
251 int shift_thermal[3] = {0, m_spaceDim * 2 + 2, m_spaceDim * 2 + 3};
252
253 for (int n = 0; n < 3; ++n)
254 {
255 pFields[0]->BwdTrans(m_outFields[shift_thermal[n]],
257 if (pFields[0]->GetWaveSpace())
258 {
259 pFields[0]->HomogeneousBwdTrans(npoints,
262 }
263 }
264 }
265 }
266 else
267 {
268 ASSERTL0(false, "Restart file is expectd. It can be generated by \
269 FieldConvert in serial. For example: \
270 FieldConvert -m bodyFittedVelocity:bnd=0:bfcOut=1 \
271 mesh.xml session.xml baseflow.fld bfvFile.fld");
272 }
273}
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
Definition: FieldIO.cpp:224
std::vector< Array< OneD, NekDouble > > m_outFieldsVels
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_bfcsDir
std::vector< Array< OneD, NekDouble > > m_curFieldsVels
std::vector< Array< OneD, NekDouble > > m_curFieldsThermalVars
std::vector< Array< OneD, NekDouble > > m_curFieldsVels_Car
std::vector< Array< OneD, NekDouble > > m_outFieldsThermalVars
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::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:322
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:50

References ASSERTL0, Nektar::LibUtilities::FieldIO::CreateForFile(), eCompressible, m_bfcsDir, m_bfsVars, m_bodyFittedCooriateFile, m_curFieldsThermalVars, m_curFieldsVels, m_curFieldsVels_Car, m_initialized, m_nVars, Nektar::SolverUtils::FilterFieldConvert::m_outFields, m_outFieldsThermalVars, m_outFieldsVels, m_problemType, Nektar::SolverUtils::Filter::m_session, m_spaceDim, and Nektar::SolverUtils::FilterFieldConvert::v_Initialise().

◆ v_PrepareOutput()

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

◆ v_ProcessSample()

void Nektar::SolverUtils::FilterBodyFittedVelocity::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 338 of file FilterBodyFittedVelocity.cpp.

342{
343 // Use shift_vel to get the current u,v,w in Cartesian coordinate
344 // cfs_2D:
345 // rho,rhou,rhov,E,u,v,p,T,s,a,Mach,Sensor,distanceToWall,u_bfc,v_bfc
346 // inc_2D: u,v,p
347 int shift_vel = (m_problemType == eCompressible) ? (m_spaceDim + 2) : 0;
348
349 // Get current u,v,w in Cartesian coordinate
350 for (int n = 0; n < m_spaceDim; ++n)
351 {
352 pFields[0]->BwdTrans(fieldcoeffs[shift_vel + n],
354 if (pFields[0]->GetWaveSpace())
355 {
356 pFields[0]->HomogeneousBwdTrans(pFields[0]->GetNpoints(),
359 }
360 }
361
362 // Project the velocity into the body-fitted coordinate system
363 int npoints = pFields[0]->GetTotPoints();
364 Array<OneD, NekDouble> vel_tmp(npoints);
365
366 for (int i = 0; i < m_spaceDim; ++i) // loop for bfc velocity
367 {
368 Vmath::Zero(npoints, m_curFieldsVels[i], 1);
369
370 for (int j = 0; j < m_spaceDim; ++j)
371 {
372 Vmath::Vmul(npoints, m_curFieldsVels_Car[j], 1, m_bfcsDir[i][j], 1,
373 vel_tmp, 1);
374 Vmath::Vadd(npoints, vel_tmp, 1, m_curFieldsVels[i], 1,
375 m_curFieldsVels[i], 1);
376 }
377 }
378
379 // Get max/min/original for the u/v/w_bfc field
380 for (int n = 0; n < m_spaceDim; ++n)
381 {
382 size_t length = m_outFieldsVels[n].size();
383
384 if (m_filterType == eMax)
385 {
386 // Compute max
387 for (int i = 0; i < length; ++i)
388 {
389 if (m_curFieldsVels[n][i] > m_outFieldsVels[n][i])
390 {
391 m_outFieldsVels[n][i] = m_curFieldsVels[n][i];
392 }
393 }
394 }
395 else if (m_filterType == eMin)
396 {
397 // Compute min
398 for (int i = 0; i < length; ++i)
399 {
400 if (m_curFieldsVels[n][i] < m_outFieldsVels[n][i])
401 {
402 m_outFieldsVels[n][i] = m_curFieldsVels[n][i];
403 }
404 }
405 }
406 else
407 {
408 // Original field
410 }
411 }
412
413 // Forward transform and put into m_outFields
414 for (int n = 0; n < m_spaceDim; ++n)
415 {
416 pFields[0]->FwdTransLocalElmt(m_outFieldsVels[n],
417 m_outFields[m_nVars + 1 + n]);
418 }
419
420 // Process the thermal variables for compressible flows
422 {
423 // Get current rho,p,T fields
424 // shift_thermal for [rho, p, T] in fieldcoeffs
425 int shift_thermal[3] = {0, m_spaceDim * 2 + 2, m_spaceDim * 2 + 3};
426
427 for (int n = 0; n < 3; ++n)
428 {
429 pFields[0]->BwdTrans(fieldcoeffs[shift_thermal[n]],
431 if (pFields[0]->GetWaveSpace())
432 {
433 pFields[0]->HomogeneousBwdTrans(pFields[0]->GetNpoints(),
436 }
437 }
438
439 // Get max/min/original for the rho/p/T fields
440 for (int n = 0; n < 3; ++n)
441 {
442 size_t length = m_outFieldsThermalVars[n].size();
443
444 if (m_filterType == eMax)
445 {
446 // Compute max
447 for (int i = 0; i < length; ++i)
448 {
449 if (m_curFieldsThermalVars[n][i] >
451 {
454 }
455 }
456 }
457 else if (m_filterType == eMin)
458 {
459 // Compute min
460 for (int i = 0; i < length; ++i)
461 {
462 if (m_curFieldsThermalVars[n][i] <
464 {
467 }
468 }
469 }
470 else
471 {
472 // Original field
474 }
475 }
476
477 // Forward transform and put into m_outFields
478 for (int n = 0; n < 3; ++n)
479 {
480 pFields[0]->FwdTransLocalElmt(m_outFieldsThermalVars[n],
481 m_outFields[shift_thermal[n]]);
482 }
483 }
484}
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 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 Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

References eCompressible, eMax, eMin, m_bfcsDir, m_curFieldsThermalVars, m_curFieldsVels, m_curFieldsVels_Car, m_filterType, m_nVars, Nektar::SolverUtils::FilterFieldConvert::m_outFields, m_outFieldsThermalVars, m_outFieldsVels, m_problemType, m_spaceDim, Vmath::Vadd(), Vmath::Vmul(), and Vmath::Zero().

Friends And Related Function Documentation

◆ MemoryManager< FilterBodyFittedVelocity >

friend class MemoryManager< FilterBodyFittedVelocity >
friend

Definition at line 1 of file FilterBodyFittedVelocity.h.

Member Data Documentation

◆ className

std::string Nektar::SolverUtils::FilterBodyFittedVelocity::className
static
Initial value:
=
"BodyFittedVelocity", FilterBodyFittedVelocity::create)
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)
Creates an instance of this class.
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:39

Name of the class.

Definition at line 62 of file FilterBodyFittedVelocity.h.

◆ m_bfcsDir

Array<OneD, Array<OneD, Array<OneD, NekDouble> > > Nektar::SolverUtils::FilterBodyFittedVelocity::m_bfcsDir
protected

Definition at line 95 of file FilterBodyFittedVelocity.h.

Referenced by v_Initialise(), and v_ProcessSample().

◆ m_bfsVars

std::vector<std::string> Nektar::SolverUtils::FilterBodyFittedVelocity::m_bfsVars
protected

Definition at line 94 of file FilterBodyFittedVelocity.h.

Referenced by v_Initialise().

◆ m_bodyFittedCooriateFile

std::string Nektar::SolverUtils::FilterBodyFittedVelocity::m_bodyFittedCooriateFile
protected

Definition at line 85 of file FilterBodyFittedVelocity.h.

Referenced by FilterBodyFittedVelocity(), and v_Initialise().

◆ m_curFieldsThermalVars

std::vector<Array<OneD, NekDouble> > Nektar::SolverUtils::FilterBodyFittedVelocity::m_curFieldsThermalVars
protected

Definition at line 101 of file FilterBodyFittedVelocity.h.

Referenced by v_Initialise(), and v_ProcessSample().

◆ m_curFieldsVels

std::vector<Array<OneD, NekDouble> > Nektar::SolverUtils::FilterBodyFittedVelocity::m_curFieldsVels
protected

Definition at line 98 of file FilterBodyFittedVelocity.h.

Referenced by v_Initialise(), and v_ProcessSample().

◆ m_curFieldsVels_Car

std::vector<Array<OneD, NekDouble> > Nektar::SolverUtils::FilterBodyFittedVelocity::m_curFieldsVels_Car
protected

Definition at line 97 of file FilterBodyFittedVelocity.h.

Referenced by v_Initialise(), and v_ProcessSample().

◆ m_filterType

FilterType Nektar::SolverUtils::FilterBodyFittedVelocity::m_filterType
protected

◆ m_initialized

bool Nektar::SolverUtils::FilterBodyFittedVelocity::m_initialized
private

Definition at line 135 of file FilterBodyFittedVelocity.h.

Referenced by FilterBodyFittedVelocity(), and v_Initialise().

◆ m_nAddFields

unsigned int Nektar::SolverUtils::FilterBodyFittedVelocity::m_nAddFields
protected

Definition at line 91 of file FilterBodyFittedVelocity.h.

Referenced by v_FillVariablesName().

◆ m_nFields

unsigned int Nektar::SolverUtils::FilterBodyFittedVelocity::m_nFields
protected

Definition at line 90 of file FilterBodyFittedVelocity.h.

Referenced by v_FillVariablesName().

◆ m_nVars

unsigned int Nektar::SolverUtils::FilterBodyFittedVelocity::m_nVars
protected

Definition at line 92 of file FilterBodyFittedVelocity.h.

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

◆ m_outFieldsThermalVars

std::vector<Array<OneD, NekDouble> > Nektar::SolverUtils::FilterBodyFittedVelocity::m_outFieldsThermalVars
protected

Definition at line 102 of file FilterBodyFittedVelocity.h.

Referenced by v_Initialise(), and v_ProcessSample().

◆ m_outFieldsVels

std::vector<Array<OneD, NekDouble> > Nektar::SolverUtils::FilterBodyFittedVelocity::m_outFieldsVels
protected

Definition at line 99 of file FilterBodyFittedVelocity.h.

Referenced by v_Initialise(), and v_ProcessSample().

◆ m_problemType

ProblemType Nektar::SolverUtils::FilterBodyFittedVelocity::m_problemType
protected

Definition at line 88 of file FilterBodyFittedVelocity.h.

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

◆ m_spaceDim

unsigned int Nektar::SolverUtils::FilterBodyFittedVelocity::m_spaceDim
protected

Definition at line 89 of file FilterBodyFittedVelocity.h.

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