Nektar++
EvaluatePoints.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: EvaluatePoints.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_EVALUATEPOINTS_H
36#define NEKTAR_SOLVERUTILS_FILTERS_EVALUATEPOINTS_H
37
44
46
47namespace Nektar::SolverUtils
48{
49
50class MobilePoint;
51class EvaluatePoints;
52class StationaryPoints;
53
54typedef std::shared_ptr<MobilePoint> MobilePointSharedPtr;
55typedef std::shared_ptr<StationaryPoints> StationaryPointsSharedPtr;
56
58{
59public:
60 virtual ~StationaryPoints() = default;
62 std::string filename, bool verbose,
63 std::map<std::string, NekDouble> &params)
64 {
65 v_OutputData(filename, verbose, params);
66 }
67 inline void ReSize(int Np)
68 {
69 v_ReSize(Np);
70 }
71 inline void AssignPoint(int id, int pid,
72 const Array<OneD, NekDouble> &gcoords)
73 {
74 v_AssignPoint(id, pid, gcoords);
75 }
76 inline void GetCoordsByPID(int pid, Array<OneD, NekDouble> &gcoords)
77 {
78 v_GetCoords(GlobalToLocal(pid), gcoords);
79 }
80 inline void SetCoordsByPID(int pid, const Array<OneD, NekDouble> &gcoords)
81 {
82 v_SetCoords(GlobalToLocal(pid), gcoords);
83 }
84 inline void GetPhysicsByPID(int pid, Array<OneD, NekDouble> &data)
85 {
86 v_GetPhysics(GlobalToLocal(pid), data);
87 }
88 inline void SetPhysicsByPID(int pid, const Array<OneD, NekDouble> &data)
89 {
90 v_SetPhysics(GlobalToLocal(pid), data);
91 }
92 inline void GetCoords(int pid, Array<OneD, NekDouble> &gcoords)
93 {
94 v_GetCoords(pid, gcoords);
95 }
96 inline void SetCoords(int pid, const Array<OneD, NekDouble> &gcoords)
97 {
98 v_SetCoords(pid, gcoords);
99 }
100 inline void GetPhysics(int pid, Array<OneD, NekDouble> &data)
101 {
102 v_GetPhysics(pid, data);
103 }
104 inline void SetPhysics(int pid, const Array<OneD, NekDouble> &data)
105 {
106 v_SetPhysics(pid, data);
107 }
108 inline void TimeAdvance(int order)
109 {
110 v_TimeAdvance(order);
111 }
112 inline int GetTotPoints()
113 {
114 return m_totPts;
115 }
116 inline int GetDim()
117 {
118 return m_dim;
119 }
120 inline int LocalToGlobal(int id)
121 {
122 return m_localIDToGlobal[id];
123 }
124 inline int GlobalToLocal(int id)
125 {
126 return m_globalIDToLocal[id];
127 }
128
129protected:
130 virtual void v_OutputData(std::string filename, bool verbose,
131 std::map<std::string, NekDouble> &params) = 0;
132 virtual void v_TimeAdvance(int order) = 0;
133 virtual void v_GetCoords(int pid, Array<OneD, NekDouble> &gcoords) = 0;
134 virtual void v_SetCoords(int pid,
135 const Array<OneD, NekDouble> &gcoords) = 0;
136 virtual void v_GetPhysics(int pid, Array<OneD, NekDouble> &data) = 0;
137 virtual void v_SetPhysics(int pid, const Array<OneD, NekDouble> &data) = 0;
138 virtual void v_ReSize(int Np) = 0;
139 virtual void v_AssignPoint(int id, int pid,
140 const Array<OneD, NekDouble> &gcoords);
141 int m_dim;
143 std::map<int, int> m_localIDToGlobal; // size should be m_totPts
144 std::map<int, int> m_globalIDToLocal; // size should be m_totPts
145};
146
148{
149public:
150 friend class MemoryManager<MobilePoint>;
153 int m_eId;
156
157 /// Creates an instance of this class
158 static MobilePointSharedPtr create(int dim, int sRank,
159 const Array<OneD, NekDouble> global)
160 {
163 return p;
164 }
165
166 MobilePoint(int dim, int sRank, const Array<OneD, NekDouble> global)
167 : m_sRank(sRank)
168 {
170 Vmath::Vcopy(dim, global, 1, m_global, 1);
172 }
173
174 void SetLocalCoords(int eId, const Array<OneD, NekDouble> &local)
175 {
176 m_eId = eId;
177 Vmath::Vcopy(m_local.size(), local, 1, m_local, 1);
178 }
179
181 {
182 Vmath::Vcopy(m_global.size(), global, 1, m_global, 1);
183 }
184
185 void SetData(int nPhys, const Array<OneD, NekDouble> data)
186 {
187 if (m_data.size() < nPhys)
188 {
190 }
191 Vmath::Vcopy(nPhys, data, 1, m_data, 1);
192 }
193
195 {
196 return m_data;
197 }
198
200 {
201 return m_global;
202 }
203
205 {
206 return m_local;
207 }
208};
209
211{
212public:
214 virtual ~EvaluatePoints();
218 const MultiRegions::ExpListSharedPtr &pField,
219 std::vector<Array<OneD, NekDouble>> &PhysicsData, NekDouble time);
221 std::map<int, std::set<int>> &callbackUpdateMobCoords);
223 std::map<int, std::set<int>> &callbackUpdateMobCoords);
228 const NekDouble &time, std::vector<std::string> &defaultValues);
230
231protected:
233 std::map<int, Array<OneD, NekDouble>> &globalCoords);
234 void GatherMobilePhysics(std::map<int, Array<OneD, NekDouble>> &revData);
237 std::set<int> &notfound);
240 std::set<int> &notfound);
241 void SyncColumnComm();
243 const NekDouble time);
244 void Pack2Int(const int &a, const int &b, double &d);
245 void unPack2Int(int &a, int &b, const double &d);
246 void Pack3Int(const int &a, const int &b, const int &c, double &d);
247 void unPack3Int(int &a, int &b, int &c, const double &d);
248 void Pack2Short(const int &a, const int &b, int &c);
249 void unPack2Short(int &a, int &b, const int &c);
250 virtual void v_ModifyVelocity(Array<OneD, NekDouble> gcoords,
260 std::map<int, MobilePointSharedPtr> m_mobilePts;
261 std::map<int, LibUtilities::EquationSharedPtr> m_defaultValueFunction;
262 // m_recvMobInfo[t] = n, in the current stationary pts, there are n mobile
263 // points from thread t; all threads; global
264 // m_recvStatInfo[t] = n, in the
265 // current mobile pts, there are n stationary points from thread t;
266 // m_columnRank == 0 only; local
267 std::map<int, int> m_recvMobInfo;
268 std::map<int, int> m_recvStatInfo;
269};
270
271} // namespace Nektar::SolverUtils
272
273#endif /* NEKTAR_SOLVERUTILS_FILTERS_EVALUATEPOINTS_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.
void PartitionLocalPoints(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::set< int > &notfound)
SOLVER_UTILS_EXPORT void EvaluateMobilePhys(const MultiRegions::ExpListSharedPtr &pField, std::vector< Array< OneD, NekDouble > > &PhysicsData, NekDouble time)
SOLVER_UTILS_EXPORT void Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time, std::vector< std::string > &defaultValues)
NekDouble GetDefaultValue(const int i, const Array< OneD, const NekDouble > x, const NekDouble time)
SOLVER_UTILS_EXPORT void PassMobilePhysToStatic(std::map< int, std::set< int > > &callbackUpdateMobCoords)
void Pack2Int(const int &a, const int &b, double &d)
SOLVER_UTILS_EXPORT void PartitionMobilePts(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
void unPack2Short(int &a, int &b, const int &c)
SOLVER_UTILS_EXPORT void SetUpCommInfo()
SOLVER_UTILS_EXPORT void CopyMobilePtsToStatic()
SOLVER_UTILS_EXPORT void CopyStaticPtsToMobile()
void Pack2Short(const int &a, const int &b, int &c)
StationaryPointsSharedPtr m_staticPts
virtual void v_ModifyVelocity(Array< OneD, NekDouble > gcoords, NekDouble time, Array< OneD, NekDouble > vel)
void PartitionExchangeNonlocalPoints(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::set< int > &notfound)
void GatherMobilePhysics(std::map< int, Array< OneD, NekDouble > > &revData)
void Pack3Int(const int &a, const int &b, const int &c, double &d)
void unPack2Int(int &a, int &b, const double &d)
void UpdateMobileCoords(std::map< int, Array< OneD, NekDouble > > &globalCoords)
std::map< int, MobilePointSharedPtr > m_mobilePts
LibUtilities::CommSharedPtr m_comm
std::map< int, LibUtilities::EquationSharedPtr > m_defaultValueFunction
SOLVER_UTILS_EXPORT void PassStaticCoordsToMobile(std::map< int, std::set< int > > &callbackUpdateMobCoords)
void unPack3Int(int &a, int &b, int &c, const double &d)
Array< OneD, NekDouble > m_data
Array< OneD, NekDouble > m_global
static MobilePointSharedPtr create(int dim, int sRank, const Array< OneD, NekDouble > global)
Creates an instance of this class.
void SetLocalCoords(int eId, const Array< OneD, NekDouble > &local)
MobilePoint(int dim, int sRank, const Array< OneD, NekDouble > global)
Array< OneD, NekDouble > GetData()
Array< OneD, NekDouble > m_local
Array< OneD, NekDouble > GetLocalCoords()
void SetGlobalCoords(const Array< OneD, NekDouble > &global)
Array< OneD, NekDouble > GetGlobalCoords()
void SetData(int nPhys, const Array< OneD, NekDouble > data)
void SetCoordsByPID(int pid, const Array< OneD, NekDouble > &gcoords)
void GetCoordsByPID(int pid, Array< OneD, NekDouble > &gcoords)
virtual void v_GetPhysics(int pid, Array< OneD, NekDouble > &data)=0
virtual void v_SetPhysics(int pid, const Array< OneD, NekDouble > &data)=0
void GetPhysicsByPID(int pid, Array< OneD, NekDouble > &data)
void SetPhysics(int pid, const Array< OneD, NekDouble > &data)
virtual void v_SetCoords(int pid, const Array< OneD, NekDouble > &gcoords)=0
virtual void v_TimeAdvance(int order)=0
void AssignPoint(int id, int pid, const Array< OneD, NekDouble > &gcoords)
SOLVER_UTILS_EXPORT void OutputData(std::string filename, bool verbose, std::map< std::string, NekDouble > &params)
void SetPhysicsByPID(int pid, const Array< OneD, NekDouble > &data)
void GetPhysics(int pid, Array< OneD, NekDouble > &data)
virtual void v_GetCoords(int pid, Array< OneD, NekDouble > &gcoords)=0
virtual void v_OutputData(std::string filename, bool verbose, std::map< std::string, NekDouble > &params)=0
void GetCoords(int pid, Array< OneD, NekDouble > &gcoords)
void SetCoords(int pid, const Array< OneD, NekDouble > &gcoords)
virtual void v_AssignPoint(int id, int pid, const Array< OneD, NekDouble > &gcoords)
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< StationaryPoints > StationaryPointsSharedPtr
std::shared_ptr< MobilePoint > MobilePointSharedPtr
std::vector< double > d(NPUPPER *NPUPPER)
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825