1 |
jmc |
1.1 |
C $Header: $ |
2 |
|
|
C $Name: $ |
3 |
|
|
|
4 |
|
|
# include "GAD_OPTIONS.h" |
5 |
|
|
|
6 |
|
|
SUBROUTINE GAD_PQM_ADV_Y(meth,bi,bj,kk, |
7 |
|
|
I calc_CFL,delT,vvel,vfac,fbar, |
8 |
|
|
O flux,myThid ) |
9 |
|
|
C |================================================================| |
10 |
|
|
C | PQM_ADV_Y: evaluate grid-cell advective flux in Y. | |
11 |
|
|
C | Lagrangian-type Piecewise Quartic Method (PQM). | |
12 |
|
|
C |================================================================| |
13 |
|
|
|
14 |
|
|
implicit none |
15 |
|
|
|
16 |
|
|
C =============================================== global variables |
17 |
|
|
# include "SIZE.h" |
18 |
|
|
# include "GRID.h" |
19 |
|
|
# include "GAD.h" |
20 |
|
|
|
21 |
|
|
C ================================================================ |
22 |
|
|
C meth :: advection method. |
23 |
|
|
C bi,bj :: tile indexing. |
24 |
|
|
C kk :: r-index. |
25 |
|
|
C calc_CFL :: TRUE to calc. CFL from vel. |
26 |
|
|
C delT :: time-step. |
27 |
|
|
C vvel :: vel.-comp in y-direction. |
28 |
|
|
C vfac :: vel.-flux in y-direction. |
29 |
|
|
C fbar :: grid-cell values. |
30 |
|
|
C flux :: adv.-flux in y-direction. |
31 |
|
|
C myThid :: thread number. |
32 |
|
|
C ================================================================ |
33 |
|
|
integer meth |
34 |
|
|
integer bi,bj,kk |
35 |
|
|
logical calc_CFL |
36 |
|
|
_RL delT |
37 |
|
|
_RL vvel(1-OLx:sNx+OLx, |
38 |
|
|
& 1-OLy:sNy+OLy) |
39 |
|
|
_RL vfac(1-OLx:sNx+OLx, |
40 |
|
|
& 1-OLy:sNy+OLy) |
41 |
|
|
_RL fbar(1-OLx:sNx+OLx, |
42 |
|
|
& 1-OLy:sNy+OLy) |
43 |
|
|
_RL flux(1-OLx:sNx+OLx, |
44 |
|
|
& 1-OLy:sNy+OLy) |
45 |
|
|
integer myThid |
46 |
|
|
|
47 |
|
|
C ================================================================ |
48 |
|
|
C ix,iy,ir :: grid indexing. |
49 |
|
|
C floc :: row of grid-cell values. |
50 |
|
|
C mloc :: row of grid-cell mask values. |
51 |
|
|
C fhat :: row of poly. coeff. |
52 |
|
|
C - FHAT(:,I) = PQM coeff. |
53 |
|
|
C edge :: row of edge-wise values/slopes. |
54 |
|
|
C - EDGE(1,:) = VALUE. |
55 |
|
|
C - EDGE(2,:) = DF/DY. |
56 |
|
|
C ohat :: row of oscl. coeff. |
57 |
|
|
C - OHAT(1,:) = D^1F/DS^1. |
58 |
|
|
C - OHAT(2,:) = D^2F/DS^2. |
59 |
|
|
C ================================================================ |
60 |
|
|
integer ix,iy |
61 |
|
|
_RL mloc(1-OLy:sNy+OLy) |
62 |
|
|
_RL floc(1-OLy:sNy+OLy) |
63 |
|
|
_RL fhat(1:5, |
64 |
|
|
& 1-OLy:sNy+OLy) |
65 |
|
|
_RL edge(1:2, |
66 |
|
|
& 1-OLy:sNy+OLy) |
67 |
|
|
_RL ohat(1:2, |
68 |
|
|
& 1-OLy:sNy+OLy) |
69 |
|
|
_RL vsum |
70 |
|
|
|
71 |
|
|
do ix = 1-OLx+0, sNx+OLx-0 |
72 |
|
|
C ==================== zero stencil "ghost" cells along boundaries |
73 |
|
|
flux(ix, +1-OLy+0) = 0. _d 0 |
74 |
|
|
flux(ix, +1-OLy+1) = 0. _d 0 |
75 |
|
|
flux(ix, +1-OLy+2) = 0. _d 0 |
76 |
|
|
flux(ix, +1-OLy+3) = 0. _d 0 |
77 |
|
|
flux(ix,sNy+OLy-0) = 0. _d 0 |
78 |
|
|
flux(ix,sNy+OLy-1) = 0. _d 0 |
79 |
|
|
flux(ix,sNy+OLy-2) = 0. _d 0 |
80 |
|
|
end do |
81 |
|
|
|
82 |
|
|
C ================================================================ |
83 |
|
|
C (1): copy a single row of data onto contiguous storage, treat |
84 |
|
|
C as a set of one-dimensional problems. |
85 |
|
|
C (2): calc. "oscillation-indicators" for each grid-cell if ad- |
86 |
|
|
C vection scheme is WENO-class. |
87 |
|
|
C (3): calc. edge-centred values/slopes by high-order interpol- |
88 |
|
|
C ation. |
89 |
|
|
C (4): calc. cell-centred polynomial profiles with appropriate |
90 |
|
|
C slope-limiting. |
91 |
|
|
C (5): calc. fluxes using a local, semi-lagrangian integration. |
92 |
|
|
C ================================================================ |
93 |
|
|
|
94 |
|
|
do ix = 1-OLx+0, sNx+OLx-0 |
95 |
|
|
|
96 |
|
|
vsum = 0.0 _d 0 |
97 |
|
|
do iy = 1-OLy+0, sNy+OLy-0 |
98 |
|
|
C ================================== quick break on zero transport |
99 |
|
|
vsum = vsum |
100 |
|
|
& + abs(vfac(ix,iy)) |
101 |
|
|
C ================================== make local unit-stride copies |
102 |
|
|
floc(iy) = fbar (ix,iy) |
103 |
|
|
mloc(iy) = |
104 |
|
|
& maskC(ix,iy,kk,bi,bj) |
105 |
|
|
end do |
106 |
|
|
|
107 |
|
|
if (vsum .gt. 0. _d 0) then |
108 |
|
|
|
109 |
|
|
C ==================== reconstruct derivatives for WENO indicators |
110 |
|
|
if (meth.eq.ENUM_PQM_WENO_LIMIT) then |
111 |
|
|
CALL GAD_OSC_HAT_Y(bi,bj,kk,ix, |
112 |
|
|
& mloc,floc, |
113 |
|
|
& ohat,myThid) |
114 |
|
|
end if |
115 |
|
|
|
116 |
|
|
C ==================== reconstruct 5th--order accurate edge values |
117 |
|
|
CALL GAD_PQM_P5E_Y(bi,bj,kk,ix, |
118 |
|
|
& mloc,floc, |
119 |
|
|
& edge,myThid) |
120 |
|
|
|
121 |
|
|
C ==================== reconstruct coeff. for grid-cell poynomials |
122 |
|
|
CALL GAD_PQM_HAT_Y(bi,bj,kk,ix, |
123 |
|
|
& meth, |
124 |
|
|
& mloc,floc, |
125 |
|
|
& edge,ohat, |
126 |
|
|
& fhat,myThid) |
127 |
|
|
|
128 |
|
|
C ==================== evaluate integral fluxes on grid-cell edges |
129 |
|
|
CALL GAD_PQM_FLX_Y(bi,bj,kk,ix, |
130 |
|
|
& calc_CFL, |
131 |
|
|
& delT,vvel, |
132 |
|
|
& vfac,fhat, |
133 |
|
|
& flux,myThid) |
134 |
|
|
|
135 |
|
|
else |
136 |
|
|
|
137 |
|
|
do iy = 1-OLy+4, sNy+OLy-3 |
138 |
|
|
C ================================== "null" flux on zero transport |
139 |
|
|
flux(ix,iy) = 0.0 _d 0 |
140 |
|
|
end do |
141 |
|
|
|
142 |
|
|
end if |
143 |
|
|
|
144 |
|
|
end do |
145 |
|
|
|
146 |
|
|
return |
147 |
|
|
|
148 |
|
|
c end subroutine GAD_PQM_ADV_Y |
149 |
|
|
end |