MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_rho_roe.t
Go to the documentation of this file.
1!> Module containing Roe solver for scalar advection
5
6 implicit none
7 private
8
9 public :: rho_roe_init
10
11contains
12
13 subroutine rho_roe_init()
15
16 nworkroe = 1
17
18 phys_average => rho_average
19 phys_get_eigenjump => rho_get_eigenjump
20 phys_rtimes => rho_rtimes
21 end subroutine rho_roe_init
22
23 subroutine rho_average(wL, wR, x, ix^L, idim, wroe, workroe)
25 integer, intent(in) :: ix^l, idim
26 double precision, intent(in) :: wl(ixg^t, nw), wr(ixg^t, nw)
27 double precision, intent(inout) :: wroe(ixg^t, nw)
28 double precision, intent(inout) :: workroe(ixg^t, nworkroe)
29 double precision, intent(in) :: x(ixg^t, 1:^nd)
30
31 wroe(ix^s, rho_)=half*(wl(ix^s, rho_)+wr(ix^s, rho_))
32 end subroutine rho_average
33
34 subroutine rho_get_eigenjump(wL, wR, wC, x, ix^L, il, idim, smalla, a, jump, workroe)
35
36 ! Calculate the characteristic speed a and the jump in the
37 ! characteristic variable in the idim direction within ixL.
38 ! For a scalar equation the characteristic and conservative variables coincide
39 ! The characteristic speed is just the velocity, but it should be averaged
40 ! for the cell interfaces
41
43
44 integer, intent(in) :: ix^l, il, idim
45 double precision, dimension(ixG^T, nw) :: wl, wr, wc
46 double precision, dimension(ixG^T) :: smalla, a, jump, v
47 double precision, dimension(ixG^T, nworkroe) :: workroe
48 double precision, intent(in) :: x(ixg^t, 1:^nd)
49 integer :: jx^l, ixc^l
50
51 jx^l=ix^l+kr(idim,^d);
52 ixcmin^d=ixmin^d; ixcmax^d=jxmax^d;
53
54 ! No entropy fix
55 smalla(ix^s)= -one
56 ! The velocity is independent of w in the transport equation,
57 ! but it may depend on the location
58 call rho_get_v(wl, x, ixg^ll, ixc^l, idim, v, .true.)
59
60 a(ix^s)=(v(jx^s)+v(ix^s))/2
61
62 jump(ix^s)=wr(ix^s, rho_)-wl(ix^s, rho_)
63
64 end subroutine rho_get_eigenjump
65
66 subroutine rho_rtimes(q, w, ix^L, iw, il, idim, rq, workroe)
67
68 ! Multiply q by R(il, iw), where R is the right eigenvalue matrix at wC.
69 ! For a scalar equation the R matrix is unity
70
72 integer, intent(in) :: ix^l, iw, il, idim
73 double precision, intent(in) :: w(ixg^t, nw), q(ixg^t)
74 double precision, intent(inout) :: rq(ixg^t)
75 double precision, intent(inout) :: workroe(ixg^t, nworkroe)
76
77 rq(ix^s)=q(ix^s)
78
79 end subroutine rho_rtimes
80
81end module mod_rho_roe
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, dimension(3, 3) kr
Kronecker delta tensor.
double precision, dimension(:), allocatable, parameter d
procedure(sub_rtimes), pointer phys_rtimes
procedure(sub_get_eigenjump), pointer phys_get_eigenjump
procedure(sub_average), pointer phys_average
Module containing the physics routines for scalar advection.
Definition mod_rho_phys.t:2
subroutine, public rho_get_v(w, x, ixil, ixol, idim, v, centered)
integer, public, protected rho_
Definition mod_rho_phys.t:7
Module containing Roe solver for scalar advection.
Definition mod_rho_roe.t:2
subroutine, public rho_roe_init()
Definition mod_rho_roe.t:14