Nektar++
FilterLagrangianPoints.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: FilterLagrangianPoints.h
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: base class for physics evaluation.
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef NEKTAR_SOLVERUTILS_FILTERS_FILTERLAGRANGIANPOINTS_H
36#define NEKTAR_SOLVERUTILS_FILTERS_FILTERLAGRANGIANPOINTS_H
37
43
44namespace Nektar::SolverUtils
45{
46
48
50{
51public:
53
54 /// Creates an instance of this class
55 static StationaryPointsSharedPtr create(int rank, int dim, int intOrder,
56 int idOffset, NekDouble dt,
57 const std::vector<int> &Np,
58 const std::vector<NekDouble> &Box,
59 std::vector<std::string> extraVars)
60 {
63 rank, dim, intOrder, idOffset, dt, Np, Box, extraVars);
64 return p;
65 }
66
67 StatLagrangianPoints(int rank, int dim, int intOrder, int idOffset,
68 NekDouble dt, const std::vector<int> &Np,
69 const std::vector<NekDouble> &Box,
70 std::vector<std::string> extraVars);
71
73 {
74 }
75
76protected:
77 void v_OutputData(std::string filename, bool verbose,
78 std::map<std::string, NekDouble> &params) override;
79 void v_TimeAdvance(int order) override;
80 void v_GetCoords(int pid, Array<OneD, NekDouble> &gcoords) override;
81 void v_SetCoords(int pid, const Array<OneD, NekDouble> &gcoords) override;
82 void v_GetPhysics(int pid, Array<OneD, NekDouble> &data) override;
83 void v_SetPhysics(int pid, const Array<OneD, NekDouble> &data) override;
84 void v_ReSize(int Np) override;
85 void v_AssignPoint(int id, int pid,
86 const Array<OneD, NekDouble> &gcoords) override;
90 std::vector<std::string> m_extraPhysVars;
93 std::vector<int> m_N;
95};
96
98{
99public:
101
104 const std::shared_ptr<EquationSystem> &pEquation,
105 const std::map<std::string, std::string> &pParams)
106 {
109 pSession, pEquation, pParams);
110 return p;
111 }
112
113 static std::string className;
114
117 const std::shared_ptr<EquationSystem> &pEquation,
118 const std::map<std::string, std::string> &pParams);
120
121protected:
122 void ExtraPhysicsVars(std::vector<std::string> &extraVars);
125 const NekDouble &time) override;
126
127 void v_Update(
129 const NekDouble &time) override;
130
131 void v_Finalise(
133 const NekDouble &time) override;
134 bool v_IsTimeDependent() override;
136 Array<OneD, NekDouble> vel) override;
137 void GetPhysicsData(
139 std::vector<Array<OneD, NekDouble>> &PhysicsData);
140 void OutputSamplePoints(NekDouble time);
141 void OutputStatPoints(NekDouble time);
142 std::vector<std::string> m_defaultValues;
144 {
145 std::map<int, LibUtilities::EquationSharedPtr> m_frameVelFunction;
146 std::map<int, LibUtilities::EquationSharedPtr> m_frameDispFunction;
147 std::map<int, NekDouble> m_frameVel;
148 std::map<int, NekDouble> m_frameDisp;
149 void Update(NekDouble time)
150 {
151 for (auto it : m_frameVelFunction)
152 {
153 m_frameVel[it.first] = it.second->Evaluate(0, 0, 0, time);
154 }
155 for (auto it : m_frameDispFunction)
156 {
157 m_frameDisp[it.first] = it.second->Evaluate(0, 0, 0, time);
158 }
159 }
161 std::vector<NekDouble> m_box;
162 std::map<std::string, std::string> v_params;
163 unsigned int m_index = 0;
164 unsigned int m_outputFrequency;
165 std::string m_outputFile;
169 std::set<int> m_samplePointIDs;
171};
172
173} // namespace Nektar::SolverUtils
174
175#endif /* NEKTAR_SOLVERUTILS_FILTERS_FILTERLAGRANGIANPOINTS_H */
#define SOLVER_UTILS_EXPORT
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
static FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_ptr< EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
void ExtraPhysicsVars(std::vector< std::string > &extraVars)
SOLVER_UTILS_EXPORT ~FilterLagrangianPoints() override
struct Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame m_frame
std::map< std::string, std::string > v_params
void v_Finalise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
SOLVER_UTILS_EXPORT void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
SOLVER_UTILS_EXPORT FilterLagrangianPoints(const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_ptr< EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
void v_ModifyVelocity(Array< OneD, NekDouble > gcoords, NekDouble time, Array< OneD, NekDouble > vel) override
void GetPhysicsData(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::vector< Array< OneD, NekDouble > > &PhysicsData)
static StationaryPointsSharedPtr create(int rank, int dim, int intOrder, int idOffset, NekDouble dt, const std::vector< int > &Np, const std::vector< NekDouble > &Box, std::vector< std::string > extraVars)
Creates an instance of this class.
StatLagrangianPoints(int rank, int dim, int intOrder, int idOffset, NekDouble dt, const std::vector< int > &Np, const std::vector< NekDouble > &Box, std::vector< std::string > extraVars)
void v_GetCoords(int pid, Array< OneD, NekDouble > &gcoords) override
void v_SetPhysics(int pid, const Array< OneD, NekDouble > &data) override
void v_OutputData(std::string filename, bool verbose, std::map< std::string, NekDouble > &params) override
void v_SetCoords(int pid, const Array< OneD, NekDouble > &gcoords) override
void v_GetPhysics(int pid, Array< OneD, NekDouble > &data) override
Array< OneD, Array< OneD, NekDouble > > m_extraPhysics
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_velocity
void v_AssignPoint(int id, int pid, const Array< OneD, NekDouble > &gcoords) override
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_coords
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< StationaryPoints > StationaryPointsSharedPtr
std::shared_ptr< Filter > FilterSharedPtr
A shared pointer to a Driver object.
Definition: Filter.h:51
double NekDouble
std::map< int, LibUtilities::EquationSharedPtr > m_frameVelFunction
std::map< int, LibUtilities::EquationSharedPtr > m_frameDispFunction