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

#include <FilterEnergy.h>

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

Public Member Functions

SOLVER_UTILS_EXPORT FilterEnergy (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
 
SOLVER_UTILS_EXPORT ~FilterEnergy () 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 SolverUtils::FilterSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const ParamMap &pParams)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of the class. More...
 

Protected Member Functions

SOLVER_UTILS_EXPORT void v_Initialise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pField, const NekDouble &time) override
 
SOLVER_UTILS_EXPORT void v_Update (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pField, const NekDouble &time) override
 
SOLVER_UTILS_EXPORT void v_Finalise (const Array< OneD, const MultiRegions::ExpListSharedPtr > &pField, const NekDouble &time) override
 
SOLVER_UTILS_EXPORT 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
 

Private Attributes

unsigned int m_index
 
unsigned int m_outputFrequency
 
std::ofstream m_outFile
 
bool m_homogeneous
 
NekDouble m_homogeneousLength
 
NekDouble m_area
 
LibUtilities::CommSharedPtr m_comm
 
Array< OneD, unsigned int > m_planes
 

Additional Inherited Members

- Public Types inherited from Nektar::SolverUtils::Filter
typedef std::map< std::string, std::string > ParamMap
 
- Protected Attributes inherited from Nektar::SolverUtils::Filter
LibUtilities::SessionReaderSharedPtr m_session
 
const std::weak_ptr< EquationSystemm_equ
 

Detailed Description

Definition at line 42 of file FilterEnergy.h.

Constructor & Destructor Documentation

◆ FilterEnergy()

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

Definition at line 48 of file FilterEnergy.cpp.

51 : Filter(pSession, pEquation), m_index(-1), m_homogeneous(false), m_planes()
52{
53 std::string outName;
54
55 // OutputFile
56 auto it = pParams.find("OutputFile");
57 if (it == pParams.end())
58 {
59 outName = m_session->GetSessionName();
60 }
61 else
62 {
63 ASSERTL0(it->second.length() > 0, "Missing parameter 'OutputFile'.");
64 outName = it->second;
65 }
66 outName += ".eny";
67
68 m_comm = pSession->GetComm();
69 if (m_comm->GetRank() == 0)
70 {
71 m_outFile.open(outName.c_str());
72 ASSERTL0(m_outFile.good(), "Unable to open: '" + outName + "'");
73 m_outFile.setf(ios::scientific, ios::floatfield);
74 m_outFile << "# Time Kinetic energy "
75 << "Enstrophy" << endl
76 << "# ---------------------------------------------"
77 << "--------------" << endl;
78 }
79 pSession->LoadParameter("LZ", m_homogeneousLength, 0.0);
80
81 // OutputFrequency
82 it = pParams.find("OutputFrequency");
83 ASSERTL0(it != pParams.end(), "Missing parameter 'OutputFrequency'.");
84 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
85 m_outputFrequency = round(equ.Evaluate());
86}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
LibUtilities::CommSharedPtr m_comm
Definition: FilterEnergy.h:85
Array< OneD, unsigned int > m_planes
Definition: FilterEnergy.h:86
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_comm, m_homogeneousLength, m_outFile, m_outputFrequency, and Nektar::SolverUtils::Filter::m_session.

◆ ~FilterEnergy()

Nektar::SolverUtils::FilterEnergy::~FilterEnergy ( )
override

Definition at line 88 of file FilterEnergy.cpp.

89{
90}

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 46 of file FilterEnergy.h.

50 {
53 pParams);
54 return p;
55 }
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::SolverUtils::FilterEnergy::v_Finalise ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pField,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 261 of file FilterEnergy.cpp.

265{
266 m_outFile.close();
267}

References m_outFile.

◆ v_Initialise()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 92 of file FilterEnergy.cpp.

95{
96 m_index = -1;
98
99 ASSERTL0(pFields[0]->GetExpType() != MultiRegions::e1D,
100 "1D expansion not supported for energy filter");
101
102 ASSERTL0(pFields[0]->GetExpType() != MultiRegions::e2D,
103 "2D expansion not supported for energy filter");
104
105 ASSERTL0(pFields[0]->GetExpType() != MultiRegions::e3DH2D,
106 "Homogeneous 2D expansion not supported for energy filter");
107
108 if (pFields[0]->GetExpType() == MultiRegions::e3DH1D)
109 {
110 m_homogeneous = true;
111 }
112
113 // Calculate area/volume of domain.
114 if (m_homogeneous)
115 {
116 m_planes = pFields[0]->GetZIDs();
117 areaField = pFields[0]->GetPlane(0);
118 }
119 else
120 {
121 areaField = pFields[0];
122 }
123
124 Array<OneD, NekDouble> inarray(areaField->GetNpoints(), 1.0);
125 m_area = areaField->Integral(inarray);
126
127 if (m_homogeneous)
128 {
130 }
131
132 // Output values at initial time.
133 v_Update(pFields, time);
134}
SOLVER_UTILS_EXPORT void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pField, const NekDouble &time) override
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.

References ASSERTL0, Nektar::MultiRegions::e1D, Nektar::MultiRegions::e2D, Nektar::MultiRegions::e3DH1D, Nektar::MultiRegions::e3DH2D, m_area, m_homogeneous, m_homogeneousLength, m_index, m_planes, and v_Update().

◆ v_IsTimeDependent()

bool Nektar::SolverUtils::FilterEnergy::v_IsTimeDependent ( )
overrideprotectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 269 of file FilterEnergy.cpp.

270{
271 return true;
272}

◆ v_Update()

void Nektar::SolverUtils::FilterEnergy::v_Update ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pField,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 136 of file FilterEnergy.cpp.

139{
140 int i, nPoints = pFields[0]->GetNpoints();
141
142 m_index++;
143
144 if (m_index % m_outputFrequency > 0)
145 {
146 return;
147 }
148
149 // Lock equation system pointer
150 auto equ = m_equ.lock();
151 ASSERTL0(equ, "Weak pointer expired");
152
153 auto fluidEqu = std::dynamic_pointer_cast<FluidInterface>(equ);
154 ASSERTL0(fluidEqu, "Energy filter is incompatible with this solver.");
155
156 // Store physical values in an array
157 Array<OneD, Array<OneD, NekDouble>> physfields(pFields.size());
158 for (i = 0; i < pFields.size(); ++i)
159 {
160 physfields[i] = pFields[i]->GetPhys();
161 }
162
163 // Calculate kinetic energy.
164 NekDouble Ek = 0.0;
165 Array<OneD, NekDouble> tmp(nPoints, 0.0);
166 Array<OneD, NekDouble> density;
167 Array<OneD, Array<OneD, NekDouble>> u(3);
168 for (i = 0; i < 3; ++i)
169 {
170 u[i] = Array<OneD, NekDouble>(nPoints);
171 }
172 fluidEqu->GetVelocity(physfields, u);
173
174 for (i = 0; i < 3; ++i)
175 {
176 if (m_homogeneous && pFields[i]->GetWaveSpace())
177 {
178 pFields[i]->HomogeneousBwdTrans(nPoints, u[i], u[i]);
179 }
180
181 Vmath::Vvtvp(nPoints, u[i], 1, u[i], 1, tmp, 1, tmp, 1);
182 }
183
184 if (!fluidEqu->HasConstantDensity())
185 {
186 density = Array<OneD, NekDouble>(nPoints);
187 fluidEqu->GetDensity(physfields, density);
188 Vmath::Vmul(nPoints, density, 1, tmp, 1, tmp, 1);
189 }
190
191 if (m_homogeneous)
192 {
193 Array<OneD, NekDouble> tmp2(nPoints, 0.0);
194 pFields[0]->HomogeneousFwdTrans(nPoints, tmp, tmp2);
195 Ek = pFields[0]->GetPlane(0)->Integral(tmp2) * m_homogeneousLength;
196 }
197 else
198 {
199 Ek = pFields[0]->Integral(tmp);
200 }
201
202 Ek /= 2.0 * m_area;
203
204 if (m_comm->GetRank() == 0)
205 {
206 m_outFile << setw(17) << setprecision(8) << time << setw(22)
207 << setprecision(11) << Ek;
208 }
209
210 bool waveSpace[3] = {pFields[0]->GetWaveSpace(), pFields[1]->GetWaveSpace(),
211 pFields[2]->GetWaveSpace()};
212
213 if (m_homogeneous)
214 {
215 for (i = 0; i < 3; ++i)
216 {
217 pFields[i]->SetWaveSpace(false);
218 }
219 }
220
221 // First calculate vorticity field.
222 Array<OneD, NekDouble> tmp2(nPoints), tmp3(nPoints);
223 Vmath::Zero(nPoints, tmp, 1);
224 for (i = 0; i < 3; ++i)
225 {
226 int f1 = (i + 2) % 3, c2 = f1;
227 int c1 = (i + 1) % 3, f2 = c1;
228 pFields[f1]->PhysDeriv(c1, u[f1], tmp2);
229 pFields[f2]->PhysDeriv(c2, u[f2], tmp3);
230 Vmath::Vsub(nPoints, tmp2, 1, tmp3, 1, tmp2, 1);
231 Vmath::Vvtvp(nPoints, tmp2, 1, tmp2, 1, tmp, 1, tmp, 1);
232 }
233
234 if (!fluidEqu->HasConstantDensity())
235 {
236 Vmath::Vmul(nPoints, density, 1, tmp, 1, tmp, 1);
237 }
238
239 if (m_homogeneous)
240 {
241 for (i = 0; i < 3; ++i)
242 {
243 pFields[i]->SetWaveSpace(waveSpace[i]);
244 }
245 pFields[0]->HomogeneousFwdTrans(nPoints, tmp, tmp);
246 Ek = pFields[0]->GetPlane(0)->Integral(tmp) * m_homogeneousLength;
247 }
248 else
249 {
250 Ek = pFields[0]->Integral(tmp);
251 }
252
253 Ek /= 2.0 * m_area;
254
255 if (m_comm->GetRank() == 0)
256 {
257 m_outFile << setw(22) << setprecision(11) << Ek << endl;
258 }
259}
const std::weak_ptr< EquationSystem > m_equ
Definition: Filter.h:84
double NekDouble
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 Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
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 ASSERTL0, m_area, m_comm, Nektar::SolverUtils::Filter::m_equ, m_homogeneous, m_homogeneousLength, m_index, m_outFile, m_outputFrequency, Vmath::Vmul(), Vmath::Vsub(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by v_Initialise().

Member Data Documentation

◆ className

std::string Nektar::SolverUtils::FilterEnergy::className
static
Initial value:
=
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static SolverUtils::FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const ParamMap &pParams)
Creates an instance of this class.
Definition: FilterEnergy.h:46
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:39

Name of the class.

Definition at line 58 of file FilterEnergy.h.

◆ m_area

NekDouble Nektar::SolverUtils::FilterEnergy::m_area
private

Definition at line 84 of file FilterEnergy.h.

Referenced by v_Initialise(), and v_Update().

◆ m_comm

LibUtilities::CommSharedPtr Nektar::SolverUtils::FilterEnergy::m_comm
private

Definition at line 85 of file FilterEnergy.h.

Referenced by FilterEnergy(), and v_Update().

◆ m_homogeneous

bool Nektar::SolverUtils::FilterEnergy::m_homogeneous
private

Definition at line 82 of file FilterEnergy.h.

Referenced by v_Initialise(), and v_Update().

◆ m_homogeneousLength

NekDouble Nektar::SolverUtils::FilterEnergy::m_homogeneousLength
private

Definition at line 83 of file FilterEnergy.h.

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

◆ m_index

unsigned int Nektar::SolverUtils::FilterEnergy::m_index
private

Definition at line 79 of file FilterEnergy.h.

Referenced by v_Initialise(), and v_Update().

◆ m_outFile

std::ofstream Nektar::SolverUtils::FilterEnergy::m_outFile
private

Definition at line 81 of file FilterEnergy.h.

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

◆ m_outputFrequency

unsigned int Nektar::SolverUtils::FilterEnergy::m_outputFrequency
private

Definition at line 80 of file FilterEnergy.h.

Referenced by FilterEnergy(), and v_Update().

◆ m_planes

Array<OneD, unsigned int> Nektar::SolverUtils::FilterEnergy::m_planes
private

Definition at line 86 of file FilterEnergy.h.

Referenced by v_Initialise().