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

#include <FilterMean.h>

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

Public Member Functions

SOLVER_UTILS_EXPORT FilterMean (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
 
SOLVER_UTILS_EXPORT ~FilterMean () 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::string m_outputFile
 
std::ofstream m_outputStream
 
bool m_homogeneous
 
NekDouble m_homogeneousLength
 
NekDouble m_area
 
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 FilterMean.h.

Constructor & Destructor Documentation

◆ FilterMean()

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

Definition at line 48 of file FilterMean.cpp.

51 : Filter(pSession, pEquation), m_index(-1), m_homogeneous(false), m_planes()
52{
53 // OutputFile
54 auto it = pParams.find("OutputFile");
55 if (it == pParams.end())
56 {
57 m_outputFile = m_session->GetSessionName();
58 }
59 else
60 {
61 ASSERTL0(it->second.length() > 0, "Missing parameter 'OutputFile'.");
62 m_outputFile = it->second;
63 }
64 m_outputFile += ".avg";
65
66 // OutputFrequency
67 it = pParams.find("OutputFrequency");
68 ASSERTL0(it != pParams.end(), "Missing parameter 'OutputFrequency'.");
69 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
70 m_outputFrequency = round(equ.Evaluate());
71
72 pSession->LoadParameter("LZ", m_homogeneousLength, 0.0);
73}
#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
Array< OneD, unsigned int > m_planes
Definition: FilterMean.h:86

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

◆ ~FilterMean()

Nektar::SolverUtils::FilterMean::~FilterMean ( )
override

Definition at line 75 of file FilterMean.cpp.

76{
77}

Member Function Documentation

◆ create()

static SolverUtils::FilterSharedPtr Nektar::SolverUtils::FilterMean::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 FilterMean.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::FilterMean::v_Finalise ( const Array< OneD, const MultiRegions::ExpListSharedPtr > &  pField,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::Filter.

Definition at line 199 of file FilterMean.cpp.

203{
204 if (pFields[0]->GetComm()->GetRank() == 0)
205 {
206 m_outputStream.close();
207 }
208}

References m_outputStream.

◆ v_Initialise()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 79 of file FilterMean.cpp.

82{
84 int spacedim = 2;
85
86 ASSERTL0(pFields[0]->GetExpType() != MultiRegions::e1D,
87 "1D expansion not supported for mean filter");
88
89 ASSERTL0(pFields[0]->GetExpType() != MultiRegions::e3DH2D,
90 "Homogeneous 2D expansion not supported for mean filter");
91
92 // Lock equation system pointer
93 auto equ = m_equ.lock();
94 ASSERTL0(equ, "Weak pointer expired");
95
96 auto fluidEqu = std::dynamic_pointer_cast<FluidInterface>(equ);
97 ASSERTL0(fluidEqu, "Mean filter is incompatible with this solver.");
98
99 if (pFields[0]->GetExpType() == MultiRegions::e3DH1D)
100 {
101 m_homogeneous = true;
102 spacedim = 3;
103 }
104
105 // Calculate area/volume of domain.
106 if (m_homogeneous)
107 {
108 m_planes = pFields[0]->GetZIDs();
109 areaField = pFields[0]->GetPlane(0);
110 }
111 else
112 {
113 areaField = pFields[0];
114 }
115
116 Array<OneD, NekDouble> inarray(areaField->GetNpoints(), 1.0);
117 m_area = areaField->Integral(inarray);
118
119 if (m_homogeneous)
120 {
122 }
123
124 // Open OutputFile
125 std::string volname[3] = {"length", "area", "volume"};
126 LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
127 if (vComm->GetRank() == 0)
128 {
129 m_outputStream.open(m_outputFile.c_str());
131 "Unable to open: '" + m_outputFile + "'");
132 m_outputStream.setf(ios::scientific, ios::floatfield);
133 m_outputStream << "# Time";
134 for (int i = 0; i < pFields.size(); ++i)
135 {
136 m_outputStream << setw(22) << equ->GetVariable(i);
137 }
138 m_outputStream << setw(22) << volname[spacedim - 1] << " " << m_area;
139 m_outputStream << endl;
140 }
141
142 // Output values at initial time.
143 m_index = 0;
144 v_Update(pFields, time);
145}
const std::weak_ptr< EquationSystem > m_equ
Definition: Filter.h:84
SOLVER_UTILS_EXPORT void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pField, const NekDouble &time) override
Definition: FilterMean.cpp:147
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.

References ASSERTL0, Nektar::MultiRegions::e1D, Nektar::MultiRegions::e3DH1D, Nektar::MultiRegions::e3DH2D, m_area, Nektar::SolverUtils::Filter::m_equ, m_homogeneous, m_homogeneousLength, m_index, m_outputFile, m_outputStream, m_planes, and v_Update().

◆ v_IsTimeDependent()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 210 of file FilterMean.cpp.

211{
212 return true;
213}

◆ v_Update()

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

Implements Nektar::SolverUtils::Filter.

Definition at line 147 of file FilterMean.cpp.

150{
151 // Only output every m_outputFrequency
152 if ((m_index++) % m_outputFrequency)
153 {
154 return;
155 }
156
157 int i;
158 LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
159
160 // Calculate average values.
161 Array<OneD, NekDouble> avg(pFields.size());
162 for (i = 0; i < pFields.size(); ++i)
163 {
164 avg[i] = 0.0;
165 }
166
167 if (m_homogeneous)
168 {
169 for (i = 0; i < pFields.size(); ++i)
170 {
171 avg[i] = pFields[0]->GetPlane(0)->Integral(pFields[i]->GetPhys()) *
173 }
174 }
175 else
176 {
177 for (i = 0; i < pFields.size(); ++i)
178 {
179 avg[i] = pFields[0]->Integral(pFields[i]->GetPhys());
180 }
181 }
182
183 for (i = 0; i < pFields.size(); ++i)
184 {
185 avg[i] /= m_area;
186 }
187
188 if (vComm->GetRank() == 0)
189 {
190 m_outputStream << setw(17) << setprecision(8) << time;
191 for (int i = 0; i < pFields.size(); ++i)
192 {
193 m_outputStream << setw(22) << setprecision(11) << avg[i];
194 }
195 m_outputStream << endl;
196 }
197}

References m_area, m_homogeneous, m_homogeneousLength, m_index, m_outputFrequency, and m_outputStream.

Referenced by v_Initialise().

Member Data Documentation

◆ className

std::string Nektar::SolverUtils::FilterMean::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: FilterMean.h:46
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:39

Name of the class.

Definition at line 58 of file FilterMean.h.

◆ m_area

NekDouble Nektar::SolverUtils::FilterMean::m_area
private

Definition at line 85 of file FilterMean.h.

Referenced by v_Initialise(), and v_Update().

◆ m_homogeneous

bool Nektar::SolverUtils::FilterMean::m_homogeneous
private

Definition at line 83 of file FilterMean.h.

Referenced by v_Initialise(), and v_Update().

◆ m_homogeneousLength

NekDouble Nektar::SolverUtils::FilterMean::m_homogeneousLength
private

Definition at line 84 of file FilterMean.h.

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

◆ m_index

unsigned int Nektar::SolverUtils::FilterMean::m_index
private

Definition at line 79 of file FilterMean.h.

Referenced by v_Initialise(), and v_Update().

◆ m_outputFile

std::string Nektar::SolverUtils::FilterMean::m_outputFile
private

Definition at line 81 of file FilterMean.h.

Referenced by FilterMean(), and v_Initialise().

◆ m_outputFrequency

unsigned int Nektar::SolverUtils::FilterMean::m_outputFrequency
private

Definition at line 80 of file FilterMean.h.

Referenced by FilterMean(), and v_Update().

◆ m_outputStream

std::ofstream Nektar::SolverUtils::FilterMean::m_outputStream
private

Definition at line 82 of file FilterMean.h.

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

◆ m_planes

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

Definition at line 86 of file FilterMean.h.

Referenced by v_Initialise().