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

Annotation of /MITgcm/pkg/generic_advdiff/gad_pqm_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: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 jmc 1.1 C $Header: $
2     C $Name: $
3    
4     # include "GAD_OPTIONS.h"
5    
6     C-- File gad_pqm_fun.F: Routines to form PQM grid-cell polynomial.
7     C-- Contents
8     C-- o QUADROOT
9     C-- o GAD_PPM_FUN_NULL
10     C-- o GAD_PPM_FUN_MONO
11    
12     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
13    
14     LOGICAL FUNCTION QUADROOT(aa,bb,cc,xx)
15     C |================================================================|
16     C | QUADROOT: find roots of quadratic ax**2 + bx + c = 0. |
17     C |================================================================|
18    
19     implicit none
20    
21     C ====================================================== arguments
22     _RL aa,bb,cc
23     _RL xx(1:2)
24    
25     C ====================================================== variables
26     _RL sq,a0,b0
27    
28     a0 = abs(aa)
29     b0 = abs(bb)
30    
31     sq = bb * bb - 4. _d 0 * aa * cc
32    
33     if (a0 .gt. 0. _d 0) then
34    
35     if (sq .ge. 0. _d 0) then
36    
37     QUADROOT = .TRUE.
38    
39     sq = sqrt(sq)
40    
41     xx(1) = - bb + sq
42     xx(2) = - bb - sq
43    
44     aa = 0.5 _d 0 / aa
45    
46     xx(1) = xx(1) * aa
47     xx(2) = xx(2) * aa
48    
49     else
50    
51     QUADROOT = .FALSE.
52    
53     end if
54    
55     else
56    
57     if (b0 .gt. 0. _d 0) then
58    
59     QUADROOT = .TRUE.
60    
61     xx(1) = - cc / bb
62     xx(2) = - cc / bb
63    
64     else
65    
66     QUADROOT = .FALSE.
67    
68     end if
69    
70     end if
71    
72     return
73    
74     c end function QUADROOT
75     end
76    
77     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
78    
79     SUBROUTINE GAD_PQM_FUN_NULL(
80     I ff00,
81     I fell,ferr,
82     I dell,derr,
83     O fhat,mono)
84     C |================================================================|
85     C | PQM_FUN_NULL: form PQM grid-cell polynomial. |
86     C | Piecewise Quartic Method (PQM) - unlimited variant. |
87     C |================================================================|
88    
89     implicit none
90    
91     C ====================================================== arguments
92     _RL ff00
93     _RL fell,ferr
94     _RL dell,derr
95     _RL fhat(+1:+5)
96     integer mono
97    
98     mono = +0
99    
100     C ============================================== unlimited profile
101     fhat(1) =
102     & + (30. _d 0 / 16. _d 0) * ff00
103     & - ( 7. _d 0 / 16. _d 0) *(ferr+fell)
104     & + ( 1. _d 0 / 16. _d 0) *(derr-dell)
105     fhat(2) =
106     & + ( 3. _d 0 / 4. _d 0) *(ferr-fell)
107     & - ( 1. _d 0 / 4. _d 0) *(derr+dell)
108     fhat(3) =
109     & - (30. _d 0 / 8. _d 0) * ff00
110     & + (15. _d 0 / 8. _d 0) *(ferr+fell)
111     & - ( 3. _d 0 / 8. _d 0) *(derr-dell)
112     fhat(4) =
113     & - ( 1. _d 0 / 4. _d 0) *(ferr-fell-derr-dell)
114     fhat(5) =
115     & + (30. _d 0 / 16. _d 0) * ff00
116     & - (15. _d 0 / 16. _d 0) *(ferr+fell)
117     & + ( 5. _d 0 / 16. _d 0) *(derr-dell)
118    
119     return
120    
121     c end subroutine GAD_PQM_FUN_NULL
122     end
123    
124     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
125    
126     SUBROUTINE GAD_PQM_FUN_MONO(
127     I ff00,
128     I ffll,ffrr,
129     I fell,ferr,
130     I dell,derr,
131     I dfds,
132     O fhat,mono)
133     C |================================================================|
134     C | PQM_FUN_MONO: form PQM grid-cell polynomial. |
135     C | Piecewise Quartic Method (PQM) - monotonic variant. |
136     C |================================================================|
137    
138     implicit none
139    
140     C ====================================================== arguments
141     _RL ff00,ffll,ffrr
142     _RL fell,ferr
143     _RL dell,derr
144     _RL dfds(-1:+1)
145     _RL fhat(+1:+5)
146     integer mono
147    
148     C ====================================================== functions
149     logical QUADROOT
150    
151     C ====================================================== variables
152     integer bind
153     _RL aval,bval,cval
154     _RL iflx(+1:+2)
155     _RL dflx(+1:+2)
156    
157     mono = +0
158    
159     C ============================================== "flatten" extrema
160     if((ffrr-ff00) *
161     & (ff00-ffll) .le. 0. _d 0) then
162    
163     mono = +1
164    
165     fhat(1) = ff00
166    
167     fhat(2) = 0. _d 0
168     fhat(3) = 0. _d 0
169     fhat(4) = 0. _d 0
170     fhat(5) = 0. _d 0
171    
172     return
173    
174     end if
175    
176     C ============================================== limit edge values
177     if((ffll - fell) *
178     & (fell - ff00) .le. 0. _d 0) then
179    
180     mono = +1
181    
182     fell = ff00 - dfds(0)
183    
184     end if
185    
186     if((ffrr - ferr) *
187     & (ferr - ff00) .le. 0. _d 0) then
188    
189     mono = +1
190    
191     ferr = ff00 + dfds(0)
192    
193     end if
194    
195     C ============================================== limit edge slopes
196     if((dell * dfds(-1)) .lt. 0. _d 0) then
197    
198     mono = +1
199    
200     dell = dfds(-1)
201    
202     end if
203    
204     if((derr * dfds(+1)) .lt. 0. _d 0) then
205    
206     mono = +1
207    
208     derr = dfds(+1)
209    
210     end if
211    
212     C ============================================== limit cell values
213     fhat(1) =
214     & + (30. _d 0 / 16. _d 0) * ff00
215     & - ( 7. _d 0 / 16. _d 0) *(ferr+fell)
216     & + ( 1. _d 0 / 16. _d 0) *(derr-dell)
217     fhat(2) =
218     & + ( 3. _d 0 / 4. _d 0) *(ferr-fell)
219     & - ( 1. _d 0 / 4. _d 0) *(derr+dell)
220     fhat(3) =
221     & - (30. _d 0 / 8. _d 0) * ff00
222     & + (15. _d 0 / 8. _d 0) *(ferr+fell)
223     & - ( 3. _d 0 / 8. _d 0) *(derr-dell)
224     fhat(4) =
225     & - ( 1. _d 0 / 4. _d 0) *(ferr-fell-derr-dell)
226     fhat(5) =
227     & + (30. _d 0 / 16. _d 0) * ff00
228     & - (15. _d 0 / 16. _d 0) *(ferr+fell)
229     & + ( 5. _d 0 / 16. _d 0) *(derr-dell)
230    
231     C ============================= calc. inflexion via 2nd-derivative
232     aval = 12. _d 0 * fhat(5)
233     bval = 6. _d 0 * fhat(4)
234     cval = 2. _d 0 * fhat(3)
235    
236     if ( QUADROOT(aval,bval,cval,iflx) ) then
237    
238     bind = +0
239    
240     if ( ( iflx(1) .gt. -1. _d 0 )
241     & .and. ( iflx(1) .lt. +1. _d 0 ) ) then
242    
243     C ============================= check for non-monotonic inflection
244     dflx(1) = fhat(2)
245     & + iflx(1) * fhat(3) * 2. _d 0
246     & +(iflx(1) ** 2) * fhat(4) * 3. _d 0
247     & +(iflx(1) ** 3) * fhat(5) * 4. _d 0
248    
249     if (dflx(1)*dfds(+0) .lt. 0. _d 0) then
250    
251     if (abs(dell)
252     & .lt. abs(derr) ) then
253    
254     bind = -1
255    
256     else
257    
258     bind = +1
259    
260     end if
261    
262     end if
263    
264     end if
265    
266     if ( ( iflx(2) .gt. -1. _d 0 )
267     & .and. ( iflx(2) .lt. +1. _d 0 ) ) then
268    
269     C ============================= check for non-monotonic inflection
270     dflx(2) = fhat(2)
271     & + iflx(2) * fhat(3) * 2. _d 0
272     & +(iflx(2) ** 2) * fhat(4) * 3. _d 0
273     & +(iflx(2) ** 3) * fhat(5) * 4. _d 0
274    
275     if (dflx(2)*dfds(+0) .lt. 0. _d 0) then
276    
277     if (abs(dell)
278     & .lt. abs(derr) ) then
279    
280     bind = -1
281    
282     else
283    
284     bind = +1
285    
286     end if
287    
288     end if
289    
290     end if
291    
292     C ============================= pop non-monotone inflexion to edge
293    
294     if (bind .eq. -1) then
295    
296     C ============================= pop inflection points onto -1 edge
297     mono = +2
298    
299     derr =
300     & -( 5. _d 0 / 1. _d 0) * ff00
301     & +( 3. _d 0 / 1. _d 0) * ferr
302     & +( 2. _d 0 / 1. _d 0) * fell
303     dell =
304     & +( 5. _d 0 / 3. _d 0) * ff00
305     & -( 1. _d 0 / 3. _d 0) * ferr
306     & -( 4. _d 0 / 3. _d 0) * fell
307    
308     if (dell*dfds(-1) .lt. 0. _d 0) then
309    
310     dell = 0. _d 0
311    
312     ferr =
313     & +( 5. _d 0 / 1. _d 0) * ff00
314     & -( 4. _d 0 / 1. _d 0) * fell
315     derr =
316     & +(10. _d 0 / 1. _d 0) * ff00
317     & -(10. _d 0 / 1. _d 0) * fell
318    
319     end if
320    
321     if (derr*dfds(+1) .lt. 0. _d 0) then
322    
323     derr = 0. _d 0
324    
325     fell =
326     & +( 5. _d 0 / 2. _d 0) * ff00
327     & -( 3. _d 0 / 2. _d 0) * ferr
328     dell =
329     & -( 5. _d 0 / 3. _d 0) * ff00
330     & +( 5. _d 0 / 3. _d 0) * ferr
331    
332     end if
333    
334     end if
335    
336     if (bind .eq. +1) then
337    
338     C ============================= pop inflection points onto +1 edge
339     mono = +2
340    
341     derr =
342     & -( 5. _d 0 / 3. _d 0) * ff00
343     & +( 4. _d 0 / 3. _d 0) * ferr
344     & +( 1. _d 0 / 3. _d 0) * fell
345     dell =
346     & +( 5. _d 0 / 1. _d 0) * ff00
347     & -( 2. _d 0 / 1. _d 0) * ferr
348     & -( 3. _d 0 / 1. _d 0) * fell
349    
350     if (dell*dfds(-1) .lt. 0. _d 0) then
351    
352     dell = 0. _d 0
353    
354     ferr =
355     & +( 5. _d 0 / 2. _d 0) * ff00
356     & -( 3. _d 0 / 2. _d 0) * fell
357     derr =
358     & +( 5. _d 0 / 3. _d 0) * ff00
359     & -( 5. _d 0 / 3. _d 0) * fell
360    
361     end if
362    
363     if (derr*dfds(+1) .lt. 0. _d 0) then
364    
365     derr = 0. _d 0
366    
367     fell =
368     & +( 5. _d 0 / 1. _d 0) * ff00
369     & -( 4. _d 0 / 1. _d 0) * ferr
370     dell =
371     & -(10. _d 0 / 1. _d 0) * ff00
372     & +(10. _d 0 / 1. _d 0) * ferr
373    
374     end if
375    
376     end if
377    
378     end if
379    
380     C ============================= re-assemble coefficients on demand
381     if (mono .eq. +2) then
382    
383     fhat(1) =
384     & + (30. _d 0 / 16. _d 0) * ff00
385     & - ( 7. _d 0 / 16. _d 0) *(ferr+fell)
386     & + ( 1. _d 0 / 16. _d 0) *(derr-dell)
387     fhat(2) =
388     & + ( 3. _d 0 / 4. _d 0) *(ferr-fell)
389     & - ( 1. _d 0 / 4. _d 0) *(derr+dell)
390     fhat(3) =
391     & - (30. _d 0 / 8. _d 0) * ff00
392     & + (15. _d 0 / 8. _d 0) *(ferr+fell)
393     & - ( 3. _d 0 / 8. _d 0) *(derr-dell)
394     fhat(4) =
395     & - ( 1. _d 0 / 4. _d 0) *(ferr-fell-derr-dell)
396     fhat(5) =
397     & + (30. _d 0 / 16. _d 0) * ff00
398     & - (15. _d 0 / 16. _d 0) *(ferr+fell)
399     & + ( 5. _d 0 / 16. _d 0) *(derr-dell)
400    
401     end if
402    
403     return
404    
405     c end subroutine GAD_PQM_FUN_MONO
406     end

  ViewVC Help
Powered by ViewVC 1.1.22