Nektar++
FilterOffsetPhase.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File FilterOffsetPhase.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// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: Outputs phase fields during time-stepping.
33//
34///////////////////////////////////////////////////////////////////////////////
35
37
38namespace Nektar
39{
42 "OffsetPhase", FilterOffsetPhase::create);
43
46 const std::weak_ptr<SolverUtils::EquationSystem> &pEquation,
47 const ParamMap &pParams)
48 : Filter(pSession, pEquation)
49{
50 if (pParams.find("OutputFile") == pParams.end())
51 {
52 m_outputFile = m_session->GetSessionName();
53 }
54 else
55 {
56 ASSERTL0(!(pParams.find("OutputFile")->second.empty()),
57 "Missing parameter 'OutputFile'.");
58 m_outputFile = pParams.find("OutputFile")->second;
59 }
60 ASSERTL0(pParams.find("OutputFrequency") != pParams.end(),
61 "Missing parameter 'OutputFrequency'.");
62 LibUtilities::Equation equ(m_session->GetInterpreter(),
63 pParams.find("OutputFrequency")->second);
64 m_outputFrequency = floor(equ.Evaluate());
65
66 if (pParams.find("MeanV") != pParams.end())
67 {
68 LibUtilities::Equation equ(m_session->GetInterpreter(),
69 pParams.find("MeanV")->second);
70 meanV = equ.Evaluate();
71 }
72
73 m_outputIndex = 0;
74 m_index = 0;
76}
77
79{
80}
81
84 [[maybe_unused]] const NekDouble &time)
85{
86 m_index = 0;
87 m_outputIndex = 0;
88
89 oldV = pFields[0]->GetPhys(); // Allocate storage for past data
90
91 v_Update(pFields, 0.0);
92}
93
96 const NekDouble &time)
97{
98 if (m_index++ % m_outputFrequency > 0)
99 {
100 return;
101 }
102
103 std::stringstream vOutputFilename;
104 vOutputFilename << m_outputFile << "_" << m_outputIndex << ".chk";
105
106 // Get FieldDefinitions
107 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef =
108 pFields[0]->GetFieldDefinitions();
109 std::vector<std::vector<NekDouble>> FieldData(FieldDef.size());
110
111 Array<OneD, NekDouble> currentV = pFields[0]->GetPhys();
112 Array<OneD, NekDouble> phase(currentV.size());
113 Array<OneD, NekDouble> phase_coeff(pFields[0]->GetNcoeffs());
114
115 // Calculate phase angle in physical space:
116 for (int i = 0; i < currentV.size(); ++i)
117 {
118 phase[i] = atan2(currentV[i] - meanV,
119 oldV[i] - meanV); // subtract mean -57.2986
120 }
121
122 // Convert to coefficient space:
123 pFields[0]->FwdTransLocalElmt(phase, phase_coeff);
124
125 // for each item in FieldDef
126 // push "phase" to FieldDef[i]->m_fields
127 // call AppendFieldData(FieldDef[i], FieldData[i], data)
128 for (int i = 0; i < FieldDef.size(); ++i)
129 {
130 FieldDef[i]->m_fields.push_back("phase");
131 pFields[0]->AppendFieldData(FieldDef[i], FieldData[i], phase_coeff);
132 }
133
134 // Save currentV
135 oldV = currentV;
136
137 // Update time in field info if required
138 LibUtilities::FieldMetaDataMap fieldMetaDataMap;
139 fieldMetaDataMap["Time"] = boost::lexical_cast<std::string>(time);
140
141 m_fld->Write(vOutputFilename.str(), FieldDef, FieldData, fieldMetaDataMap);
143}
144
147 &pFields,
148 [[maybe_unused]] const NekDouble &time)
149{
150}
151
153{
154 return true;
155}
156} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
static std::string className
Name of the class.
Array< OneD, NekDouble > oldV
static SolverUtils::FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const ParamMap &pParams)
Creates an instance of this class.
void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
FilterOffsetPhase(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const ParamMap &pParams)
LibUtilities::FieldIOSharedPtr m_fld
void v_Finalise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
static std::shared_ptr< FieldIO > CreateDefault(const LibUtilities::SessionReaderSharedPtr session)
Returns an object for the default FieldIO method.
Definition: FieldIO.cpp:195
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
std::map< std::string, std::string > ParamMap
Definition: Filter.h:65
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:50
std::shared_ptr< SessionReader > SessionReaderSharedPtr
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:39
double NekDouble