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 
37 #include <boost/core/ignore_unused.hpp>
38 
41 
42 using namespace std;
43 
44 namespace Nektar
45 {
46 namespace SolverUtils
47 {
48 std::string FilterMean::className =
50  FilterMean::create);
51 
52 FilterMean::FilterMean(const LibUtilities::SessionReaderSharedPtr &pSession,
53  const std::weak_ptr<EquationSystem> &pEquation,
54  const ParamMap &pParams)
55  : Filter(pSession, pEquation), m_index(-1), m_homogeneous(false), m_planes()
56 {
57  // OutputFile
58  auto it = pParams.find("OutputFile");
59  if (it == pParams.end())
60  {
61  m_outputFile = m_session->GetSessionName();
62  }
63  else
64  {
65  ASSERTL0(it->second.length() > 0, "Missing parameter 'OutputFile'.");
66  m_outputFile = it->second;
67  }
68  m_outputFile += ".avg";
69 
70  // OutputFrequency
71  it = pParams.find("OutputFrequency");
72  ASSERTL0(it != pParams.end(), "Missing parameter 'OutputFrequency'.");
73  LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
74  m_outputFrequency = round(equ.Evaluate());
75 
76  pSession->LoadParameter("LZ", m_homogeneousLength, 0.0);
77 }
78 
80 {
81 }
82 
85  const NekDouble &time)
86 {
88  int spacedim = 2;
89 
90  ASSERTL0(pFields[0]->GetExpType() != MultiRegions::e1D,
91  "1D expansion not supported for mean filter");
92 
93  ASSERTL0(pFields[0]->GetExpType() != MultiRegions::e3DH2D,
94  "Homogeneous 2D expansion not supported for mean filter");
95 
96  // Lock equation system pointer
97  auto equ = m_equ.lock();
98  ASSERTL0(equ, "Weak pointer expired");
99 
100  auto fluidEqu = std::dynamic_pointer_cast<FluidInterface>(equ);
101  ASSERTL0(fluidEqu, "Mean filter is incompatible with this solver.");
102 
103  if (pFields[0]->GetExpType() == MultiRegions::e3DH1D)
104  {
105  m_homogeneous = true;
106  spacedim = 3;
107  }
108 
109  // Calculate area/volume of domain.
110  if (m_homogeneous)
111  {
112  m_planes = pFields[0]->GetZIDs();
113  areaField = pFields[0]->GetPlane(0);
114  }
115  else
116  {
117  areaField = pFields[0];
118  }
119 
120  Array<OneD, NekDouble> inarray(areaField->GetNpoints(), 1.0);
121  m_area = areaField->Integral(inarray);
122 
123  if (m_homogeneous)
124  {
126  }
127 
128  // Open OutputFile
129  std::string volname[3] = {"length", "area", "volume"};
130  LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
131  if (vComm->GetRank() == 0)
132  {
133  m_outputStream.open(m_outputFile.c_str());
134  ASSERTL0(m_outputStream.good(),
135  "Unable to open: '" + m_outputFile + "'");
136  m_outputStream.setf(ios::scientific, ios::floatfield);
137  m_outputStream << "# Time";
138  for (int i = 0; i < pFields.size(); ++i)
139  m_outputStream << setw(22) << equ->GetVariable(i);
140  m_outputStream << setw(22) << volname[spacedim - 1] << " " << m_area;
141  m_outputStream << endl;
142  }
143 
144  // Output values at initial time.
145  m_index = 0;
146  v_Update(pFields, time);
147 }
148 
151  const NekDouble &time)
152 {
153  // Only output every m_outputFrequency
154  if ((m_index++) % m_outputFrequency)
155  {
156  return;
157  }
158 
159  int i;
160  LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
161 
162  // Calculate average values.
163  Array<OneD, NekDouble> avg(pFields.size());
164  for (i = 0; i < pFields.size(); ++i)
165  {
166  avg[i] = 0.0;
167  }
168 
169  if (m_homogeneous)
170  {
171  for (i = 0; i < pFields.size(); ++i)
172  avg[i] = pFields[0]->GetPlane(0)->Integral(pFields[i]->GetPhys()) *
174  }
175  else
176  {
177  for (i = 0; i < pFields.size(); ++i)
178  avg[i] = pFields[0]->Integral(pFields[i]->GetPhys());
179  }
180 
181  for (i = 0; i < pFields.size(); ++i)
182  avg[i] /= m_area;
183 
184  if (vComm->GetRank() == 0)
185  {
186  m_outputStream << setw(17) << setprecision(8) << time;
187  for (int i = 0; i < pFields.size(); ++i)
188  m_outputStream << setw(22) << setprecision(11) << avg[i];
189  m_outputStream << endl;
190  }
191 }
192 
195  const NekDouble &time)
196 {
197  boost::ignore_unused(pFields, time);
198 
199  if (pFields[0]->GetComm()->GetRank() == 0)
200  {
201  m_outputStream.close();
202  }
203 }
204 
206 {
207  return true;
208 }
209 
210 } // namespace SolverUtils
211 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:85
const std::weak_ptr< EquationSystem > m_equ
Definition: Filter.h:86
std::map< std::string, std::string > ParamMap
Definition: Filter.h:67
virtual SOLVER_UTILS_EXPORT void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pField, const NekDouble &time) override
Definition: FilterMean.cpp:149
virtual SOLVER_UTILS_EXPORT void v_Finalise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pField, const NekDouble &time) override
Definition: FilterMean.cpp:193
virtual SOLVER_UTILS_EXPORT bool v_IsTimeDependent() override
Definition: FilterMean.cpp:205
virtual SOLVER_UTILS_EXPORT void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pField, const NekDouble &time) override
Definition: FilterMean.cpp:83
Array< OneD, unsigned int > m_planes
Definition: FilterMean.h:88
SOLVER_UTILS_EXPORT ~FilterMean()
Definition: FilterMean.cpp:79
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:54
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:41
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble