Nektar++
PulseWavePropagation.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: PulseWavePropagation.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: Pulse Wave Propagation solve routines based on the weak
33// formulation (1):
34//
35///////////////////////////////////////////////////////////////////////////////
36
37#include <iostream>
38
41
42using namespace std;
43
44namespace Nektar
45{
46
49 "PulseWavePropagation", PulseWavePropagation::create,
50 "Pulse Wave Propagation equation.");
51/**
52 * @class PulseWavePropagation
53 *
54 * Set up the routines based on the weak formulation from
55 * "Computational Modelling of 1D blood flow with variable
56 * mechanical properties" by S. J. Sherwin et al. The weak
57 * formulation (1) reads:
58 * \f$ \sum_{e=1}^{N_{el}} \left[ \left( \frac{\partial \mathbf{U}^{\delta}
59 * }{\partial t} , \mathbf{\psi}^{\delta} \right)_{\Omega_e} - \left(
60 * \frac{\partial \mathbf{F(\mathbf{U})}^{\delta} }
61 * {\partial x}, \mathbf{\psi}^{\delta} \right)_{\Omega_e} + \left[
62 * \mathbf{\psi}^{\delta} \cdot \{ \mathbf{F}^u -
63 * \mathbf{F}(\mathbf{U}^{\delta}) \} \right]_{x_e^l}^{x_eû} \right] = 0 \f$
64 */
68 : PulseWaveSystem(pSession, pGraph)
69{
70}
71
72void PulseWavePropagation::v_InitObject([[maybe_unused]] bool DeclareField)
73{
74 // Will set up an array of vessels/fields in PulseWaveSystem::v_InitObject
75 // so set DeclareField to false so that the fields are not set up in
76 // EquationSystem unnecessarily. Note the number of fields in Equation
77 // system is related to the number of variables. The number of vessels is
78 // therefore held in PulwWaveSystem.
80
81 if (m_session->DefinesSolverInfo("PressureArea"))
82 {
84 m_session->GetSolverInfo("PressureArea"), m_vessels, m_session);
85 }
86 else
87 {
89 "Beta", m_vessels, m_session);
90 }
91
93 {
96 }
97 else
98 {
99 ASSERTL0(false, "Implicit Pulse Wave Propagation not set up.");
100 }
101
102 // Create advection object
103 string advName;
104 string riemName;
105 switch (m_upwindTypePulse)
106 {
107 case eUpwindPulse:
108 {
109 advName = "WeakDG";
110 riemName = "UpwindPulse";
111 }
112 break;
113 default:
114 {
115 ASSERTL0(false, "populate switch statement for upwind flux");
116 }
117 break;
118 }
123 riemName, m_session);
124 m_riemannSolver->SetScalar("A0", &PulseWavePropagation::GetA0, this);
125 m_riemannSolver->SetScalar("beta", &PulseWavePropagation::GetBeta, this);
126 m_riemannSolver->SetScalar("alpha", &PulseWavePropagation::GetAlpha, this);
127 m_riemannSolver->SetScalar("N", &PulseWavePropagation::GetN, this);
128 m_riemannSolver->SetParam("rho", &PulseWavePropagation::GetRho, this);
130 this);
131
132 m_advObject->SetRiemannSolver(m_riemannSolver);
133 m_advObject->InitObject(m_session, m_fields);
134}
135
137{
138}
139
140/**
141 * Computes the right hand side of (1). The RHS is everything
142 * except the term that contains the time derivative
143 * \f$\frac{\partial \mathbf{U}}{\partial t}\f$. In case of a
144 * Discontinuous Galerkin projection, m_advObject->Advect
145 * will be called
146 *
147 */
149 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
150 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
151{
152 size_t i;
153
155
156 // Dummy array for WeakDG advection
158
159 // Output array for advection
161
162 size_t cnt = 0;
163
164 // Set up Inflow and Outflow boundary conditions.
165 SetPulseWaveBoundaryConditions(inarray, outarray, time);
166
167 // Set up any interface conditions and write into boundary condition
169
170 // do advection evaluation in all domains
171 for (size_t omega = 0; omega < m_nDomains; ++omega)
172 {
174 m_currentDomain = omega;
175 size_t nq = m_vessels[omega * m_nVariables]->GetTotPoints();
176
177 timer.Start();
178 for (i = 0; i < m_nVariables; ++i)
179 {
180 physarray[i] = inarray[i] + cnt;
181 out[i] = outarray[i] + cnt;
182 }
183
184 for (i = 0; i < m_nVariables; ++i)
185 {
186 m_fields[i] = m_vessels[omega * m_nVariables + i];
187 }
188
189 m_advObject->Advect(m_nVariables, m_fields, advVel, physarray, out,
190 time);
191 for (i = 0; i < m_nVariables; ++i)
192 {
193 Vmath::Neg(nq, out[i], 1);
194 }
195 timer.Stop();
196 timer.AccumulateRegion("PulseWavePropagation:_DoOdeRHS", 1);
197 cnt += nq;
198 }
199}
200
202 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
204 [[maybe_unused]] const NekDouble time)
205{
206 // Just copy over array
207 if (inarray != outarray)
208 {
209 for (size_t i = 0; i < m_nVariables; ++i)
210 {
211 Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
212 }
213 }
214}
215
216/**
217 * Does the projection between ... space and the ... space. Also checks for
218 *Q-inflow boundary conditions at the inflow of the current arterial segment and
219 *applies the Q-inflow if specified
220 */
222 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
223 [[maybe_unused]] Array<OneD, Array<OneD, NekDouble>> &outarray,
224 const NekDouble time)
225
226{
227 size_t omega;
228
230
231 size_t offset = 0;
232
233 // This will be moved to the RCR boundary condition once factory is setup
234 if (time == 0)
235 {
237
238 for (omega = 0; omega < m_nDomains; ++omega)
239 {
240 vessel[0] = m_vessels[2 * omega];
241 vessel[1] = m_vessels[2 * omega + 1];
242
243 for (size_t j = 0; j < 2; ++j)
244 {
245 std::string BCType;
246
247 if (j < vessel[0]->GetBndConditions().size())
248 {
249 BCType = vessel[0]->GetBndConditions()[j]->GetUserDefined();
250 }
251
252 // if no condition given define it to be NoUserDefined
253 if (BCType.empty() || BCType == "Interface")
254 {
255 BCType = "NoUserDefined";
256 }
257
260
261 // turn on time dependent BCs
262 if (BCType == "Q-inflow")
263 {
264 vessel[0]->GetBndConditions()[j]->SetIsTimeDependent(true);
265 }
266 else if (BCType == "A-inflow")
267 {
268 vessel[0]->GetBndConditions()[j]->SetIsTimeDependent(true);
269 }
270 else if (BCType == "U-inflow")
271 {
272 vessel[1]->GetBndConditions()[j]->SetIsTimeDependent(true);
273 }
274 else if (BCType == "RCR-terminal")
275 {
276 vessel[0]->GetBndConditions()[j]->SetIsTimeDependent(true);
277 }
278 }
279 }
280 }
281
283
284 // Loop over all vessels and set boundary conditions
286 for (omega = 0; omega < m_nDomains; ++omega)
287 {
288 timer.Start();
289 for (size_t n = 0; n < 2; ++n)
290 {
291 m_Boundary[2 * omega + n]->DoBoundary(
292 inarray, m_A_0, m_beta, m_alpha, time, omega, offset, n);
293 }
294
295 offset += m_vessels[2 * omega]->GetTotPoints();
296 timer.Stop();
297 timer.AccumulateRegion("PulseWavePropagation:_SetBCs", 1);
298 }
299}
300
301/**
302 * Calculates the second term of the weak form (1): \f$
303 * \left( \frac{\partial \mathbf{F(\mathbf{U})}^{\delta}
304 * }{\partial x}, \mathbf{\psi}^{\delta} \right)_{\Omega_e}
305 * \f$
306 * The variables of the system are $\mathbf{U} = [A,u]^T$
307 * physfield[0] = A physfield[1] = u
308 * flux[0] = F[0] = A*u flux[1] = F[1] = u^2/2 + p/rho
309 */
311 const Array<OneD, Array<OneD, NekDouble>> &physfield,
313{
314 size_t nq = m_vessels[m_currentDomain * m_nVariables]->GetTotPoints();
315 NekDouble domain = m_currentDomain;
317 Array<OneD, NekDouble> dAUdx(nq);
318 NekDouble viscoelasticGradient = 0.0;
319
321
322 for (size_t j = 0; j < nq; ++j)
323 {
324 timer.Start();
325 flux[0][0][j] = physfield[0][j] * physfield[1][j];
326 timer.Stop();
327 timer.AccumulateRegion("PulseWavePropagation:GetFluxVector-flux", 3);
328 }
329
330 // d/dx of AU, for the viscoelastic tube law and extra fields
331 m_fields[0]->PhysDeriv(flux[0][0], dAUdx);
332
333 for (size_t j = 0; j < nq; ++j)
334 {
335 if ((j == 0) || (j == nq - 1))
336 {
337 viscoelasticGradient = dAUdx[j];
338 }
339 else
340 {
341 viscoelasticGradient = (dAUdx[j] + dAUdx[j + 1]) / 2;
342 }
343
344 m_pressureArea->GetPressure(m_pressure[domain][j], m_beta[domain][j],
345 physfield[0][j], m_A_0[domain][j],
346 viscoelasticGradient, m_gamma[domain][j],
347 m_alpha[domain][j]);
348
349 flux[1][0][j] = physfield[1][j] * physfield[1][j] / 2 +
350 m_pressure[domain][j] / m_rho;
351 }
352
353 m_session->MatchSolverInfo("OutputExtraFields", "True", extraFields, true);
354
355 if (extraFields)
356 {
357 /*
358 Calculates wave speed and characteristic variables.
359
360 Ideally this should be moved to PulseWaveSystem, but it's easiest to
361 implement here.
362 */
363 size_t counter = 0;
364
365 m_PWV[domain] = Array<OneD, NekDouble>(nq);
366 m_W1[domain] = Array<OneD, NekDouble>(nq);
367 m_W2[domain] = Array<OneD, NekDouble>(nq);
368
369 for (size_t j = 0; j < nq; ++j)
370 {
371 m_pressureArea->GetC(m_PWV[domain][j], m_beta[domain][j],
372 physfield[0][counter + j], m_A_0[domain][j],
373 m_alpha[domain][j]);
374 m_pressureArea->GetW1(m_W1[domain][j], physfield[1][counter + j],
375 m_beta[domain][j], physfield[0][counter + j],
376 m_A_0[domain][j], m_alpha[domain][j]);
377 m_pressureArea->GetW2(m_W2[domain][j], physfield[1][counter + j],
378 m_beta[domain][j], physfield[0][counter + j],
379 m_A_0[domain][j], m_alpha[domain][j]);
380 }
381
382 counter += nq;
383 }
384}
385
387{
389}
390
392{
394}
395
397{
399}
400
402{
404}
405
407{
408 return m_rho;
409}
410
412{
413 return m_nDomains;
414}
415
416/**
417 * Print summary routine, calls virtual routine reimplemented in
418 * UnsteadySystem
419 */
421{
423}
424
425} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
Definition: Timer.cpp:70
Array< OneD, NekDouble > & GetAlpha()
void v_InitObject(bool DeclareField=false) override
Init object for UnsteadySystem class.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
PulseWavePropagation(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
SolverUtils::AdvectionSharedPtr m_advObject
void SetPulseWaveBoundaryConditions(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Array< OneD, PulseWaveBoundarySharedPtr > m_Boundary
Array< OneD, NekDouble > & GetN()
Array< OneD, NekDouble > & GetA0()
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
DG Pulse Wave Propagation routines:
void v_GenerateSummary(SolverUtils::SummaryList &s) override
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Array< OneD, NekDouble > & GetBeta()
static std::string className
Name of class.
static EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Base class for unsteady solvers.
Array< OneD, Array< OneD, NekDouble > > m_A_0
Array< OneD, Array< OneD, NekDouble > > m_beta_trace
Array< OneD, Array< OneD, NekDouble > > m_W2
Array< OneD, Array< OneD, NekDouble > > m_trace_fwd_normal
PulseWavePressureAreaSharedPtr m_pressureArea
Array< OneD, Array< OneD, NekDouble > > m_W1
void v_InitObject(bool DeclareField=false) override
Array< OneD, Array< OneD, NekDouble > > m_alpha_trace
UpwindTypePulse m_upwindTypePulse
Array< OneD, Array< OneD, NekDouble > > m_A_0_trace
Array< OneD, Array< OneD, NekDouble > > m_pressure
Array< OneD, Array< OneD, NekDouble > > m_gamma
Array< OneD, Array< OneD, NekDouble > > m_alpha
Array< OneD, Array< OneD, NekDouble > > m_PWV
Array< OneD, MultiRegions::ExpListSharedPtr > m_vessels
void EnforceInterfaceConditions(const Array< OneD, const Array< OneD, NekDouble > > &fields)
Array< OneD, Array< OneD, NekDouble > > m_beta
int m_spacedim
Spatial dimension (>= expansion dim).
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
EquationSystemFactory & GetEquationSystemFactory()
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
PressureAreaFactory & GetPressureAreaFactory()
@ eUpwindPulse
simple upwinding scheme
BoundaryFactory & GetBoundaryFactory()
double NekDouble
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825