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