aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_SlopeLimiter.F90
blob: ffdfc3ad6319428e9d0f0352e053c0297254dc53 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
 /*@@
   @file      GRHydro_SlopeLimiter.F90
   @date      Sat Jan 26 02:00:32 2002
   @author    
   @desc 
   The routine for the more complex slope limiters. 
   See Toros book for most of these.
   @enddesc 
 @@*/

#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"
#include "cctk_Functions.h"

 /*@@
   @routine    slopelimiter
   @date       Sat Jan 26 02:01:15 2002
   @author     Luca Baiotti, Ian Hawke
   @desc 
   Given the slope delta and the ratio of the two local slopes,
   limit the slope so the reconstruction is TVD.
   @enddesc 
   @calls     
   @calledby   
   @history 
   Original minmod from GR3D, author Mark Miller. Other limiters from 
   codes of Hawke, Nikiforakis and Toro.
   @endhistory 

@@*/

subroutine slopelimiter(r, delta)
  
  USE GRHydro_Scalars

  implicit none
  
  DECLARE_CCTK_FUNCTIONS

  CCTK_REAL :: phi, r, denor, phir, delta
  
  denor = 1.d0 + r
  phir = 2.d0/denor
  
  if (MC1) then
    
    if (r > 0.d0) then 
      phi = min( 2.d0*r/(1.d0+r), phir )
    else
      phi = 0.d0
    end if
    
!!$  else if (MC2) then
!!$    
!!$    if (r > 3.d0) then 
!!$      phi = 2*phir
!!$    else if (r > 1.d0/3.d0) then
!!$      phi = 1
!!$    else if (r > 0) then
!!$      phi = 2*r*phir
!!$    else
!!$      phi = 0.d0
!!$    end if
    
  else if (MINMOD2) then
    
    if (r > 1.d0) then 
      phi = min(1.d0, phir )
    else if (r > 0.d0) then
      phi = r
    else
      phi = 0.d0
    end if
    
  else if (MINMOD3) then
    
    if (r > 0.d0) then 
      phi = min(1.d0, 4.d0/(1.d0+r) )
    else
      phi = 0.d0
    end if
    
  else if (SUPERBEE) then
    
    if (r > 1.d0) then
      phi = min(2.d0, r, 2.d0 / (1.d0 + r))
    else if (r > 0.5d0) then
      phi = 1.d0
    else if (r > 0.d0) then
      phi = 2.d0 * r
    else
      phi = 0.d0
    end if
    
  else 
    
    call CCTK_WARN(0, "Type of limiter not recognized")
    
  end if

  delta = delta * phi

end subroutine slopelimiter