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

Contents of /MITgcm/pkg/generic_advdiff/gad_pqm_hat_r.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 SUBROUTINE GAD_PQM_HAT_R(bi,bj,ix,iy,
7 & method,mask,fbar,edge,
8 & ohat,fhat,myThid)
9 C |================================================================|
10 C | PQM_HAT_R: reconstruct grid-cell PQM polynomials. |
11 C |================================================================|
12
13 implicit none
14
15 C =============================================== global variables
16 # include "SIZE.h"
17 # include "GRID.h"
18 # include "GAD.h"
19
20 C ====================================================== arguments
21 C bi,bj :: tile indexing.
22 C ix,iy :: x-,y-indexing.
23 C method :: advection scheme.
24 C mask :: col. of cell-wise mask values.
25 C fbar :: col. of cell-wise values.
26 C edge :: col. of edge-wise values/slopes.
27 C - EDGE(1,:) = VALUE.
28 C - EDGE(2,:) = DF/DR.
29 C ohat :: col. of oscl. coeff.
30 C - OHAT(1,:) = D^1F/DS^1.
31 C - OHAT(2,:) = D^2F/DS^2.
32 C fhat :: col. of poly. coeff.
33 C - FHAT(:,I) = PQM coeff.
34 C myThid :: thread number.
35 C ================================================================
36 integer bi,bj
37 integer ix,iy
38 integer method
39 _RL mask(1-3:Nr+3)
40 _RL fbar(1-3:Nr+3)
41 _RL edge(1:2,
42 & 1-0:Nr+1)
43 _RL ohat(1:2,
44 & 1-3:Nr+3)
45 _RL fhat(1:5,
46 & 1-0:Nr+0)
47 integer myThid
48
49 C ====================================================== variables
50 C ii,ir :: local r-indexing.
51 C ff00 :: centre-biased cell value.
52 C ffll,ffrr :: left-, and right-biased cell values.
53 C xhat :: local coord. scaling.
54 C fell,ferr :: left-, and right-biased edge values.
55 C dell,derr :: left-, and right-biased edge slopes.
56 C dfds :: linear slope estimates.
57 C uhat :: "NULL" limited poly. coeff.
58 C lhat :: "MONO" limited poly. coeff.
59 C scal :: "WENO" weights.
60 C fmag,fdel :: local perturbation indicators.
61 C ================================================================
62 integer ii,ir
63 _RL ff00
64 _RL ffll,ffrr
65 _RL rhat
66 _RL fell,ferr
67 _RL dell,derr
68 _RL dfds(-1:+1)
69 _RL uhat(+1:+5)
70 _RL lhat(+1:+5)
71 _RL scal(+1:+2)
72 _RL fmag,fdel
73 integer mono
74
75 do ir = +1, Nr
76
77 C =============================== scale to local grid-cell co-ords
78 rhat = drF(ir) * .5 _d 0
79
80 C =============================== assemble cell mean + edge values
81 ff00 = fbar(ir+0)
82 ffll = ff00
83 & + mask(ir-1)*(fbar(ir-1)-ff00)
84 ffrr = ff00
85 & + mask(ir+1)*(fbar(ir+1)-ff00)
86
87 fell = edge(+1,ir-0)
88 ferr = edge(+1,ir+1)
89
90 dell = edge(+2,ir-0)
91 derr = edge(+2,ir+1)
92
93 dell = dell * rhat
94 derr = derr * rhat
95
96 c select case(method)
97 c case(ENUM_PQM_NULL_LIMIT)
98 if ( method.eq.ENUM_PQM_NULL_LIMIT ) then
99 C =============================== "NULL" limited grid-cell profile
100 CALL GAD_PQM_FUN_NULL ( ff00,
101 & fell,ferr,dell,derr,lhat,mono)
102
103 c case(ENUM_PQM_MONO_LIMIT)
104 elseif ( method.eq.ENUM_PQM_MONO_LIMIT ) then
105 C =============================== "MONO" limited grid-cell profile
106 CALL GAD_PLM_FUN_U(ffll,ff00,ffrr,dfds)
107
108 CALL GAD_PQM_FUN_MONO ( ff00,ffll,ffrr,
109 & fell,ferr,dell,derr,dfds,lhat,
110 & mono)
111
112 c case(ENUM_PQM_WENO_LIMIT)
113 elseif ( method.eq.ENUM_PQM_WENO_LIMIT ) then
114 C =============================== "WENO" limited grid-cell profile
115 CALL GAD_PLM_FUN_U(ffll,ff00,ffrr,dfds)
116
117 CALL GAD_PQM_FUN_NULL ( ff00,
118 & fell,ferr,dell,derr,uhat,mono)
119
120 CALL GAD_PQM_FUN_MONO ( ff00,ffll,ffrr,
121 & fell,ferr,dell,derr,dfds,lhat,
122 & mono)
123
124 if ( mono .gt. 0) then
125
126 C =============================== only apply WENO if it is worth it
127 fdel = abs(ffrr-ff00)+abs(ff00-ffll)
128 fmag = abs(ffll)+abs(ff00)+abs(ffrr)
129
130 if (fdel .gt. 1. _d -6 * fmag) then
131
132 C =============================== calc. WENO oscillation weighting
133 CALL GAD_OSC_MUL_R(ir,+2,mask,
134 & ohat,scal)
135
136 do ii = +1, +5
137 C =============================== blend limited/un-limited profile
138 lhat(ii) = scal(1) * uhat(ii)
139 & + scal(2) * lhat(ii)
140 end do
141
142 end if
143
144 end if
145
146 c end select
147 endif
148
149 do ii = +1, +5
150 C =============================== copy polynomial onto output data
151 fhat(ii,ir) = lhat(ii)
152 end do
153
154 end do
155
156 return
157
158 c end subroutine GAD_PQM_HAT_R
159 end

  ViewVC Help
Powered by ViewVC 1.1.22