MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_rho_roe.t
Go to the documentation of this file.
1 !> Module containing Roe solver for scalar advection
2 module mod_rho_roe
4  use mod_rho_phys
5 
6  implicit none
7  private
8 
9  public :: rho_roe_init
10 
11 contains
12 
13  subroutine rho_roe_init()
14  use mod_physics_roe
15 
16  nworkroe = 1
17 
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 
81 end 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.
integer, 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)
Definition: mod_rho_phys.t:123
integer, public, protected rho_
Definition: mod_rho_phys.t:7
Module containing Roe solver for scalar advection.
Definition: mod_rho_roe.t:2
subroutine rho_get_eigenjump(wL, wR, wC, x, ixL, il, idim, smalla, a, jump, workroe)
Definition: mod_rho_roe.t:35
subroutine rho_rtimes(q, w, ixL, iw, il, idim, rq, workroe)
Definition: mod_rho_roe.t:67
subroutine rho_average(wL, wR, x, ixL, idim, wroe, workroe)
Definition: mod_rho_roe.t:24
subroutine, public rho_roe_init()
Definition: mod_rho_roe.t:14