Nektar++
FilterMean.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: FilterMean.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Output mean.
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <iomanip>
36
39
40using namespace std;
41
42namespace Nektar::SolverUtils
43{
44std::string FilterMean::className =
47
49 const std::weak_ptr<EquationSystem> &pEquation,
50 const ParamMap &pParams)
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}
74
76{
77}
78
81 const NekDouble &time)
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}
146
149 const NekDouble &time)
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}
198
201 &pFields,
202 [[maybe_unused]] const NekDouble &time)
203{
204 if (pFields[0]->GetComm()->GetRank() == 0)
205 {
206 m_outputStream.close();
207 }
208}
209
211{
212 return true;
213}
214
215} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:83
const std::weak_ptr< EquationSystem > m_equ
Definition: Filter.h:84
std::map< std::string, std::string > ParamMap
Definition: Filter.h:65
SOLVER_UTILS_EXPORT void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pField, const NekDouble &time) override
Definition: FilterMean.cpp:147
SOLVER_UTILS_EXPORT void v_Finalise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pField, const NekDouble &time) override
Definition: FilterMean.cpp:199
SOLVER_UTILS_EXPORT bool v_IsTimeDependent() override
Definition: FilterMean.cpp:210
SOLVER_UTILS_EXPORT void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pField, const NekDouble &time) override
Definition: FilterMean.cpp:79
Array< OneD, unsigned int > m_planes
Definition: FilterMean.h:86
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
SOLVER_UTILS_EXPORT ~FilterMean() override
Definition: FilterMean.cpp:75
SOLVER_UTILS_EXPORT FilterMean(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
Definition: FilterMean.cpp:48
static std::string className
Name of the class.
Definition: FilterMean.h:58
std::shared_ptr< SessionReader > SessionReaderSharedPtr
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.
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:39
double NekDouble