Nektar++
ForcingMovingReferenceFrame.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ForcingMovingReferenceFrame.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: Allows for a moving frame of reference, through adding c * du/dx
32// to the body force, where c is the frame velocity vector
33//
34///////////////////////////////////////////////////////////////////////////////
35
36#ifndef NEKTAR_SOLVERUTILS_FORCINGMOVINGREFERENCEFRAME
37#define NEKTAR_SOLVERUTILS_FORCINGMOVINGREFERENCEFRAME
38
39#include <string>
40
49#include <boost/numeric/ublas/io.hpp>
50#include <boost/numeric/ublas/matrix.hpp>
51#include <boost/numeric/ublas/vector.hpp>
52#include <cmath>
53
54namespace Nektar::SolverUtils
55{
56
57namespace bn = boost::numeric;
58
59/***
60 * Solve the body's motion using Newmark-Beta method
61 * M ddx + C dx + K x = F
62 * In discrete form
63 * CoeffMatrix dx = rhs
64 ***/
66{
67public:
72 Array<OneD, NekDouble> K, std::set<int> DirDoFs);
74 std::map<int, NekDouble> motionPrescribed);
79 std::map<int, NekDouble> motionPrescribed);
80 int m_rows;
82 std::vector<int> m_index;
88};
89
91{
92public:
95 void SetAngle(const Array<OneD, NekDouble> theta);
96 void BodyToInerital(const int dim, const Array<OneD, NekDouble> &body,
97 Array<OneD, NekDouble> &inertial);
98 void IneritalToBody(const int dim, const Array<OneD, NekDouble> &inertial,
100
101private:
103};
104
106{
107
108public:
110
111 /// Creates an instance of this class
114 const std::weak_ptr<EquationSystem> &pEquation,
116 const unsigned int &pNumForcingFields, const TiXmlElement *pForce)
117 {
120 pSession, pEquation);
121 p->InitObject(pFields, pNumForcingFields, pForce);
122 return p;
123 }
124
125 /// Name of the class
126 static std::string classNameBody;
127
128protected:
131 const unsigned int &pNumForcingFields,
132 const TiXmlElement *pForce) override;
133
136 const Array<OneD, Array<OneD, NekDouble>> &inarray,
138 const NekDouble &time) override;
139
142 const Array<OneD, Array<OneD, NekDouble>> &inarray,
144 const NekDouble &time) override;
145
146private:
147 // name of the function for linear and angular velocities in the session
148 // file
149 // pivot point
153
154 // prescribed functions in the session file
155 std::map<int, LibUtilities::EquationSharedPtr> m_frameVelFunction;
156 std::ofstream m_outputStream;
158
159 // a boolean switch indicating for which direction the velocities are
160 // available. The available velocites could be different from the
161 // precscribed one because of the rotation which result in change of basis
162 // vector of local frame to the inertial frame.
165 bool m_hasRotation; // m_hasOmega[0] || m_hasOmega[1] || m_hasOmega[2]
166
167 // frame linear velocities in local translating-rotating frame
169
170 // frame angular velocities in local translating-rotating frame
172 // coordinate vector
174
177
183 unsigned int m_index;
184 unsigned int m_outputFrequency;
185
186 struct
187 {
190 std::set<int> dirDoFs;
192 // fluid force filter
194 // externel force
195 std::map<int, LibUtilities::EquationSharedPtr> extForceFunction;
202
205 const std::weak_ptr<EquationSystem> &pEquation);
206
207 ~ForcingMovingReferenceFrame(void) override;
208
209 void UpdatePrescribed(const NekDouble &time,
210 std::map<int, NekDouble> &Dirs);
213 const NekDouble &time);
215 const int step);
217 void SetInitialConditions(std::map<int, NekDouble> &Dirs);
218
219 void addRotation(int npoints,
220 const Array<OneD, Array<OneD, NekDouble>> &inarray0,
221 NekDouble angVelScale,
222 const Array<OneD, Array<OneD, NekDouble>> &inarray1,
223 Array<OneD, Array<OneD, NekDouble>> &outarray);
224 void InitBodySolver(const TiXmlElement *pForce);
226 const Array<OneD, NekDouble> &forcebody,
227 std::map<int, NekDouble> &Dirs);
228 void LoadParameters(const TiXmlElement *pForce);
229 NekDouble EvaluateExpression(std::string expression);
230 void InitialiseFilter(
233 const TiXmlElement *pForce);
235};
236
237} // namespace Nektar::SolverUtils
238
239#endif
#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.
Defines a forcing term to be explicitly applied.
Definition: Forcing.h:71
SOLVER_UTILS_EXPORT void v_PreApply(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time) override
Compute the moving frame velocity at given time.
void addRotation(int npoints, const Array< OneD, Array< OneD, NekDouble > > &inarray0, NekDouble angVelScale, const Array< OneD, Array< OneD, NekDouble > > &inarray1, Array< OneD, Array< OneD, NekDouble > > &outarray)
outarray = inarray0 + angVelScale Omega x inarray1
std::map< int, LibUtilities::EquationSharedPtr > m_frameVelFunction
void SolveBodyMotion(Array< OneD, Array< OneD, NekDouble > > &bodyVel, const Array< OneD, NekDouble > &forcebody, std::map< int, NekDouble > &Dirs)
void UpdateFluidInterface(Array< OneD, Array< OneD, NekDouble > > &bodyVel, const int step)
void UpdatePrescribed(const NekDouble &time, std::map< int, NekDouble > &Dirs)
ForcingMovingReferenceFrame(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation)
SOLVER_UTILS_EXPORT void v_Apply(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time) override
Adds the body force, -Omega X u.
std::map< int, LibUtilities::EquationSharedPtr > extForceFunction
SOLVER_UTILS_EXPORT void v_InitObject(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce) override
Initialise the forcing module.
void UpdateFrameVelocity(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
Updates the forcing array with the current required forcing.
void UpdateBoundaryConditions(NekDouble time)
Set velocity boundary condition at the next time step.
static SOLVER_UTILS_EXPORT ForcingSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce)
Creates an instance of this class.
struct Nektar::SolverUtils::ForcingMovingReferenceFrame::@2 m_body
void InitialiseFilter(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pForce)
void SetAngle(const Array< OneD, NekDouble > theta)
void BodyToInerital(const int dim, const Array< OneD, NekDouble > &body, Array< OneD, NekDouble > &inertial)
void IneritalToBody(const int dim, const Array< OneD, NekDouble > &inertial, Array< OneD, NekDouble > &body)
void SetNewmarkBeta(NekDouble beta, NekDouble gamma, NekDouble dt, Array< OneD, NekDouble > M, Array< OneD, NekDouble > C, Array< OneD, NekDouble > K, std::set< int > DirDoFs)
Array< OneD, Array< OneD, NekDouble > > m_K
void SolveFree(Array< OneD, Array< OneD, NekDouble > > u, Array< OneD, NekDouble > force)
void SolvePrescribed(Array< OneD, Array< OneD, NekDouble > > u, std::map< int, NekDouble > motionPrescribed)
Array< OneD, Array< OneD, NekDouble > > m_Matrix
Array< OneD, Array< OneD, NekDouble > > m_C
void Solve(Array< OneD, Array< OneD, NekDouble > > u, Array< OneD, NekDouble > force, std::map< int, NekDouble > motionPrescribed)
Array< OneD, Array< OneD, NekDouble > > m_M
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
std::shared_ptr< FilterAeroForces > FilterAeroForcesSharedPtr
SOLVER_UTILS_EXPORT typedef std::shared_ptr< Forcing > ForcingSharedPtr
A shared pointer to an EquationSystem object.
Definition: Forcing.h:53
std::int32_t NekInt
double NekDouble