Skip to content

Commit 4f59f04

Browse files
authored
Merge pull request #3 from yolhan83/develop
Add second order 1D model
2 parents 4c57e89 + 0848953 commit 4f59f04

File tree

4 files changed

+72
-4
lines changed

4 files changed

+72
-4
lines changed

README.md

Lines changed: 0 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -61,7 +61,6 @@ pkg> add https://github.com/your-repo/BloodFlowTrixi.jl
6161

6262
**short term**
6363
- Add second order 1D model.
64-
- Design prim variables for 1D and 2D models.
6564
- Add proper tests for 1D and 2D models.
6665
- Add 3D representations of the solutions for 1D and 2D models.
6766
- Design easy to use interfaces for users to define their own initial and boundary conditions and source terms.

docs/src/math.md

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -56,10 +56,10 @@ Le modèle 2D est dérivé à partir d’une **intégration radiale des équatio
5656

5757
2. **Conservation de la quantité de mouvement (composante radiale et axiale)** :
5858
```math
59-
∂_t (Q_{Rθ}) + ∂_θ \left( \frac{Q_{Rθ}^2}{2 A^2} + A P \right) + ∂_s \left( \frac{Q_{Rθ} Q_s}{A} \right) = \frac{2 R}{3} C \sin θ \frac{Q_s^2}{A} + \frac{2 R k Q_{Rθ}}{A} + ∂_θ (A P)
59+
∂_t (Q_{Rθ}) + ∂_θ \left( \frac{Q_{Rθ}^2}{2 A^2} + A P \right) + ∂_s \left( \frac{Q_{Rθ} Q_s}{A} \right) = \frac{2 R}{3} C \sin θ \frac{Q_s^2}{A} + \frac{2 R k Q_{Rθ}}{A} + P∂_θ (A)
6060
```
6161
```math
62-
∂_t (Q_s) + ∂_θ \left( \frac{Q_s Q_{Rθ}}{A^2} \right) + ∂_s \left( \frac{Q_s^2}{A} - \frac{Q_{Rθ}^2}{2 A^2} + A P \right) = - \frac{2 R}{3} C \sin θ \frac{Q_{Rθ} Q_s}{A^2} + \frac{k R Q_s}{A} + ∂_s (A P)
62+
∂_t (Q_s) + ∂_θ \left( \frac{Q_s Q_{Rθ}}{A^2} \right) + ∂_s \left( \frac{Q_s^2}{A} - \frac{Q_{Rθ}^2}{2 A^2} + A P \right) = - \frac{2 R}{3} C \sin θ \frac{Q_{Rθ} Q_s}{A^2} + \frac{k R Q_s}{A} + P∂_s (A)
6363
```
6464

6565
### Énergie et relation d’entropie du modèle 2D

src/1DModel/1dmodel.jl

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -148,4 +148,5 @@ end
148148
include("./variables.jl")
149149
include("./bc1d.jl")
150150
include("./Test_Cases/pressure_in.jl")
151-
include("./Test_Cases/convergence_test.jl")
151+
include("./Test_Cases/convergence_test.jl")
152+
include("Ord2/1dmodelord2.jl")

src/1DModel/Ord2/1dmodelord2.jl

Lines changed: 68 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,68 @@
1+
struct BloodFlowEquations1DOrd2{T <:Real,E} <: Trixi.AbstractEquationsParabolic{1, 4, GradientVariablesConservative}
2+
nu ::T
3+
model1d ::E
4+
end
5+
Trixi.varnames(mapin,eq::BloodFlowTrixi.BloodFlowEquations1DOrd2) = Trixi.varnames(mapin,eq.model1d)
6+
7+
function Trixi.flux(u,gradients,orientation::Int,eq_parab ::BloodFlowEquations1DOrd2)
8+
dudx = gradients
9+
a,Q,_,A0 = u
10+
A = a+A0
11+
val = 3*eq_parab.nu * (-(dudx[1] + dudx[4])*Q/A + dudx[2])
12+
return SVector(0.0,val,0,0)
13+
end
14+
15+
function source_term_simple(u, x, t, eq::BloodFlowEquations1DOrd2)
16+
res = source_term_simple(u, x, t, eq.model1d)
17+
k = friction(u,x,eq.model1d)
18+
R = radius(u,eq.model1d)
19+
return SVector(res[1],res[2]/(1-R*k/(4*eq.nu)),res[3],res[4])
20+
end
21+
22+
@inline function boundary_condition_pressure_in(flux_inner, u_inner,
23+
orientation_or_normal,direction,
24+
x, t,
25+
operator_type::Trixi.Gradient,
26+
equations_parabolic::BloodFlowEquations1DOrd2)
27+
return boundary_condition_pressure_in(u_inner,orientation_or_normal,direction,x,t,flux_lax_friedrichs,equations_parabolic.model1d)
28+
end
29+
@inline function boundary_condition_pressure_in(flux_inner, u_inner,
30+
orientation_or_normal,direction,
31+
x, t,
32+
operator_type::Trixi.Divergence,
33+
equations_parabolic::BloodFlowEquations1DOrd2)
34+
return flux_inner
35+
end
36+
# Dirichlet and Neumann boundary conditions for use with parabolic solvers in weak form.
37+
# Note that these are general, so they apply to LaplaceDiffusion in any spatial dimension.
38+
@inline function (boundary_condition::Trixi.BoundaryConditionDirichlet)(flux_inner, u_inner,
39+
normal::AbstractVector,
40+
x, t,
41+
operator_type::Trixi.Gradient,
42+
equations_parabolic::BloodFlowEquations1DOrd2)
43+
return boundary_condition.boundary_value_function(x, t, equations_parabolic)
44+
end
45+
46+
@inline function (boundary_condition::Trixi.BoundaryConditionDirichlet)(flux_inner, u_inner,
47+
normal::AbstractVector,
48+
x, t,
49+
operator_type::Trixi.Divergence,
50+
equations_parabolic::BloodFlowEquations1DOrd2)
51+
return flux_inner
52+
end
53+
54+
@inline function (boundary_condition::Trixi.BoundaryConditionNeumann)(flux_inner, u_inner,
55+
normal::AbstractVector,
56+
x, t,
57+
operator_type::Trixi.Divergence,
58+
equations_parabolic::BloodFlowEquations1DOrd2)
59+
return boundary_condition.boundary_normal_flux_function(x, t, equations_parabolic)
60+
end
61+
62+
@inline function (boundary_condition::Trixi.BoundaryConditionNeumann)(flux_inner, u_inner,
63+
normal::AbstractVector,
64+
x, t,
65+
operator_type::Trixi.Gradient,
66+
equations_parabolic::BloodFlowEquations1DOrd2)
67+
return flux_inner
68+
end

0 commit comments

Comments
 (0)