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

Contents of /MITgcm/pkg/generic_advdiff/gad_ppm_hat_x.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_PPM_HAT_X(bi,bj,kk,iy,
7 & method,mask,fbar,edge,
8 & ohat,fhat,myThid)
9 C |================================================================|
10 C | PPM_HAT_X: reconstruct grid-cell PPM 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 kk :: r-indexing.
23 C iy :: y-indexing.
24 C method :: advection scheme.
25 C mask :: row of cell-wise masking values.
26 C fbar :: row of cell-wise values.
27 C edge :: row of edge-wise values.
28 C ohat :: row of oscl. coeff.
29 C - OHAT(1,:) = D^1F/DS^1.
30 C - OHAT(2,:) = D^2F/DS^2.
31 C fhat :: row of poly. coeff.
32 C - FHAT(:,I) = PPM coeff.
33 C myThid :: thread number.
34 C ================================================================
35 integer bi,bj,kk,iy
36 integer method
37 _RL mask(1-OLx:sNx+OLx)
38 _RL fbar(1-OLx:sNx+OLx)
39 _RL edge(1-OLx:sNx+OLx)
40 _RL ohat(1:2,
41 & 1-OLx:sNx+OLx)
42 _RL fhat(1:3,
43 & 1-OLx:sNx+OLx)
44 integer myThid
45
46 C ====================================================== variables
47 C ii,ix :: local x-indexing.
48 C ff00 :: centre-biased cell mean.
49 C ffll,ffrr :: left-, and right-biased cell values.
50 C fell,ferr :: left-, and right-biased edge values.
51 C dfds :: linear slope estimates.
52 C uhat :: "NULL" limited poly. coeff.
53 C lhat :: "MONO" limited poly. coeff.
54 C scal :: "WENO" weights.
55 C fmag,fdel :: local perturbation indicators.
56 C ================================================================
57 integer ii,ix
58 _RL ff00,ffll,ffrr
59 _RL fell,ferr
60 _RL dfds(-1:+1)
61 _RL uhat(+1:+3)
62 _RL lhat(+1:+3)
63 _RL scal(+1:+2)
64 _RL fmag,fdel
65 integer mono
66
67 do ix = 1-OLx+2, sNx+OLx-2
68
69 C =============================== assemble cell mean + edge values
70 ff00 = fbar(ix+0)
71 ffll = ff00
72 & + mask(ix-1)*(fbar(ix-1)-ff00)
73 ffrr = ff00
74 & + mask(ix+1)*(fbar(ix+1)-ff00)
75
76 fell = edge(ix-0)
77 ferr = edge(ix+1)
78
79 c select case(method)
80 c case(ENUM_PPM_NULL_LIMIT)
81 if ( method.eq.ENUM_PPM_NULL_LIMIT ) then
82 C =============================== "NULL" limited grid-cell profile
83 CALL GAD_PPM_FUN_NULL ( ff00,fell,ferr,
84 & lhat,mono)
85
86 c case(ENUM_PPM_MONO_LIMIT)
87 elseif ( method.eq.ENUM_PPM_MONO_LIMIT ) then
88 C =============================== "MONO" limited grid-cell profile
89 CALL GAD_PLM_FUN_U(ffll,ff00,ffrr,dfds)
90
91 CALL GAD_PPM_FUN_MONO ( ff00,ffll,ffrr,
92 & fell,ferr,dfds,lhat,mono)
93
94 c case(ENUM_PPM_WENO_LIMIT)
95 elseif ( method.eq.ENUM_PPM_WENO_LIMIT ) then
96 C =============================== "WENO" limited grid-cell profile
97 CALL GAD_PLM_FUN_U(ffll,ff00,ffrr,dfds)
98
99 CALL GAD_PPM_FUN_NULL ( ff00,fell,ferr,
100 & uhat,mono)
101
102 CALL GAD_PPM_FUN_MONO ( ff00,ffll,ffrr,
103 & fell,ferr,dfds,lhat,mono)
104
105 if ( mono .gt. 0) then
106
107 C =============================== only apply WENO if it is worth it
108 fdel = abs(ffrr-ff00)+abs(ff00-ffll)
109 fmag = abs(ffll)+abs(ff00)+abs(ffrr)
110
111 if (fdel .gt. 1. _d -6 * fmag) then
112
113 C =============================== calc. WENO oscillation weighting
114 CALL GAD_OSC_MUL_X(ix,+2,mask,
115 & ohat,scal)
116
117 do ii = +1, +3
118 C =============================== blend limited/un-limited profile
119 lhat(ii) = scal(1) * uhat(ii)
120 & + scal(2) * lhat(ii)
121 end do
122
123 end if
124
125 end if
126
127 c end select
128 endif
129
130 do ii = +1, +3
131 C =============================== copy polynomial onto output data
132 fhat(ii,ix) = lhat(ii)
133 end do
134
135 end do
136
137 return
138
139 c end subroutine GAD_PPM_HAT_X
140 end

  ViewVC Help
Powered by ViewVC 1.1.22