/[MITgcm]/MITgcm/pkg/generic_advdiff/gad_plm_fun.F
ViewVC logotype

Annotation of /MITgcm/pkg/generic_advdiff/gad_plm_fun.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.1 - (hide annotations) (download)
Sun Mar 13 01:44:02 2016 UTC (8 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65v, checkpoint65w, checkpoint65u, HEAD
- from Darren: add PPM and PQM advection schemes (number 40-42 and 50-52)
  with 2 types of limiter (see: Engwirda & Kelley, submit. to JCP);
  Note (from Darren): unlimited PPM/PQM scheme (40 & 50) are just for
  testing and not for actual use.

1 jmc 1.1 C $Header: $
2     C $Name: $
3    
4     # include "GAD_OPTIONS.h"
5    
6     C-- File gad_plm_fun.F: Routines for monotone piecewise linear method.
7     C-- Contents
8     C-- o GAD_PLM_FUN_U
9     C-- o GAD_PLM_FUN_V
10    
11     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
12    
13     SUBROUTINE GAD_PLM_FUN_U(
14     I ffll,ff00,ffrr,
15     O dfds)
16     C |================================================================|
17     C | PLM_FUN_U: monotone piecewise linear method. |
18     C | - uniform grid-spacing variant. |
19     C |================================================================|
20    
21     implicit none
22    
23     C ====================================================== arguments
24     _RL ffll,ff00,ffrr
25     _RL dfds(-1:+1)
26    
27     C ====================================================== variables
28     _RL fell,ferr,scal
29     _RL epsil
30     PARAMETER( epsil = 1. _d -16 )
31    
32     dfds(-1) = ff00 - ffll
33     dfds(+1) = ffrr - ff00
34    
35     if (dfds(-1) * dfds(+1) .gt. 0.0 _d 0) then
36    
37     C ======================================= calc. ll//rr edge values
38     fell = 0.5 _d 0 * (ffll + ff00)
39     ferr = 0.5 _d 0 * (ff00 + ffrr)
40    
41     C ======================================= calc. centred derivative
42     dfds(+0) =
43     & 0.5 _d 0 * (ferr - fell)
44    
45     C ======================================= monotonic slope-limiting
46     scal = min(abs(dfds(-1)),
47     & abs(dfds(+1)))
48     & / max(abs(dfds(+0)), epsil)
49     c & / max(abs(dfds(+0)), epsilon(ff00))
50     scal = min(scal, 1.0 _d 0)
51    
52     dfds(+0) = scal * dfds(+0)
53    
54     else
55    
56     C ======================================= flatten if local extrema
57     dfds(+0) = 0.0 _d 0
58    
59     end if
60    
61     dfds(-1) = 0.5 _d 0 * dfds(-1)
62     dfds(+1) = 0.5 _d 0 * dfds(+1)
63    
64     return
65    
66     c end subroutine GAD_PLM_FUN_U
67     end
68    
69     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
70    
71     SUBROUTINE GAD_PLM_FUN_V(
72     I ffll,ff00,ffrr,
73     I ddll,dd00,ddrr,
74     O dfds)
75     C |================================================================|
76     C | PLM_FUN_V: monotone piecewise linear method. |
77     C | - variable grid-spacing variant. |
78     C |================================================================|
79    
80     implicit none
81    
82     C ====================================================== arguments
83     _RL ffll,ff00,ffrr
84     _RL ddll,dd00,ddrr
85     _RL dfds(-1:+1)
86    
87     C ====================================================== variables
88     _RL fell,ferr,scal
89     _RL epsil
90     PARAMETER( epsil = 1. _d -16 )
91    
92     dfds(-1) = ff00 - ffll
93     dfds(+1) = ffrr - ff00
94    
95     if (dfds(-1) * dfds(+1) .gt. 0.0 _d 0) then
96    
97     C ======================================= calc. ll//rr edge values
98     fell = (dd00 * ffll + ddll * ff00)
99     & / (ddll + dd00)
100     ferr = (ddrr * ff00 + dd00 * ffrr)
101     & / (dd00 + ddrr)
102    
103     C ======================================= calc. centred derivative
104     dfds(+0) =
105     & 0.5 _d 0 * (ferr - fell)
106    
107     C ======================================= monotonic slope-limiting
108     scal = min(abs(dfds(-1)),
109     & abs(dfds(+1)))
110     & / max(abs(dfds(+0)), epsil)
111     c & / max(abs(dfds(+0)), epsilon(ff00))
112     scal = min(scal, 1.0 _d 0)
113    
114     dfds(+0) = scal * dfds(+0)
115    
116     else
117    
118     C ======================================= flatten if local extrema
119     dfds(+0) = 0.0 _d 0
120    
121     end if
122    
123     C == !! check this
124     dfds(-1) = dfds(-1) / (ddll + dd00) * dd00
125     dfds(+1) = dfds(+1) / (dd00 + ddrr) * dd00
126    
127     return
128    
129     c end subroutine GAD_PLM_FUN_V
130     end

  ViewVC Help
Powered by ViewVC 1.1.22