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

#include <FilterAeroForcesSPM.h>

Inheritance diagram for Nektar::FilterAeroForcesSPM:
[legend]

Public Member Functions

void CalculateForces (const Array< OneD, Array< OneD, NekDouble > > &pIntVel, const Array< OneD, Array< OneD, NekDouble > > &pUpPrev, const MultiRegions::ExpListSharedPtr &pPhi, NekDouble time, NekDouble dt)
 Determine the total force on the body defined by \(\Phi\) (note that if the shape function represents more than one body, this function calculates the value of the final force after adding up the values for each body). This value must be scaled with the density to get the real force vector. More...
 
- 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 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...
 

Protected Member Functions

 FilterAeroForcesSPM (const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_ptr< SolverUtils::EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
 
 ~FilterAeroForcesSPM () override=default
 
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)
 

Protected Attributes

unsigned int m_index
 
unsigned int m_outputFrequency
 
std::string m_outputFile
 
std::ofstream m_outputStream
 
NekDouble m_startTime
 
std::vector< std::string > m_dirNames
 STL vector containing the names of the different directions. More...
 
NekDouble m_spaceDim
 Dimension of the fluid domain. More...
 
Array< OneD, NekDoublem_Forces
 Array storing the last value of the aerodynamic forces. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::Filter
LibUtilities::SessionReaderSharedPtr m_session
 
const std::weak_ptr< EquationSystemm_equ
 
bool m_updateOnInitialise = true
 

Friends

class MemoryManager< FilterAeroForcesSPM >
 

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 FilterAeroForcesSPM.h.

Constructor & Destructor Documentation

◆ FilterAeroForcesSPM()

Nektar::FilterAeroForcesSPM::FilterAeroForcesSPM ( const LibUtilities::SessionReaderSharedPtr pSession,
const std::shared_ptr< SolverUtils::EquationSystem > &  pEquation,
const std::map< std::string, std::string > &  pParams 
)
protected

Definition at line 48 of file FilterAeroForcesSPM.cpp.

52 : Filter(pSession, pEquation)
53{
54 // OutputFile
55 std::string ext = ".fce";
56 m_outputFile = Filter::SetupOutput(ext, pParams);
57
58 // OutputFrequency
59 auto it = pParams.find("OutputFrequency");
60 if (it == pParams.end())
61 {
63 }
64 else
65 {
66 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
67 m_outputFrequency = round(equ.Evaluate());
68 }
69
70 // Time after which we need to calculate the forces
71 it = pParams.find("StartTime");
72 if (it == pParams.end())
73 {
74 m_startTime = 0;
75 }
76 else
77 {
78 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
79 m_startTime = equ.Evaluate();
80 }
81}
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 Nektar::LibUtilities::Equation::Evaluate(), m_outputFile, m_outputFrequency, Nektar::SolverUtils::Filter::m_session, m_startTime, and Nektar::SolverUtils::Filter::SetupOutput().

◆ ~FilterAeroForcesSPM()

Nektar::FilterAeroForcesSPM::~FilterAeroForcesSPM ( )
overrideprotecteddefault

Member Function Documentation

◆ CalculateForces()

void Nektar::FilterAeroForcesSPM::CalculateForces ( const Array< OneD, Array< OneD, NekDouble > > &  pIntVel,
const Array< OneD, Array< OneD, NekDouble > > &  pUpPrev,
const MultiRegions::ExpListSharedPtr pPhi,
NekDouble  time,
NekDouble  dt 
)

Determine the total force on the body defined by \(\Phi\) (note that if the shape function represents more than one body, this function calculates the value of the final force after adding up the values for each body). This value must be scaled with the density to get the real force vector.

Parameters
pIntVel
pUpPrev
pPhi
time
dt

Definition at line 205 of file FilterAeroForcesSPM.cpp.

209{
210 // Only output every m_outputFrequency.
211 if ((m_index % m_outputFrequency) || (time < m_startTime))
212 {
213 return;
214 }
215
216 int nq = pIntVel[0].size();
217 Array<OneD, NekDouble> tmp(nq);
218
219 // Aerodynamic forces are computed according eq. 18a in Luo et al. (2009).
220 // Smoothed profile method for particulate flows: Error analysis and
221 // simulations. Journal of Computational Physics, 228(5)
222 for (int i = 0; i < m_spaceDim; ++i)
223 {
224 // "Scalar" field to be integrated
225 Vmath::Vsub(nq, pIntVel[i], 1, pUpPrev[i], 1, tmp, 1);
226 Vmath::Vmul(nq, pPhi->GetPhys(), 1, tmp, 1, tmp, 1);
227
228 // Integration of force throughout the domain
229 m_Forces[i] = pPhi->Integral(tmp) / dt;
230 }
231}
NekDouble m_spaceDim
Dimension of the fluid domain.
Array< OneD, NekDouble > m_Forces
Array storing the last value of the aerodynamic forces.
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 Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220

References m_Forces, m_index, m_outputFrequency, m_spaceDim, m_startTime, Vmath::Vmul(), and Vmath::Vsub().

◆ create()

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

Creates an instance of this class.

Definition at line 50 of file FilterAeroForcesSPM.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:52

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

◆ v_Finalise()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 174 of file FilterAeroForcesSPM.cpp.

177{
178 if (pFields[0]->GetComm()->GetRank() == 0)
179 {
180 m_outputStream.close();
181 }
182}

References m_outputStream.

◆ v_Initialise()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 86 of file FilterAeroForcesSPM.cpp.

89{
90 // Save space dimension
91 m_spaceDim = pFields[0]->GetGraph()->GetMeshDimension();
92
93 // Fill in the directions vector
94 m_dirNames.push_back("X");
95 if (m_spaceDim > 1)
96 {
97 m_dirNames.push_back("Y");
98 }
99 if (m_spaceDim > 2)
100 {
101 m_dirNames.push_back("Z");
102 }
103
104 // Allocate the aero-forces vector
105 m_Forces = Array<OneD, NekDouble>(m_spaceDim);
106
107 // Write header
108 LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
109 if (vComm->GetRank() == 0)
110 {
111 // Open output stream
112 bool adaptive;
113 m_session->MatchSolverInfo("Driver", "Adaptive", adaptive, false);
114 if (adaptive)
115 {
116 m_outputStream.open(m_outputFile.c_str(), std::ofstream::app);
117 }
118 else
119 {
120 m_outputStream.open(m_outputFile.c_str());
121 }
122 m_outputStream << "# Forces acting on bodies\n";
123 m_outputStream << "#";
124 m_outputStream.width(7);
125 m_outputStream << "Time";
126 for (int i = 0; i < m_spaceDim; ++i)
127 {
128 m_outputStream.width(14);
129 m_outputStream << "F_" << m_dirNames[i];
130 }
131
132 m_outputStream << std::endl;
133 }
134
135 m_index = 0;
136}
std::vector< std::string > m_dirNames
STL vector containing the names of the different directions.
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55

References m_dirNames, m_Forces, m_index, m_outputFile, m_outputStream, Nektar::SolverUtils::Filter::m_session, and m_spaceDim.

◆ v_IsTimeDependent()

bool Nektar::FilterAeroForcesSPM::v_IsTimeDependent ( )
overrideprotectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 187 of file FilterAeroForcesSPM.cpp.

188{
189 return true;
190}

◆ v_Update()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 141 of file FilterAeroForcesSPM.cpp.

144{
145 // Only output every m_outputFrequency.
146 if ((m_index++) % m_outputFrequency || (time < m_startTime))
147 {
148 return;
149 }
150
151 // Communicators to exchange results
152 LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
153
154 // Write Results
155 if (vComm->GetRank() == 0)
156 {
157 // Write time
158 m_outputStream.width(8);
159 m_outputStream << std::setprecision(6) << time;
160 // Write forces
161 for (int i = 0; i < m_spaceDim; ++i)
162 {
163 m_outputStream.width(15);
164 m_outputStream << std::setprecision(8) << m_Forces[i];
165 }
166 m_outputStream.width(10);
167 m_outputStream << std::endl;
168 }
169}

References m_Forces, m_index, m_outputFrequency, m_outputStream, m_spaceDim, and m_startTime.

Friends And Related Function Documentation

◆ MemoryManager< FilterAeroForcesSPM >

friend class MemoryManager< FilterAeroForcesSPM >
friend

Definition at line 108 of file FilterAeroForcesSPM.h.

Member Data Documentation

◆ className

std::string Nektar::FilterAeroForcesSPM::className
static
Initial value:
=
"AeroForcesSPM", FilterAeroForcesSPM::create)
static SolverUtils::FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_ptr< SolverUtils::EquationSystem > &pEquation, const std::map< std::string, std::string > &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 62 of file FilterAeroForcesSPM.h.

◆ m_dirNames

std::vector<std::string> Nektar::FilterAeroForcesSPM::m_dirNames
protected

STL vector containing the names of the different directions.

Definition at line 78 of file FilterAeroForcesSPM.h.

Referenced by v_Initialise().

◆ m_Forces

Array<OneD, NekDouble> Nektar::FilterAeroForcesSPM::m_Forces
protected

Array storing the last value of the aerodynamic forces.

Definition at line 82 of file FilterAeroForcesSPM.h.

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

◆ m_index

unsigned int Nektar::FilterAeroForcesSPM::m_index
protected

Definition at line 71 of file FilterAeroForcesSPM.h.

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

◆ m_outputFile

std::string Nektar::FilterAeroForcesSPM::m_outputFile
protected

Definition at line 73 of file FilterAeroForcesSPM.h.

Referenced by FilterAeroForcesSPM(), and v_Initialise().

◆ m_outputFrequency

unsigned int Nektar::FilterAeroForcesSPM::m_outputFrequency
protected

Definition at line 72 of file FilterAeroForcesSPM.h.

Referenced by CalculateForces(), FilterAeroForcesSPM(), and v_Update().

◆ m_outputStream

std::ofstream Nektar::FilterAeroForcesSPM::m_outputStream
protected

Definition at line 74 of file FilterAeroForcesSPM.h.

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

◆ m_spaceDim

NekDouble Nektar::FilterAeroForcesSPM::m_spaceDim
protected

Dimension of the fluid domain.

Definition at line 80 of file FilterAeroForcesSPM.h.

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

◆ m_startTime

NekDouble Nektar::FilterAeroForcesSPM::m_startTime
protected

Definition at line 76 of file FilterAeroForcesSPM.h.

Referenced by CalculateForces(), FilterAeroForcesSPM(), and v_Update().