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

Contents of /MITgcm/pkg/generic_advdiff/gad_pqm_flx_y.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_FLX_Y(bi,bj,kk,ix,
7 & calc_CFL,delT,vvel,
8 & vfac,fhat,flux,myThid)
9 C |================================================================|
10 C | PQM_FLX_Y: evaluate PQM flux on grid-cell edges. |
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 ================================================================
21 C bi,bj :: tile indexing.
22 C kk :: r-index.
23 C ix :: x-index.
24 C calc_CFL :: TRUE to calc. CFL from vel.
25 C delT :: time-step.
26 C vvel :: vel.-comp in y-direction.
27 C vfac :: vel.-flux in y-direction.
28 C fhat :: row of poly. coeff.
29 C flux :: adv.-flux in y-direction.
30 C myThid :: thread number.
31 C ================================================================
32 integer bi,bj,kk,ix
33 logical calc_CFL
34 _RL delT
35 _RL vvel(1-OLx:sNx+OLx,
36 & 1-OLy:sNy+OLy)
37 _RL vfac(1-OLx:sNx+OLx,
38 & 1-OLy:sNy+OLy)
39 _RL fhat(1:5,
40 & 1-OLy:sNy+OLy)
41 _RL flux(1-OLx:sNx+OLx,
42 & 1-OLy:sNy+OLy)
43 integer myThid
44
45 C ================================================================
46 C iy :: y-indexing.
47 C vCFL :: CFL number.
48 C intF :: upwind tracer edge-value.
49 C ss11,ss22 :: int. endpoints.
50 C ivec :: int. basis vec.
51 C ================================================================
52 integer iy
53 _RL vCFL,intF
54 _RL ss11,ss22
55 _RL ivec(1:5)
56
57 C ================================================================
58 C (1): calc. "departure-points" for each grid-cell edge by int-
59 C egrating edge position backward in time over one single
60 C time-step. This is a "single-cell" implementation: requ-
61 C ires CFL <= 1.0.
62 C (2): calc. flux as the integral of the upwind grid-cell poly-
63 C nomial over the deformation interval found in (1).
64 C ================================================================
65
66 do iy = 1-OLy+4, sNy+OLy-3
67
68 if (vvel(ix,iy) .eq. 0. _d 0) then
69
70 flux(ix,iy) = 0. _d 0
71
72 else
73
74 if (vvel(ix,iy) .gt. 0. _d 0) then
75
76 C ==================== integrate PQM profile over upwind cell IY-1
77 if ( calc_CFL ) then
78 vCFL = vvel(ix,iy) * delT
79 & * recip_dyF(ix,iy-1,bi,bj)
80 & * recip_deepFacC(kk)
81 else
82 vCFL = vvel(ix,iy)
83 end if
84
85 ss11 = +1. _d 0 - 2. _d 0 * vCFL
86 ss22 = +1. _d 0
87
88 C ==================== integrate profile over region swept by face
89 ivec(1) = ss22 - ss11
90 ivec(2) =(ss22 ** 2
91 & - ss11 ** 2)*(1. _d 0 / 2. _d 0)
92 ivec(3) =(ss22 ** 3
93 & - ss11 ** 3)*(1. _d 0 / 3. _d 0)
94 ivec(4) =(ss22 ** 4
95 & - ss11 ** 4)*(1. _d 0 / 4. _d 0)
96 ivec(5) =(ss22 ** 5
97 & - ss11 ** 5)*(1. _d 0 / 5. _d 0)
98
99 intF = ivec(1) * fhat(1,iy-1)
100 & + ivec(2) * fhat(2,iy-1)
101 & + ivec(3) * fhat(3,iy-1)
102 & + ivec(4) * fhat(4,iy-1)
103 & + ivec(5) * fhat(5,iy-1)
104
105 intF = intF / (ss22 - ss11)
106
107 else
108
109 C ==================== integrate PQM profile over upwind cell IY+0
110 if ( calc_CFL ) then
111 vCFL = vvel(ix,iy) * delT
112 & * recip_dyF(ix,iy-0,bi,bj)
113 & * recip_deepFacC(kk)
114 else
115 vCFL = vvel(ix,iy)
116 end if
117
118 ss11 = -1. _d 0 - 2. _d 0 * vCFL
119 ss22 = -1. _d 0
120
121 C ==================== integrate profile over region swept by face
122 ivec(1) = ss22 - ss11
123 ivec(2) =(ss22 ** 2
124 & - ss11 ** 2)*(1. _d 0 / 2. _d 0)
125 ivec(3) =(ss22 ** 3
126 & - ss11 ** 3)*(1. _d 0 / 3. _d 0)
127 ivec(4) =(ss22 ** 4
128 & - ss11 ** 4)*(1. _d 0 / 4. _d 0)
129 ivec(5) =(ss22 ** 5
130 & - ss11 ** 5)*(1. _d 0 / 5. _d 0)
131
132 intF = ivec(1) * fhat(1,iy-0)
133 & + ivec(2) * fhat(2,iy-0)
134 & + ivec(3) * fhat(3,iy-0)
135 & + ivec(4) * fhat(4,iy-0)
136 & + ivec(5) * fhat(5,iy-0)
137
138 intF = intF / (ss22 - ss11)
139
140 end if
141
142 C ==================== calc. flux = upwind tracer * face-transport
143 flux(ix,iy) = + vfac(ix,iy) * intF
144
145 end if
146
147 end do
148
149 return
150
151 c end subroutine GAD_PQM_FLX_Y
152 end

  ViewVC Help
Powered by ViewVC 1.1.22