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

#include <FilterAeroForcesSPM.h>

Inheritance diagram for Nektar::FilterAeroForcesSPM:
[legend]

Public Member Functions

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

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
 

Additional Inherited Members

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

Detailed Description

Definition at line 43 of file FilterAeroForcesSPM.h.

Constructor & Destructor Documentation

◆ FilterAeroForcesSPM()

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

Definition at line 49 of file FilterAeroForcesSPM.cpp.

53 : Filter(pSession, pEquation)
54{
55 // OutputFile
56 auto it = pParams.find("OutputFile");
57 if (it == pParams.end())
58 {
59 m_outputFile = m_session->GetSessionName();
60 }
61 else
62 {
63 ASSERTL0(it->second.length() > 0, "Missing parameter 'OutputFile'.");
64 m_outputFile = it->second;
65 }
66 if (!(m_outputFile.length() >= 4 &&
67 m_outputFile.substr(m_outputFile.length() - 4) == ".fce"))
68 {
69 m_outputFile += ".fce";
70 }
71
72 // OutputFrequency
73 it = pParams.find("OutputFrequency");
74 if (it == pParams.end())
75 {
77 }
78 else
79 {
80 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
81 m_outputFrequency = round(equ.Evaluate());
82 }
83
84 // Time after which we need to calculate the forces
85 it = pParams.find("StartTime");
86 if (it == pParams.end())
87 {
88 m_startTime = 0;
89 }
90 else
91 {
92 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
93 m_startTime = equ.Evaluate();
94 }
95}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:83
SOLVER_UTILS_EXPORT Filter(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation)
Definition: Filter.cpp:45

References ASSERTL0, Nektar::LibUtilities::Equation::Evaluate(), m_outputFile, m_outputFrequency, Nektar::SolverUtils::Filter::m_session, and m_startTime.

◆ ~FilterAeroForcesSPM()

Nektar::FilterAeroForcesSPM::~FilterAeroForcesSPM ( )
override

Definition at line 100 of file FilterAeroForcesSPM.cpp.

101{
102}

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 226 of file FilterAeroForcesSPM.cpp.

230{
231 // Only output every m_outputFrequency.
232 if ((m_index % m_outputFrequency) || (time < m_startTime))
233 {
234 return;
235 }
236
237 int nq = pIntVel[0].size();
238 Array<OneD, NekDouble> tmp(nq);
239
240 // Aerodynamic forces are computed according eq. 18a in Luo et al. (2009).
241 // Smoothed profile method for particulate flows: Error analysis and
242 // simulations. Journal of Computational Physics, 228(5)
243 for (int i = 0; i < m_spaceDim; ++i)
244 {
245 // "Scalar" field to be integrated
246 Vmath::Vsub(nq, pIntVel[i], 1, pUpPrev[i], 1, tmp, 1);
247 Vmath::Vmul(nq, pPhi->GetPhys(), 1, tmp, 1, tmp, 1);
248
249 // Integration of force throughout the domain
250 m_Forces[i] = pPhi->Integral(tmp) / dt;
251 }
252}
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::weak_ptr< SolverUtils::EquationSystem > &  pEquation,
const std::map< std::string, std::string > &  pParams 
)
inlinestatic

Creates an instance of this class.

Definition at line 47 of file FilterAeroForcesSPM.h.

51 {
54 pSession, pEquation, pParams);
55 return p;
56 }
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_Finalise()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 195 of file FilterAeroForcesSPM.cpp.

198{
199 if (pFields[0]->GetComm()->GetRank() == 0)
200 {
201 m_outputStream.close();
202 }
203}

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 107 of file FilterAeroForcesSPM.cpp.

110{
111 // Save space dimension
112 m_spaceDim = pFields[0]->GetGraph()->GetMeshDimension();
113
114 // Fill in the directions vector
115 m_dirNames.push_back("X");
116 if (m_spaceDim > 1)
117 {
118 m_dirNames.push_back("Y");
119 }
120 if (m_spaceDim > 2)
121 {
122 m_dirNames.push_back("Z");
123 }
124
125 // Allocate the aero-forces vector
126 m_Forces = Array<OneD, NekDouble>(m_spaceDim);
127
128 // Write header
129 LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
130 if (vComm->GetRank() == 0)
131 {
132 // Open output stream
133 bool adaptive;
134 m_session->MatchSolverInfo("Driver", "Adaptive", adaptive, false);
135 if (adaptive)
136 {
137 m_outputStream.open(m_outputFile.c_str(), ofstream::app);
138 }
139 else
140 {
141 m_outputStream.open(m_outputFile.c_str());
142 }
143 m_outputStream << "# Forces acting on bodies\n";
144 m_outputStream << "#";
145 m_outputStream.width(7);
146 m_outputStream << "Time";
147 for (int i = 0; i < m_spaceDim; ++i)
148 {
149 m_outputStream.width(14);
150 m_outputStream << "F_" << m_dirNames[i];
151 }
152
153 m_outputStream << endl;
154 }
155
156 m_index = 0;
157}
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 208 of file FilterAeroForcesSPM.cpp.

209{
210 return true;
211}

◆ 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 162 of file FilterAeroForcesSPM.cpp.

165{
166 // Only output every m_outputFrequency.
167 if ((m_index++) % m_outputFrequency || (time < m_startTime))
168 {
169 return;
170 }
171
172 // Communicators to exchange results
173 LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
174
175 // Write Results
176 if (vComm->GetRank() == 0)
177 {
178 // Write time
179 m_outputStream.width(8);
180 m_outputStream << setprecision(6) << time;
181 // Write forces
182 for (int i = 0; i < m_spaceDim; ++i)
183 {
184 m_outputStream.width(15);
185 m_outputStream << setprecision(8) << m_Forces[i];
186 }
187 m_outputStream.width(10);
188 m_outputStream << endl;
189 }
190}

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

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::weak_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.
Definition: NekFactory.hpp:197
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:39

Name of the class.

Definition at line 59 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 82 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 86 of file FilterAeroForcesSPM.h.

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

◆ m_index

unsigned int Nektar::FilterAeroForcesSPM::m_index
protected

Definition at line 75 of file FilterAeroForcesSPM.h.

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

◆ m_outputFile

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

Definition at line 77 of file FilterAeroForcesSPM.h.

Referenced by FilterAeroForcesSPM(), and v_Initialise().

◆ m_outputFrequency

unsigned int Nektar::FilterAeroForcesSPM::m_outputFrequency
protected

Definition at line 76 of file FilterAeroForcesSPM.h.

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

◆ m_outputStream

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

Definition at line 78 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 84 of file FilterAeroForcesSPM.h.

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

◆ m_startTime

NekDouble Nektar::FilterAeroForcesSPM::m_startTime
protected

Definition at line 80 of file FilterAeroForcesSPM.h.

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