Nektar++
FilterEnergy1D.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: FilterEnergy1D.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: Outputs orthogonal expansion of 1D elements.
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <boost/core/ignore_unused.hpp>
36
39
40using namespace std;
41
42namespace Nektar
43{
44namespace SolverUtils
45{
46std::string FilterEnergy1D::className =
49
50/**
51 * @brief Set up filter with output file and frequency parameters.
52 *
53 * @param pSession Current session.
54 * @param pParams Map of parameters defined in XML file.
55 */
58 const std::weak_ptr<EquationSystem> &pEquation, const ParamMap &pParams)
59 : Filter(pSession, pEquation), m_index(0)
60{
61 ASSERTL0(pSession->GetComm()->GetSize() == 1,
62 "The 1D energy filter currently only works in serial.");
63
64 std::string outName;
65
66 // OutputFile
67 auto it = pParams.find("OutputFile");
68 if (it == pParams.end())
69 {
70 outName = m_session->GetSessionName();
71 }
72 else
73 {
74 ASSERTL0(it->second.length() > 0, "Missing parameter 'OutputFile'.");
75 outName = it->second;
76 }
77 outName += ".eny";
78 m_out.open(outName.c_str());
79
80 // OutputFrequency
81 it = pParams.find("OutputFrequency");
82 if (it == pParams.end())
83 {
85 }
86 else
87 {
88 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
89 m_outputFrequency = round(equ.Evaluate());
90 }
91}
92
93/**
94 * @brief Destructor.
95 */
97{
98}
99
100/**
101 * @brief Initialize filter.
102 */
105 const NekDouble &time)
106{
107 boost::ignore_unused(time);
108 ASSERTL0(pFields[0]->GetExp(0)->GetNumBases() == 1,
109 "The Energy 1D filter is only valid in 1D.");
110}
111
112/**
113 * @brief Update filter output with the current timestep's orthogonal
114 * coefficients.
115 */
118 const NekDouble &time)
119{
120 // Only output every m_outputFrequency
121 if ((m_index++) % m_outputFrequency)
122 {
123 return;
124 }
125
126 int nElmt = pFields[0]->GetExpSize();
127
128 // Loop over all elements
129 m_out << "##" << endl;
130 m_out << "## Time = " << time << endl;
131 m_out << "##" << endl;
132
133 for (int i = 0; i < nElmt; ++i)
134 {
135 // Figure out number of modes in this expansion.
136 LocalRegions::ExpansionSharedPtr exp = pFields[0]->GetExp(i);
137 int nModes = exp->GetBasis(0)->GetNumModes();
138
139 // Set uo basis key for orthogonal basis
141 LibUtilities::BasisKey bkeyOrth(btype, nModes,
142 exp->GetBasis(0)->GetPointsKey());
143
144 // Get basis key for existing expansion
145 LibUtilities::BasisKey bkey(exp->GetBasis(0)->GetBasisType(),
146 exp->GetBasis(0)->GetNumModes(),
147 exp->GetBasis(0)->GetPointsKey());
148
149 // Find coeffs for this element in the list of all coefficients
151 pFields[0]->GetCoeffs() + pFields[0]->GetCoeff_Offset(i);
152
153 // Storage for orthogonal coefficients
154 Array<OneD, NekDouble> coeffsOrth(nModes);
155
156 // Project from coeffs -> orthogonal coeffs
157 LibUtilities::InterpCoeff1D(bkey, coeffs, bkeyOrth, coeffsOrth);
158
159 // Write coeffs to file
160 m_out << "# Element " << i << " (ID " << exp->GetGeom()->GetGlobalID()
161 << ")" << endl;
162 for (int j = 0; j < nModes; ++j)
163 {
164 m_out << coeffsOrth[j] << endl;
165 }
166 }
167 m_out << endl;
168}
169
172 const NekDouble &time)
173{
174 boost::ignore_unused(pFields, time);
175 m_out.close();
176}
177
179{
180 return true;
181}
182} // namespace SolverUtils
183} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
Describes the specification for a Basis.
Definition: Basis.h:47
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
virtual void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pField, const NekDouble &time) override
Initialize filter.
static std::string className
Name of the class.
virtual bool v_IsTimeDependent() override
SOLVER_UTILS_EXPORT FilterEnergy1D(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
Set up filter with output file and frequency parameters.
unsigned int m_outputFrequency
Output frequency.
virtual void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pField, const NekDouble &time) override
Update filter output with the current timestep's orthogonal coefficients.
unsigned int m_index
Current index counter.
SOLVER_UTILS_EXPORT ~FilterEnergy1D()
Destructor.
virtual void v_Finalise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pField, const NekDouble &time) override
std::ofstream m_out
Output file.
static FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
Creates an instance of this class.
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:85
std::map< std::string, std::string > ParamMap
Definition: Filter.h:67
std::shared_ptr< SessionReader > SessionReaderSharedPtr
void InterpCoeff1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
Definition: InterpCoeff.cpp:44
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:44
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:41
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble