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

Contents of /MITgcm/pkg/generic_advdiff/gad_ppm_fun.F

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


Revision 1.1 - (show annotations) (download)
Sun Mar 13 01:44:03 2016 UTC (8 years, 1 month 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 C $Header: $
2 C $Name: $
3
4 # include "GAD_OPTIONS.h"
5
6 C-- File gad_ppm_fun.F: Routines to form PPM grid-cell polynomial.
7 C-- Contents
8 C-- o GAD_PPM_FUN_NULL
9 C-- o GAD_PPM_FUN_MONO
10
11 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
12
13 SUBROUTINE GAD_PPM_FUN_NULL(
14 I ff00,
15 I fell,ferr,
16 O fhat,mono)
17 C |================================================================|
18 C | PPM_FUN_NULL: form PPM grid-cell polynomial. |
19 C | Piecewise Parabolic Method (PPM), unlimited variant. |
20 C |================================================================|
21
22 implicit none
23
24 C ====================================================== arguments
25 _RL ff00
26 _RL fell,ferr
27 _RL fhat(+1:+3)
28 integer mono
29
30 mono = +0
31
32 C ============================================== unlimited profile
33 fhat( 1 ) =
34 & +(3. _d 0 / 2. _d 0) * ff00
35 & -(1. _d 0 / 4. _d 0) *(ferr+fell)
36 fhat( 2 ) =
37 & +(1. _d 0 / 2. _d 0) *(ferr-fell)
38 fhat( 3 ) =
39 & -(3. _d 0 / 2. _d 0) * ff00
40 & +(3. _d 0 / 4. _d 0) *(ferr+fell)
41
42 return
43
44 c end subroutine GAD_PPM_FUN_NULL
45 end
46
47 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
48
49 SUBROUTINE GAD_PPM_FUN_MONO(
50 I ff00,
51 I ffll,ffrr,
52 I fell,ferr,
53 I dfds,
54 O fhat,mono)
55 C |================================================================|
56 C | PPM_FUN_MONO: form PPM grid-cell polynomial. |
57 C | Piecewise Parabolic Method (PPM) - monotonic variant. |
58 C |================================================================|
59
60 implicit none
61
62 C ====================================================== arguments
63 _RL ff00
64 _RL ffll,ffrr
65 _RL fell,ferr
66 _RL dfds(-1:+1)
67 _RL fhat(+1:+3)
68 integer mono
69
70 C ====================================================== variables
71 _RL turn
72
73 mono = 0
74
75 C ============================================== "flatten" extrema
76 if((ffrr-ff00) *
77 & (ff00-ffll) .le. 0. _d 0) then
78
79 mono = +1
80
81 fhat(1) = ff00
82
83 fhat(2) = 0. _d 0
84 fhat(3) = 0. _d 0
85
86 return
87
88 end if
89
90 C ============================================== limit edge values
91 if((ffll-fell) *
92 & (fell-ff00) .le. 0. _d 0) then
93
94 mono = +1
95
96 fell = ff00 - dfds(0)
97
98 end if
99
100 if((ffrr-ferr) *
101 & (ferr-ff00) .le. 0. _d 0) then
102
103 mono = +1
104
105 ferr = ff00 + dfds(0)
106
107 end if
108
109 C ============================================== limit cell values
110 fhat( 1 ) =
111 & +(3. _d 0 / 2. _d 0) * ff00
112 & -(1. _d 0 / 4. _d 0) *(ferr+fell)
113 fhat( 2 ) =
114 & +(1. _d 0 / 2. _d 0) *(ferr-fell)
115 fhat( 3 ) =
116 & -(3. _d 0 / 2. _d 0) * ff00
117 & +(3. _d 0 / 4. _d 0) *(ferr+fell)
118
119 if (abs(fhat(3)) .gt.
120 & abs(fhat(2))*.5 _d 0) then
121
122 turn = -0.5 _d 0 * fhat(2)
123 & / fhat(3)
124
125 if ((turn .ge. -1. _d 0)
126 & .and.(turn .le. +0. _d 0)) then
127
128 mono = +2
129
130 C ====================================== push TURN onto lower edge
131 ferr = +3. _d 0 * ff00
132 & -2. _d 0 * fell
133
134 end if
135
136 if ((turn .gt. +0. _d 0)
137 & .and.(turn .le. +1. _d 0)) then
138
139 mono = +2
140
141 C ====================================== push TURN onto upper edge
142 fell = +3. _d 0 * ff00
143 & -2. _d 0 * ferr
144
145 end if
146
147 end if
148
149 if (mono .gt. +1) then
150
151 C ====================================== re-calc. coeff. on demand
152 fhat( 1 ) =
153 & +(3. _d 0 / 2. _d 0) * ff00
154 & -(1. _d 0 / 4. _d 0) *(ferr+fell)
155 fhat( 2 ) =
156 & +(1. _d 0 / 2. _d 0) *(ferr-fell)
157 fhat( 3 ) =
158 & -(3. _d 0 / 2. _d 0) * ff00
159 & +(3. _d 0 / 4. _d 0) *(ferr+fell)
160
161 end if
162
163 return
164
165 c end subroutine GAD_PPM_FUN_MONO
166 end

  ViewVC Help
Powered by ViewVC 1.1.22