MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_physics_roe.t
Go to the documentation of this file.
2 use mod_comm_lib, only: mpistop
3
4 implicit none
5 public
6
7 procedure(sub_average), pointer :: phys_average => null()
8 procedure(sub_get_eigenjump), pointer :: phys_get_eigenjump => null()
9 procedure(sub_rtimes), pointer :: phys_rtimes => null()
10
11 integer :: nworkroe = -1
12
13 abstract interface
14 subroutine sub_average(wL, wR, x, ix^L, idim, wroe, workroe)
16 import
17 integer, intent(in) :: ix^L, idim
18 double precision, intent(in) :: wL(ixG^T, nw), wR(ixG^T, nw)
19 double precision, intent(inout) :: wroe(ixG^T, nw)
20 double precision, intent(inout) :: workroe(ixG^T, nworkroe)
21 double precision, intent(in) :: x(ixG^T, 1:^ND)
22 end subroutine sub_average
23
24 subroutine sub_get_eigenjump(wL, wR, wC, x, ix^L, il, &
25 idim, smalla, a, jump, workroe)
27 import
28 integer, intent(in) :: ix^L, il, idim
29 double precision, dimension(ixG^T, nw) :: wL, wR, wC
30 double precision, dimension(ixG^T) :: smalla, a, jump
31 double precision, dimension(ixG^T, nworkroe) :: workroe
32 double precision, intent(in) :: x(ixG^T, 1:^ND)
33 end subroutine sub_get_eigenjump
34
35 subroutine sub_rtimes(q, w, ix^L, iw, il, idim, rq, workroe)
37 import
38 integer, intent(in) :: ix^L, iw, il, idim
39 double precision, intent(in) :: w(ixG^T, nw), q(ixG^T)
40 double precision, intent(inout) :: rq(ixG^T)
41 double precision, intent(inout) :: workroe(ixG^T, nworkroe)
42 end subroutine sub_rtimes
43 end interface
44
45contains
46
47 subroutine phys_roe_check()
48 if (.not. associated(phys_average)) &
50
51 if (.not. associated(phys_get_eigenjump)) &
53
54 if (.not. associated(phys_rtimes)) &
56
57 end subroutine phys_roe_check
58
59 subroutine dummy_roe_average(wL, wR, x, ix^L, idim, wroe, workroe)
61 integer, intent(in) :: ix^L, idim
62 double precision, intent(in) :: wL(ixG^T, nw), wR(ixG^T, nw)
63 double precision, intent(inout) :: wroe(ixG^T, nw)
64 double precision, intent(inout) :: workroe(ixG^T, nworkroe)
65 double precision, intent(in) :: x(ixG^T, 1:^ND)
66
67 call mpistop("dummy_roe_average should not be called")
68 end subroutine dummy_roe_average
69
70 subroutine dummy_roe_get_eigenjump(wL, wR, wC, x, ix^L, il, &
71 idim, smalla, a, jump, workroe)
73 integer, intent(in) :: ix^L, il, idim
74 double precision, dimension(ixG^T, nw) :: wL, wR, wC
75 double precision, dimension(ixG^T) :: smalla, a, jump
76 double precision, dimension(ixG^T, nworkroe) :: workroe
77 double precision, intent(in) :: x(ixG^T, 1:^ND)
78
79 call mpistop("dummy_roe_get_eigenjump should not be called")
80 end subroutine dummy_roe_get_eigenjump
81
82 subroutine dummy_roe_rtimes(q, w, ix^L, iw, il, idim, rq, workroe)
84 integer, intent(in) :: ix^L, iw, il, idim
85 double precision, intent(in) :: w(ixG^T, nw), q(ixG^T)
86 double precision, intent(inout) :: rq(ixG^T)
87 double precision, intent(inout) :: workroe(ixG^T, nworkroe)
88
89 call mpistop("dummy_roe_rtimes should not be called")
90 end subroutine dummy_roe_rtimes
91
92end module mod_physics_roe
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
This module contains definitions of global parameters and variables and some generic functions/subrou...
subroutine phys_roe_check()
procedure(sub_rtimes), pointer phys_rtimes
subroutine dummy_roe_average(wl, wr, x, ixl, idim, wroe, workroe)
subroutine dummy_roe_get_eigenjump(wl, wr, wc, x, ixl, il, idim, smalla, a, jump, workroe)
subroutine dummy_roe_rtimes(q, w, ixl, iw, il, idim, rq, workroe)
procedure(sub_get_eigenjump), pointer phys_get_eigenjump
procedure(sub_average), pointer phys_average