Nektar++
PulseWavePressureArea.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: PulseWavePressureArea.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// 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: PulseWavePressureArea header
33//
34///////////////////////////////////////////////////////////////////////////////
35
36#ifndef NEKTAR_PULSEWAVEPRESSUREAREA_H
37#define NEKTAR_PULSEWAVEPRESSUREAREA_H
38
43
44namespace Nektar
45{
46
47class PulseWavePressureArea;
48
49typedef std::shared_ptr<PulseWavePressureArea> PulseWavePressureAreaSharedPtr;
50
52
58
60{
61public:
64
65 virtual ~PulseWavePressureArea() = default;
66
67 inline void GetPressure(NekDouble &P, const NekDouble &beta,
68 const NekDouble &A, const NekDouble &A0,
69 const NekDouble &dAUdx, const NekDouble &gamma = 0,
70 const NekDouble &alpha = 0.5)
71 {
72 v_GetPressure(P, beta, A, A0, dAUdx, gamma, alpha);
73 }
74
75 inline void GetC(NekDouble &c, const NekDouble &beta, const NekDouble &A,
76 const NekDouble &A0, const NekDouble &alpha = 0.5)
77 {
78 v_GetC(c, beta, A, A0, alpha);
79 }
80
81 inline void GetW1(NekDouble &W1, const NekDouble &u, const NekDouble &beta,
82 const NekDouble &A, const NekDouble &A0,
83 const NekDouble &alpha = 0.5)
84 {
85 v_GetW1(W1, u, beta, A, A0, alpha);
86 }
87
88 inline void GetW2(NekDouble &W2, const NekDouble &u, const NekDouble &beta,
89 const NekDouble &A, const NekDouble &A0,
90 const NekDouble &alpha = 0.5)
91 {
92 v_GetW2(W2, u, beta, A, A0, alpha);
93 }
94
95 inline void GetAFromChars(NekDouble &A, const NekDouble &W1,
96 const NekDouble &W2, const NekDouble &beta,
97 const NekDouble &A0, const NekDouble &alpha = 0.5)
98 {
99 v_GetAFromChars(A, W1, W2, beta, A0, alpha);
100 }
101
102 inline void GetUFromChars(NekDouble &u, const NekDouble &W1,
103 const NekDouble &W2)
104 {
105 v_GetUFromChars(u, W1, W2);
106 }
107
108 inline void GetCharIntegral(NekDouble &I, const NekDouble &beta,
109 const NekDouble &A, const NekDouble &A0,
110 const NekDouble &alpha = 0.5)
111 {
112 v_GetCharIntegral(I, beta, A, A0, alpha);
113 }
114
116 const Array<OneD, NekDouble> &Au,
117 const Array<OneD, NekDouble> &uu,
119 const Array<OneD, NekDouble> &A0,
120 const Array<OneD, NekDouble> &alpha,
121 const std::string &type)
122 {
123 v_GetJacobianInverse(invJ, Au, uu, beta, A0, alpha, type);
124 }
125
126protected:
129
132
133 virtual void v_GetPressure(NekDouble &P, const NekDouble &beta,
134 const NekDouble &A, const NekDouble &A0,
135 const NekDouble &dAUdx,
136 const NekDouble &gamma = 0,
137 const NekDouble &alpha = 0.5) = 0;
138
139 virtual void v_GetC(NekDouble &c, const NekDouble &beta, const NekDouble &A,
140 const NekDouble &A0, const NekDouble &alpha = 0.5) = 0;
141
142 virtual void v_GetW1(NekDouble &W1, const NekDouble &u,
143 const NekDouble &beta, const NekDouble &A,
144 const NekDouble &A0, const NekDouble &alpha = 0.5) = 0;
145
146 virtual void v_GetW2(NekDouble &W2, const NekDouble &u,
147 const NekDouble &beta, const NekDouble &A,
148 const NekDouble &A0, const NekDouble &alpha = 0.5) = 0;
149
150 virtual void v_GetAFromChars(NekDouble &A, const NekDouble &W1,
151 const NekDouble &W2, const NekDouble &beta,
152 const NekDouble &A0,
153 const NekDouble &alpha = 0.5) = 0;
154
155 virtual void v_GetUFromChars(NekDouble &u, const NekDouble &W1,
156 const NekDouble &W2) = 0;
157
158 virtual void v_GetCharIntegral(NekDouble &I, const NekDouble &beta,
159 const NekDouble &A, const NekDouble &A0,
160 const NekDouble &alpha = 0.5) = 0;
161
163 const Array<OneD, NekDouble> &Au,
164 const Array<OneD, NekDouble> &uu,
166 const Array<OneD, NekDouble> &A0,
167 const Array<OneD, NekDouble> &alpha,
168 const std::string &type) = 0;
169
170private:
171};
172
173} // namespace Nektar
174#endif
Provides a generic Factory class.
void GetCharIntegral(NekDouble &I, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
void GetW2(NekDouble &W2, const NekDouble &u, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
void GetC(NekDouble &c, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
virtual void v_GetW1(NekDouble &W1, const NekDouble &u, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)=0
virtual void v_GetJacobianInverse(NekMatrix< NekDouble > &invJ, const Array< OneD, NekDouble > &Au, const Array< OneD, NekDouble > &uu, const Array< OneD, NekDouble > &beta, const Array< OneD, NekDouble > &A0, const Array< OneD, NekDouble > &alpha, const std::string &type)=0
virtual void v_GetCharIntegral(NekDouble &I, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)=0
void GetW1(NekDouble &W1, const NekDouble &u, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)
virtual void v_GetUFromChars(NekDouble &u, const NekDouble &W1, const NekDouble &W2)=0
Array< OneD, MultiRegions::ExpListSharedPtr > m_vessels
virtual void v_GetPressure(NekDouble &P, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &dAUdx, const NekDouble &gamma=0, const NekDouble &alpha=0.5)=0
void GetJacobianInverse(NekMatrix< NekDouble > &invJ, const Array< OneD, NekDouble > &Au, const Array< OneD, NekDouble > &uu, const Array< OneD, NekDouble > &beta, const Array< OneD, NekDouble > &A0, const Array< OneD, NekDouble > &alpha, const std::string &type)
LibUtilities::SessionReaderSharedPtr m_session
virtual ~PulseWavePressureArea()=default
virtual void v_GetC(NekDouble &c, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)=0
void GetUFromChars(NekDouble &u, const NekDouble &W1, const NekDouble &W2)
virtual void v_GetW2(NekDouble &W2, const NekDouble &u, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &alpha=0.5)=0
void GetAFromChars(NekDouble &A, const NekDouble &W1, const NekDouble &W2, const NekDouble &beta, const NekDouble &A0, const NekDouble &alpha=0.5)
PulseWavePressureArea(Array< OneD, MultiRegions::ExpListSharedPtr > &pVessel, const LibUtilities::SessionReaderSharedPtr &pSession)
virtual void v_GetAFromChars(NekDouble &A, const NekDouble &W1, const NekDouble &W2, const NekDouble &beta, const NekDouble &A0, const NekDouble &alpha=0.5)=0
void GetPressure(NekDouble &P, const NekDouble &beta, const NekDouble &A, const NekDouble &A0, const NekDouble &dAUdx, const NekDouble &gamma=0, const NekDouble &alpha=0.5)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
@ P
Monomial polynomials .
Definition: BasisType.h:62
std::shared_ptr< PulseWavePressureArea > PulseWavePressureAreaSharedPtr
LibUtilities::NekFactory< std::string, PulseWavePressureArea, Array< OneD, MultiRegions::ExpListSharedPtr > &, const LibUtilities::SessionReaderSharedPtr & > PressureAreaFactory
PressureAreaFactory & GetPressureAreaFactory()
static PulseWavePressureAreaSharedPtr NullPulseWavePressureAreaSharedPtr
double NekDouble