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

Annotation of /MITgcm/pkg/generic_advdiff/gad_osc_hat_r.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_osc_hat_r.F: Routines ???
7     C-- Contents
8     C-- o GAD_OSC_LOC_R
9     C-- o GAD_OSC_HAT_R
10    
11     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
12    
13     SUBROUTINE GAD_OSC_LOC_R(ir,mask,fbar,ohat)
14    
15     implicit none
16    
17     C =============================================== global variables
18     # include "SIZE.h"
19    
20     C ====================================================== arguments
21     integer ir
22     _RL mask(1-3:Nr+3)
23     _RL fbar(1-3:Nr+3)
24     _RL ohat(1:2,
25     & 1-3:Nr+3)
26    
27     C ====================================================== variables
28     _RL floc(-1:+1)
29    
30     if (ir.gt.+1-3 .and.
31     & ir.lt.Nr+3) then
32    
33     C ================ mask local stencil: expand from centre outwards
34    
35     floc(+0) = fbar(+0+ir)
36    
37     floc(-1) = floc(+0) +
38     & mask(ir-1)*(fbar(ir-1)-floc(+0))
39     floc(+1) = floc(+0) +
40     & mask(ir+1)*(fbar(ir+1)-floc(+0))
41    
42     C ================ calc. 1st & 2nd derivatives over masked stencil
43    
44     ohat(+1,ir) = floc(+1)*0.25 _d 0
45     & - floc(-1)*0.25 _d 0
46    
47     ohat(+2,ir) = floc(+1)*0.25 _d 0
48     & - floc(+0)*0.50 _d 0
49     & + floc(-1)*0.25 _d 0
50    
51     else
52    
53     if (ir.eq.+1-3) then
54    
55     C ================ mask local stencil: expand from centre outwards
56    
57     floc(+0) = fbar(+0+ir)
58    
59     floc(+1) = floc(+0) +
60     & mask(ir+1)*(fbar(ir+1)-floc(+0))
61    
62     C ================ calc. 1st & 2nd derivatives over masked stencil
63    
64     ohat(+1,ir) = floc(+1)*0.50 _d 0
65     & - floc(+0)*0.50 _d 0
66    
67     ohat(+2,ir) = 0. _d 0
68    
69     end if
70    
71     if (ir.eq.Nr+3) then
72    
73     C ================ mask local stencil: expand from centre outwards
74    
75     floc(+0) = fbar(+0+ir)
76    
77     floc(-1) = floc(+0) +
78     & mask(ir-1)*(fbar(ir-1)-floc(+0))
79    
80     C ================ calc. 1st & 2nd derivatives over masked stencil
81    
82     ohat(+1,ir) = floc(+0)*0.50 _d 0
83     & - floc(-1)*0.50 _d 0
84    
85     ohat(+2,ir) = 0. _d 0
86    
87     end if
88    
89     end if
90    
91     return
92    
93     c end subroutine GAD_OSC_LOC_R
94     end
95    
96     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
97    
98     SUBROUTINE GAD_OSC_HAT_R(bi,bj,ix,iy,
99     & mask,fbar,ohat,myThid)
100     C |================================================================|
101     C | OSC_HAT_R: compute WENO oscillation derivatives in R. |
102     C |================================================================|
103    
104     implicit none
105    
106     C =============================================== global variables
107     # include "SIZE.h"
108    
109     C ====================================================== arguments
110     integer bi,bj,ix,iy
111     _RL mask(1-3:Nr+3)
112     _RL fbar(1-3:Nr+3)
113     _RL ohat(1:2,
114     & 1-3:Nr+3)
115     integer myThid
116    
117     C ====================================================== variables
118     integer ir
119    
120     C ================================ derivatives for WENO indicators
121     do ir = +1-3, Nr+3
122    
123     CALL GAD_OSC_LOC_R(ir,mask,fbar,ohat)
124    
125     end do
126    
127     return
128    
129     c end subroutine GAD_OSC_HAT_R
130     end

  ViewVC Help
Powered by ViewVC 1.1.22