FV3DYCORE  Version 2.0.0
a2b_edge.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the FV3 dynamical core.
5 !*
6 !* The FV3 dynamical core is free software: you can redistribute it
7 !* and/or modify it under the terms of the
8 !* GNU Lesser General Public License as published by the
9 !* Free Software Foundation, either version 3 of the License, or
10 !* (at your option) any later version.
11 !*
12 !* The FV3 dynamical core is distributed in the hope that it will be
13 !* useful, but WITHOUT ANYWARRANTY; without even the implied warranty
14 !* of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
15 !* See the GNU General Public License for more details.
16 !*
17 !* You should have received a copy of the GNU Lesser General Public
18 !* License along with the FV3 dynamical core.
19 !* If not, see <http://www.gnu.org/licenses/>.
20 !***********************************************************************
21 
23 
25 !
26 ! Modules Included:
27 ! <table>
28 ! <tr>
29 ! <th>Module Name</th>
30 ! <th>Functions Included</th>
31 ! </tr>
32 ! <tr>
33 ! <td>fv_arrays_mod</td>
34 ! <td>fv_grid_type, R_GRID</td>
35 ! </tr>
36 ! <tr>
37 ! <td>fv_grid_utils_mod</td>
38 ! <td>great_circle_dist, van2 </td>
39 ! </tr>
40 ! </table>
41 
43 #ifdef VAN2
44  use fv_grid_utils_mod, only: van2
45 #endif
46 
48 
49  implicit none
50 
51  real, parameter:: r3 = 1./3.
52 !----------------------------
53 ! 4-pt Lagrange interpolation
54 !----------------------------
55  real, parameter:: a1 = 0.5625
56  real, parameter:: a2 = -0.0625
57 !----------------------
58 ! PPM volume mean form:
59 !----------------------
60  real, parameter:: b1 = 7./12.
61  real, parameter:: b2 = -1./12.
62 
63  private
64  public :: a2b_ord2, a2b_ord4
65 
66 contains
67 
68 #ifndef USE_OLD_ALGORITHM
69  subroutine a2b_ord4(qin, qout, gridstruct, npx, npy, is, ie, js, je, ng, replace)
70  integer, intent(IN):: npx, npy, is, ie, js, je, ng
71  real, intent(INOUT):: qin(is-ng:ie+ng,js-ng:je+ng)
72  real, intent(INOUT):: qout(is-ng:ie+ng,js-ng:je+ng)
73  type(fv_grid_type), intent(IN), target :: gridstruct
74  logical, optional, intent(IN):: replace
75 ! local: compact 4-pt cubic
76  real, parameter:: c1 = 2./3.
77  real, parameter:: c2 = -1./6.
78 ! Parabolic spline
79 ! real, parameter:: c1 = 0.75
80 ! real, parameter:: c2 = -0.25
81 
82  real qx(is:ie+1,js-ng:je+ng)
83  real qy(is-ng:ie+ng,js:je+1)
84  real qxx(is-ng:ie+ng,js-ng:je+ng)
85  real qyy(is-ng:ie+ng,js-ng:je+ng)
86  real g_in, g_ou
87  real:: p0(2)
88  real:: q1(is-1:ie+1), q2(js-1:je+1)
89  integer:: i, j, is1, js1, is2, js2, ie1, je1
90 
91  real, pointer, dimension(:,:,:) :: grid, agrid
92  real, pointer, dimension(:,:) :: dxa, dya
93  real(kind=R_GRID), pointer, dimension(:) :: edge_w, edge_e, edge_s, edge_n
94 
95  edge_w => gridstruct%edge_w
96  edge_e => gridstruct%edge_e
97  edge_s => gridstruct%edge_s
98  edge_n => gridstruct%edge_n
99 
100  grid => gridstruct%grid
101  agrid => gridstruct%agrid
102  dxa => gridstruct%dxa
103  dya => gridstruct%dya
104 
105  if (gridstruct%grid_type < 3) then
106 
107  is1 = max(1,is-1)
108  js1 = max(1,js-1)
109  is2 = max(2,is)
110  js2 = max(2,js)
111 
112  ie1 = min(npx-1,ie+1)
113  je1 = min(npy-1,je+1)
114 
115 ! Corners:
116 ! 3-way extrapolation
117  if (gridstruct%bounded_domain) then
118 
119  do j=js-2,je+2
120  do i=is,ie+1
121  qx(i,j) = b2*(qin(i-2,j)+qin(i+1,j)) + b1*(qin(i-1,j)+qin(i,j))
122  enddo
123  enddo
124 
125 
126  else
127 
128  if ( gridstruct%sw_corner ) then
129  p0(1:2) = grid(1,1,1:2)
130  qout(1,1) = (extrap_corner(p0, agrid(1,1,1:2), agrid( 2, 2,1:2), qin(1,1), qin( 2, 2)) + &
131  extrap_corner(p0, agrid(0,1,1:2), agrid(-1, 2,1:2), qin(0,1), qin(-1, 2)) + &
132  extrap_corner(p0, agrid(1,0,1:2), agrid( 2,-1,1:2), qin(1,0), qin( 2,-1)))*r3
133 
134  endif
135  if ( gridstruct%se_corner ) then
136  p0(1:2) = grid(npx,1,1:2)
137  qout(npx,1) = (extrap_corner(p0, agrid(npx-1,1,1:2), agrid(npx-2, 2,1:2), qin(npx-1,1), qin(npx-2, 2)) + &
138  extrap_corner(p0, agrid(npx-1,0,1:2), agrid(npx-2,-1,1:2), qin(npx-1,0), qin(npx-2,-1)) + &
139  extrap_corner(p0, agrid(npx ,1,1:2), agrid(npx+1, 2,1:2), qin(npx ,1), qin(npx+1, 2)))*r3
140  endif
141  if ( gridstruct%ne_corner ) then
142  p0(1:2) = grid(npx,npy,1:2)
143  qout(npx,npy) = (extrap_corner(p0, agrid(npx-1,npy-1,1:2), agrid(npx-2,npy-2,1:2), qin(npx-1,npy-1), qin(npx-2,npy-2)) + &
144  extrap_corner(p0, agrid(npx ,npy-1,1:2), agrid(npx+1,npy-2,1:2), qin(npx ,npy-1), qin(npx+1,npy-2)) + &
145  extrap_corner(p0, agrid(npx-1,npy ,1:2), agrid(npx-2,npy+1,1:2), qin(npx-1,npy ), qin(npx-2,npy+1)))*r3
146  endif
147  if ( gridstruct%nw_corner ) then
148  p0(1:2) = grid(1,npy,1:2)
149  qout(1,npy) = (extrap_corner(p0, agrid(1,npy-1,1:2), agrid( 2,npy-2,1:2), qin(1,npy-1), qin( 2,npy-2)) + &
150  extrap_corner(p0, agrid(0,npy-1,1:2), agrid(-1,npy-2,1:2), qin(0,npy-1), qin(-1,npy-2)) + &
151  extrap_corner(p0, agrid(1,npy, 1:2), agrid( 2,npy+1,1:2), qin(1,npy ), qin( 2,npy+1)))*r3
152  endif
153 
154 !------------
155 ! X-Interior:
156 !------------
157  do j=max(1,js-2),min(npy-1,je+2)
158  do i=max(3,is), min(npx-2,ie+1)
159  qx(i,j) = b2*(qin(i-2,j)+qin(i+1,j)) + b1*(qin(i-1,j)+qin(i,j))
160  enddo
161  enddo
162 
163  ! *** West Edges:
164  if ( is==1 ) then
165  do j=js1, je1
166  q2(j) = (qin(0,j)*dxa(1,j) + qin(1,j)*dxa(0,j))/(dxa(0,j) + dxa(1,j))
167  enddo
168  do j=js2, je1
169  qout(1,j) = edge_w(j)*q2(j-1) + (1.-edge_w(j))*q2(j)
170  enddo
171 !
172  do j=max(1,js-2),min(npy-1,je+2)
173  g_in = dxa(2,j) / dxa(1,j)
174  g_ou = dxa(-1,j) / dxa(0,j)
175  qx(1,j) = 0.5*( ((2.+g_in)*qin(1,j)-qin( 2,j))/(1.+g_in) + &
176  ((2.+g_ou)*qin(0,j)-qin(-1,j))/(1.+g_ou) )
177  qx(2,j) = ( 3.*(g_in*qin(1,j)+qin(2,j))-(g_in*qx(1,j)+qx(3,j)) ) / (2.+2.*g_in)
178  enddo
179  endif
180 
181  ! East Edges:
182  if ( (ie+1)==npx ) then
183  do j=js1, je1
184  q2(j) = (qin(npx-1,j)*dxa(npx,j) + qin(npx,j)*dxa(npx-1,j))/(dxa(npx-1,j) + dxa(npx,j))
185  enddo
186  do j=js2, je1
187  qout(npx,j) = edge_e(j)*q2(j-1) + (1.-edge_e(j))*q2(j)
188  enddo
189 !
190  do j=max(1,js-2),min(npy-1,je+2)
191  g_in = dxa(npx-2,j) / dxa(npx-1,j)
192  g_ou = dxa(npx+1,j) / dxa(npx,j)
193  qx(npx,j) = 0.5*( ((2.+g_in)*qin(npx-1,j)-qin(npx-2,j))/(1.+g_in) + &
194  ((2.+g_ou)*qin(npx, j)-qin(npx+1,j))/(1.+g_ou) )
195  qx(npx-1,j) = (3.*(qin(npx-2,j)+g_in*qin(npx-1,j)) - (g_in*qx(npx,j)+qx(npx-2,j)))/(2.+2.*g_in)
196  enddo
197  endif
198 
199  end if
200 !------------
201 ! Y-Interior:
202 !------------
203 
204  if (gridstruct%bounded_domain) then
205 
206 
207  do j=js,je+1
208  do i=is-2,ie+2
209  qy(i,j) = b2*(qin(i,j-2)+qin(i,j+1)) + b1*(qin(i,j-1) + qin(i,j))
210  enddo
211  enddo
212 
213  else
214 
215  do j=max(3,js),min(npy-2,je+1)
216  do i=max(1,is-2), min(npx-1,ie+2)
217  qy(i,j) = b2*(qin(i,j-2)+qin(i,j+1)) + b1*(qin(i,j-1) + qin(i,j))
218  enddo
219  enddo
220 
221  ! South Edges:
222  if ( js==1 ) then
223  do i=is1, ie1
224  q1(i) = (qin(i,0)*dya(i,1) + qin(i,1)*dya(i,0))/(dya(i,0) + dya(i,1))
225  enddo
226  do i=is2, ie1
227  qout(i,1) = edge_s(i)*q1(i-1) + (1.-edge_s(i))*q1(i)
228  enddo
229 !
230  do i=max(1,is-2),min(npx-1,ie+2)
231  g_in = dya(i,2) / dya(i,1)
232  g_ou = dya(i,-1) / dya(i,0)
233  qy(i,1) = 0.5*( ((2.+g_in)*qin(i,1)-qin(i,2))/(1.+g_in) + &
234  ((2.+g_ou)*qin(i,0)-qin(i,-1))/(1.+g_ou) )
235  qy(i,2) = (3.*(g_in*qin(i,1)+qin(i,2)) - (g_in*qy(i,1)+qy(i,3)))/(2.+2.*g_in)
236  enddo
237  endif
238 
239  ! North Edges:
240  if ( (je+1)==npy ) then
241  do i=is1, ie1
242  q1(i) = (qin(i,npy-1)*dya(i,npy) + qin(i,npy)*dya(i,npy-1))/(dya(i,npy-1)+dya(i,npy))
243  enddo
244  do i=is2, ie1
245  qout(i,npy) = edge_n(i)*q1(i-1) + (1.-edge_n(i))*q1(i)
246  enddo
247 !
248  do i=max(1,is-2),min(npx-1,ie+2)
249  g_in = dya(i,npy-2) / dya(i,npy-1)
250  g_ou = dya(i,npy+1) / dya(i,npy)
251  qy(i,npy) = 0.5*( ((2.+g_in)*qin(i,npy-1)-qin(i,npy-2))/(1.+g_in) + &
252  ((2.+g_ou)*qin(i,npy )-qin(i,npy+1))/(1.+g_ou) )
253  qy(i,npy-1) = (3.*(qin(i,npy-2)+g_in*qin(i,npy-1)) - (g_in*qy(i,npy)+qy(i,npy-2)))/(2.+2.*g_in)
254  enddo
255  endif
256 
257  end if
258 !--------------------------------------
259 
260  if (gridstruct%bounded_domain) then
261 
262  do j=js, je+1
263  do i=is,ie+1
264  qxx(i,j) = a2*(qx(i,j-2)+qx(i,j+1)) + a1*(qx(i,j-1)+qx(i,j))
265  enddo
266  enddo
267 
268  do j=js,je+1
269  do i=is,ie+1
270  qyy(i,j) = a2*(qy(i-2,j)+qy(i+1,j)) + a1*(qy(i-1,j)+qy(i,j))
271  enddo
272 
273  do i=is,ie+1
274  qout(i,j) = 0.5*(qxx(i,j) + qyy(i,j)) ! averaging
275  enddo
276  enddo
277 
278 
279 
280  else
281 
282  do j=max(3,js),min(npy-2,je+1)
283  do i=max(2,is),min(npx-1,ie+1)
284  qxx(i,j) = a2*(qx(i,j-2)+qx(i,j+1)) + a1*(qx(i,j-1)+qx(i,j))
285  enddo
286  enddo
287 
288  if ( js==1 ) then
289  do i=max(2,is),min(npx-1,ie+1)
290  qxx(i,2) = c1*(qx(i,1)+qx(i,2))+c2*(qout(i,1)+qxx(i,3))
291  enddo
292  endif
293  if ( (je+1)==npy ) then
294  do i=max(2,is),min(npx-1,ie+1)
295  qxx(i,npy-1) = c1*(qx(i,npy-2)+qx(i,npy-1))+c2*(qout(i,npy)+qxx(i,npy-2))
296  enddo
297  endif
298 
299 
300  do j=max(2,js),min(npy-1,je+1)
301  do i=max(3,is),min(npx-2,ie+1)
302  qyy(i,j) = a2*(qy(i-2,j)+qy(i+1,j)) + a1*(qy(i-1,j)+qy(i,j))
303  enddo
304  if ( is==1 ) qyy(2,j) = c1*(qy(1,j)+qy(2,j))+c2*(qout(1,j)+qyy(3,j))
305  if((ie+1)==npx) qyy(npx-1,j) = c1*(qy(npx-2,j)+qy(npx-1,j))+c2*(qout(npx,j)+qyy(npx-2,j))
306 
307  do i=max(2,is),min(npx-1,ie+1)
308  qout(i,j) = 0.5*(qxx(i,j) + qyy(i,j)) ! averaging
309  enddo
310  enddo
311 
312  end if
313 
314  else ! grid_type>=3
315 !------------------------
316 ! Doubly periodic domain:
317 !------------------------
318 ! X-sweep: PPM
319  do j=js-2,je+2
320  do i=is,ie+1
321  qx(i,j) = b1*(qin(i-1,j)+qin(i,j)) + b2*(qin(i-2,j)+qin(i+1,j))
322  enddo
323  enddo
324 ! Y-sweep: PPM
325  do j=js,je+1
326  do i=is-2,ie+2
327  qy(i,j) = b1*(qin(i,j-1)+qin(i,j)) + b2*(qin(i,j-2)+qin(i,j+1))
328  enddo
329  enddo
330 
331  do j=js,je+1
332  do i=is,ie+1
333  qout(i,j) = 0.5*( a1*(qx(i,j-1)+qx(i,j ) + qy(i-1,j)+qy(i, j)) + &
334  a2*(qx(i,j-2)+qx(i,j+1) + qy(i-2,j)+qy(i+1,j)) )
335  enddo
336  enddo
337  endif
338 
339  if ( present(replace) ) then
340  if ( replace ) then
341  do j=js,je+1
342  do i=is,ie+1
343  qin(i,j) = qout(i,j)
344  enddo
345  enddo
346  endif
347  endif
348 
349  end subroutine a2b_ord4
350 
351 #else
352 
353 ! Working version:
354  subroutine a2b_ord4(qin, qout, gridstruct, npx, npy, is, ie, js, je, ng, replace)
355  integer, intent(IN):: npx, npy, is, ie, js, je, ng
356  real, intent(INOUT):: qin(is-ng:ie+ng,js-ng:je+ng)
357  real, intent(INOUT):: qout(is-ng:ie+ng,js-ng:je+ng)
358  type(fv_grid_type), intent(IN), target :: gridstruct
359  logical, optional, intent(IN):: replace
360 ! local: compact 4-pt cubic
361  real, parameter:: c1 = 2./3.
362  real, parameter:: c2 = -1./6.
363 ! Parabolic spline
364 ! real, parameter:: c1 = 0.75
365 ! real, parameter:: c2 = -0.25
366 !-----------------------------
367 ! 6-pt corner interpolation:
368 !-----------------------------
369  real, parameter:: d1 = 0.375
370  real, parameter:: d2 = -1./24.
371 
372  real qx(is:ie+1,js-ng:je+ng)
373  real qy(is-ng:ie+ng,js:je+1)
374  real qxx(is-ng:ie+ng,js-ng:je+ng)
375  real qyy(is-ng:ie+ng,js-ng:je+ng)
376  real gratio
377  real g_in, g_ou
378  real:: p0(2)
379  real q1(npx), q2(npy)
380  integer:: i, j, is1, js1, is2, js2, ie1, je1
381 
382  real, pointer, dimension(:,:,:) :: grid, agrid
383  real, pointer, dimension(:,:) :: dxa, dya
384  real, pointer, dimension(:) :: edge_w, edge_e, edge_s, edge_n
385 
386  edge_w => gridstruct%edge_w
387  edge_e => gridstruct%edge_e
388  edge_s => gridstruct%edge_s
389  edge_n => gridstruct%edge_n
390 
391 
392  grid => gridstruct%grid
393  agrid => gridstruct%agrid
394  dxa => gridstruct%dxa
395  dya => gridstruct%dya
396 
397  if (gridstruct%grid_type < 3) then
398 
399  is1 = max(1,is-1)
400  js1 = max(1,js-1)
401  is2 = max(2,is)
402  js2 = max(2,js)
403 
404  ie1 = min(npx-1,ie+1)
405  je1 = min(npy-1,je+1)
406 
407 ! Corners:
408 #ifdef USE_3PT
409  if ( gridstruct%sw_corner ) qout(1, 1) = r3*(qin(1, 1)+qin(1, 0)+qin(0, 1))
410  if ( gridstruct%se_corner ) qout(npx, 1) = r3*(qin(npx-1, 1)+qin(npx-1, 0)+qin(npx, 1))
411  if ( gridstruct%ne_corner ) qout(npx,npy) = r3*(qin(npx-1,npy-1)+qin(npx,npy-1)+qin(npx-1,npy))
412  if ( gridstruct%nw_corner ) qout(1, npy) = r3*(qin(1, npy-1)+qin(0, npy-1)+qin(1, npy))
413 #else
414 
415 #ifdef SYMM_GRID
416 ! Symmetrical 6-point formular:
417  if ( gridstruct%sw_corner ) then
418  qout(1,1) = d1*(qin(1, 0) + qin( 0,1) + qin(1,1)) + &
419  d2*(qin(2,-1) + qin(-1,2) + qin(2,2))
420  endif
421  if ( gridstruct%se_corner ) then
422  qout(npx,1) = d1*(qin(npx-1, 0) + qin(npx-1,1) + qin(npx, 1)) + &
423  d2*(qin(npx-2,-1) + qin(npx-2,2) + qin(npx+1,2))
424  endif
425  if ( gridstruct%ne_corner ) then
426  qout(npx,npy) = d1*(qin(npx-1,npy-1) + qin(npx, npy-1) + qin(npx-1,npy)) + &
427  d2*(qin(npx-2,npy-2) + qin(npx+1,npy-2) + qin(npx-2,npy+1))
428  endif
429  if ( gridstruct%nw_corner ) then
430  qout(1,npy) = d1*(qin( 0,npy-1) + qin(1,npy-1) + qin(1,npy)) + &
431  d2*(qin(-1,npy-2) + qin(2,npy-2) + qin(2,npy+1))
432  endif
433 #else
434 ! 3-way extrapolation
435  if ( gridstruct%sw_corner ) then
436  p0(1:2) = grid(1,1,1:2)
437  qout(1,1) = (extrap_corner(p0, agrid(1,1,1:2), agrid( 2, 2,1:2), qin(1,1), qin( 2, 2)) + &
438  extrap_corner(p0, agrid(0,1,1:2), agrid(-1, 2,1:2), qin(0,1), qin(-1, 2)) + &
439  extrap_corner(p0, agrid(1,0,1:2), agrid( 2,-1,1:2), qin(1,0), qin( 2,-1)))*r3
440 
441  endif
442  if ( gridstruct%se_corner ) then
443  p0(1:2) = grid(npx,1,1:2)
444  qout(npx,1) = (extrap_corner(p0, agrid(npx-1,1,1:2), agrid(npx-2, 2,1:2), qin(npx-1,1), qin(npx-2, 2)) + &
445  extrap_corner(p0, agrid(npx-1,0,1:2), agrid(npx-2,-1,1:2), qin(npx-1,0), qin(npx-2,-1)) + &
446  extrap_corner(p0, agrid(npx ,1,1:2), agrid(npx+1, 2,1:2), qin(npx ,1), qin(npx+1, 2)))*r3
447  endif
448  if ( gridstruct%ne_corner ) then
449  p0(1:2) = grid(npx,npy,1:2)
450  qout(npx,npy) = (extrap_corner(p0, agrid(npx-1,npy-1,1:2), agrid(npx-2,npy-2,1:2), qin(npx-1,npy-1), qin(npx-2,npy-2)) + &
451  extrap_corner(p0, agrid(npx ,npy-1,1:2), agrid(npx+1,npy-2,1:2), qin(npx ,npy-1), qin(npx+1,npy-2)) + &
452  extrap_corner(p0, agrid(npx-1,npy ,1:2), agrid(npx-2,npy+1,1:2), qin(npx-1,npy ), qin(npx-2,npy+1)))*r3
453  endif
454  if ( gridstruct%nw_corner ) then
455  p0(1:2) = grid(1,npy,1:2)
456  qout(1,npy) = (extrap_corner(p0, agrid(1,npy-1,1:2), agrid( 2,npy-2,1:2), qin(1,npy-1), qin( 2,npy-2)) + &
457  extrap_corner(p0, agrid(0,npy-1,1:2), agrid(-1,npy-2,1:2), qin(0,npy-1), qin(-1,npy-2)) + &
458  extrap_corner(p0, agrid(1,npy, 1:2), agrid( 2,npy+1,1:2), qin(1,npy ), qin( 2,npy+1)))*r3
459  endif
460 #endif
461 #endif
462 
463 !------------
464 ! X-Interior:
465 !------------
466  if (gridstruct%bounded_domain) then
467 
468  do j=js-2,je+2
469  do i=is, ie+1
470  qx(i,j) = b2*(qin(i-2,j)+qin(i+1,j)) + b1*(qin(i-1,j)+qin(i,j))
471  enddo
472  enddo
473 
474 
475  else
476 
477  do j=max(1,js-2),min(npy-1,je+2)
478  do i=max(3,is), min(npx-2,ie+1)
479  qx(i,j) = b2*(qin(i-2,j)+qin(i+1,j)) + b1*(qin(i-1,j)+qin(i,j))
480  enddo
481  enddo
482 
483 ! West Edges:
484  if ( is==1 ) then
485 
486  do j=max(1,js-2),min(npy-1,je+2)
487  gratio = dxa(2,j) / dxa(1,j)
488 #ifdef SYMM_GRID
489  qx(1,j) = 0.5*((2.+gratio)*(qin(0,j)+qin(1,j))-(qin(-1,j)+qin(2,j))) / (1.+gratio)
490 #else
491  g_in = gratio
492  g_ou = dxa(-1,j) / dxa(0,j)
493  qx(1,j) = 0.5*( ((2.+g_in)*qin(1,j)-qin( 2,j))/(1.+g_in) + &
494  ((2.+g_ou)*qin(0,j)-qin(-1,j))/(1.+g_ou) )
495 #endif
496  qx(2,j) = ( 3.*(gratio*qin(1,j)+qin(2,j))-(gratio*qx(1,j)+qx(3,j)) ) / (2.+2.*gratio)
497  enddo
498 
499  do j=max(3,js),min(npy-2,je+1)
500  qout(1,j) = a2*(qx(1,j-2)+qx(1,j+1)) + a1*(qx(1,j-1)+qx(1,j))
501  enddo
502 
503  if( js==1 ) qout(1, 2) = c1*(qx(1,1)+qx(1,2)) + c2*(qout(1,1)+qout(1,3))
504  if((je+1)==npy) qout(1,npy-1) = c1*(qx(1,npy-2)+qx(1,npy-1)) + c2*(qout(1,npy-2)+qout(1,npy))
505  endif
506 
507 ! East Edges:
508  if ( (ie+1)==npx ) then
509 
510  do j=max(1,js-2),min(npy-1,je+2)
511  gratio = dxa(npx-2,j) / dxa(npx-1,j)
512 #ifdef SYMM_GRID
513  qx(npx,j) = 0.5*((2.+gratio)*(qin(npx-1,j)+qin(npx,j)) &
514  - (qin(npx-2,j)+qin(npx+1,j))) / (1.+gratio )
515 #else
516  g_in = gratio
517  g_ou = dxa(npx+1,j) / dxa(npx,j)
518  qx(npx,j) = 0.5*( ((2.+g_in)*qin(npx-1,j)-qin(npx-2,j))/(1.+g_in) + &
519  ((2.+g_ou)*qin(npx, j)-qin(npx+1,j))/(1.+g_ou) )
520 #endif
521  qx(npx-1,j) = (3.*(qin(npx-2,j)+gratio*qin(npx-1,j)) - (gratio*qx(npx,j)+qx(npx-2,j)))/(2.+2.*gratio)
522  enddo
523 
524  do j=max(3,js),min(npy-2,je+1)
525  qout(npx,j) = a2*(qx(npx,j-2)+qx(npx,j+1)) + a1*(qx(npx,j-1)+qx(npx,j))
526  enddo
527 
528  if(js==1) qout(npx,2) = c1*(qx(npx,1)+qx(npx,2))+c2*(qout(npx,1)+qout(npx,3))
529  if((je+1)==npy) qout(npx,npy-1) = &
530  c1*(qx(npx,npy-2)+qx(npx,npy-1))+c2*(qout(npx,npy-2)+qout(npx,npy))
531  endif
532 
533  endif
534 !------------
535 ! Y-Interior:
536 !------------
537  if (gridstruct%bounded_domain) then
538 
539  do j=js,je+1
540  do i=is-2, ie+2
541  qy(i,j) = b2*(qin(i,j-2)+qin(i,j+1)) + b1*(qin(i,j-1)+qin(i,j))
542  enddo
543  enddo
544 
545  else
546 
547  do j=max(3,js),min(npy-2,je+1)
548  do i=max(1,is-2), min(npx-1,ie+2)
549  qy(i,j) = b2*(qin(i,j-2)+qin(i,j+1)) + b1*(qin(i,j-1)+qin(i,j))
550  enddo
551  enddo
552 
553 ! South Edges:
554  if ( js==1 ) then
555 
556  do i=max(1,is-2),min(npx-1,ie+2)
557  gratio = dya(i,2) / dya(i,1)
558 #ifdef SYMM_GRID
559  qy(i,1) = 0.5*((2.+gratio)*(qin(i,0)+qin(i,1)) &
560  - (qin(i,-1)+qin(i,2))) / (1.+gratio )
561 #else
562  g_in = gratio
563  g_ou = dya(i,-1) / dya(i,0)
564  qy(i,1) = 0.5*( ((2.+g_in)*qin(i,1)-qin(i,2))/(1.+g_in) + &
565  ((2.+g_ou)*qin(i,0)-qin(i,-1))/(1.+g_ou) )
566 #endif
567  qy(i,2) = (3.*(gratio*qin(i,1)+qin(i,2)) - (gratio*qy(i,1)+qy(i,3)))/(2.+2.*gratio)
568  enddo
569 
570  do i=max(3,is),min(npx-2,ie+1)
571  qout(i,1) = a2*(qy(i-2,1)+qy(i+1,1)) + a1*(qy(i-1,1)+qy(i,1))
572  enddo
573 
574  if( is==1 ) qout(2,1) = c1*(qy(1,1)+qy(2,1))+c2*(qout(1,1)+qout(3,1))
575  if((ie+1)==npx) qout(npx-1,1) = c1*(qy(npx-2,1)+qy(npx-1,1))+c2*(qout(npx-2,1)+qout(npx,1))
576  endif
577 
578 
579 ! North Edges:
580  if ( (je+1)==npy ) then
581  do i=max(1,is-2),min(npx-1,ie+2)
582  gratio = dya(i,npy-2) / dya(i,npy-1)
583 #ifdef SYMM_GRID
584  qy(i,npy ) = 0.5*((2.+gratio)*(qin(i,npy-1)+qin(i,npy)) &
585  - (qin(i,npy-2)+qin(i,npy+1))) / (1.+gratio)
586 #else
587  g_in = gratio
588  g_ou = dya(i,npy+1) / dya(i,npy)
589  qy(i,npy) = 0.5*( ((2.+g_in)*qin(i,npy-1)-qin(i,npy-2))/(1.+g_in) + &
590  ((2.+g_ou)*qin(i,npy )-qin(i,npy+1))/(1.+g_ou) )
591 #endif
592  qy(i,npy-1) = (3.*(qin(i,npy-2)+gratio*qin(i,npy-1)) - (gratio*qy(i,npy)+qy(i,npy-2)))/(2.+2.*gratio)
593  enddo
594 
595  do i=max(3,is),min(npx-2,ie+1)
596  qout(i,npy) = a2*(qy(i-2,npy)+qy(i+1,npy)) + a1*(qy(i-1,npy)+qy(i,npy))
597  enddo
598 
599  if( is==1 ) qout(2,npy) = c1*(qy(1,npy)+qy(2,npy))+c2*(qout(1,npy)+qout(3,npy))
600  if((ie+1)==npx) qout(npx-1,npy) = c1*(qy(npx-2,npy)+qy(npx-1,npy))+c2*(qout(npx-2,npy)+qout(npx,npy))
601  endif
602 
603  end if
604 
605  if (gridstruct%bounded_domain) then
606 
607  do j=js,je+1
608  do i=is,ie+1
609  qxx(i,j) = a2*(qx(i,j-2)+qx(i,j+1)) + a1*(qx(i,j-1)+qx(i,j))
610  enddo
611  enddo
612 
613  do j=js,je+1
614  do i=is,ie+1
615  qyy(i,j) = a2*(qy(i-2,j)+qy(i+1,j)) + a1*(qy(i-1,j)+qy(i,j))
616  enddo
617 
618  do i=is,ie+1
619  qout(i,j) = 0.5*(qxx(i,j) + qyy(i,j)) ! averaging
620  enddo
621  enddo
622 
623 
624  else
625 
626  do j=max(3,js),min(npy-2,je+1)
627  do i=max(2,is),min(npx-1,ie+1)
628  qxx(i,j) = a2*(qx(i,j-2)+qx(i,j+1)) + a1*(qx(i,j-1)+qx(i,j))
629  enddo
630  enddo
631 
632  if ( js==1 ) then
633  do i=max(2,is),min(npx-1,ie+1)
634  qxx(i,2) = c1*(qx(i,1)+qx(i,2))+c2*(qout(i,1)+qxx(i,3))
635  enddo
636  endif
637  if ( (je+1)==npy ) then
638  do i=max(2,is),min(npx-1,ie+1)
639  qxx(i,npy-1) = c1*(qx(i,npy-2)+qx(i,npy-1))+c2*(qout(i,npy)+qxx(i,npy-2))
640  enddo
641  endif
642 
643 
644  do j=max(2,js),min(npy-1,je+1)
645  do i=max(3,is),min(npx-2,ie+1)
646  qyy(i,j) = a2*(qy(i-2,j)+qy(i+1,j)) + a1*(qy(i-1,j)+qy(i,j))
647  enddo
648  if ( is==1 ) qyy(2,j) = c1*(qy(1,j)+qy(2,j))+c2*(qout(1,j)+qyy(3,j))
649  if((ie+1)==npx) qyy(npx-1,j) = c1*(qy(npx-2,j)+qy(npx-1,j))+c2*(qout(npx,j)+qyy(npx-2,j))
650 
651  do i=max(2,is),min(npx-1,ie+1)
652  qout(i,j) = 0.5*(qxx(i,j) + qyy(i,j)) ! averaging
653  enddo
654  enddo
655 
656  endif
657 
658  else ! grid_type>=3
659 !------------------------
660 ! Doubly periodic domain:
661 !------------------------
662 ! X-sweep: PPM
663  do j=js-2,je+2
664  do i=is,ie+1
665  qx(i,j) = b1*(qin(i-1,j)+qin(i,j)) + b2*(qin(i-2,j)+qin(i+1,j))
666  enddo
667  enddo
668 ! Y-sweep: PPM
669  do j=js,je+1
670  do i=is-2,ie+2
671  qy(i,j) = b1*(qin(i,j-1)+qin(i,j)) + b2*(qin(i,j-2)+qin(i,j+1))
672  enddo
673  enddo
674 
675  do j=js,je+1
676  do i=is,ie+1
677  qout(i,j) = 0.5*( a1*(qx(i,j-1)+qx(i,j ) + qy(i-1,j)+qy(i, j)) + &
678  a2*(qx(i,j-2)+qx(i,j+1) + qy(i-2,j)+qy(i+1,j)) )
679  enddo
680  enddo
681  endif
682 
683  if ( present(replace) ) then
684  if ( replace ) then
685  do j=js,je+1
686  do i=is,ie+1
687  qin(i,j) = qout(i,j)
688  enddo
689  enddo
690  endif
691  endif
692 
693  end subroutine a2b_ord4
694 #endif
695 
696 
697  subroutine a2b_ord2(qin, qout, gridstruct, npx, npy, is, ie, js, je, ng, replace)
698  integer, intent(IN ) :: npx, npy, is, ie, js, je, ng
699  real , intent(INOUT) :: qin(is-ng:ie+ng,js-ng:je+ng)
700  real , intent( OUT) :: qout(is-ng:ie+ng,js-ng:je+ng)
701  type(fv_grid_type), intent(IN), target :: gridstruct
702  logical, optional, intent(IN) :: replace
703  ! local:
704  real q1(npx), q2(npy)
705  integer :: i,j
706  integer :: is1, js1, is2, js2, ie1, je1
707 
708  real, pointer, dimension(:,:,:) :: grid, agrid
709  real, pointer, dimension(:,:) :: dxa, dya
710 
711  real(kind=R_GRID), pointer, dimension(:) :: edge_w, edge_e, edge_s, edge_n
712 
713  edge_w => gridstruct%edge_w
714  edge_e => gridstruct%edge_e
715  edge_s => gridstruct%edge_s
716  edge_n => gridstruct%edge_n
717 
718  grid => gridstruct%grid
719  agrid => gridstruct%agrid
720  dxa => gridstruct%dxa
721  dya => gridstruct%dya
722 
723  if (gridstruct%grid_type < 3) then
724 
725  if (gridstruct%bounded_domain) then
726 
727  do j=js-2,je+1+2
728  do i=is-2,ie+1+2
729  qout(i,j) = 0.25*(qin(i-1,j-1)+qin(i,j-1)+qin(i-1,j)+qin(i,j))
730  enddo
731  enddo
732 
733  else
734 
735  is1 = max(1,is-1)
736  js1 = max(1,js-1)
737  is2 = max(2,is)
738  js2 = max(2,js)
739 
740  ie1 = min(npx-1,ie+1)
741  je1 = min(npy-1,je+1)
742 
743  do j=js2,je1
744  do i=is2,ie1
745  qout(i,j) = 0.25*(qin(i-1,j-1)+qin(i,j-1)+qin(i-1,j)+qin(i,j))
746  enddo
747  enddo
748 
749 ! Fix the 4 Corners:
750  if ( gridstruct%sw_corner ) qout(1, 1) = r3*(qin(1, 1)+qin(1, 0)+qin(0, 1))
751  if ( gridstruct%se_corner ) qout(npx, 1) = r3*(qin(npx-1, 1)+qin(npx-1, 0)+qin(npx, 1))
752  if ( gridstruct%ne_corner ) qout(npx,npy) = r3*(qin(npx-1,npy-1)+qin(npx,npy-1)+qin(npx-1,npy))
753  if ( gridstruct%nw_corner ) qout(1, npy) = r3*(qin(1, npy-1)+qin(0, npy-1)+qin(1, npy))
754 
755  ! *** West Edges:
756  if ( is==1 ) then
757  do j=js1, je1
758  q2(j) = 0.5*(qin(0,j) + qin(1,j))
759  enddo
760  do j=js2, je1
761  qout(1,j) = edge_w(j)*q2(j-1) + (1.-edge_w(j))*q2(j)
762  enddo
763  endif
764 
765  ! East Edges:
766  if ( (ie+1)==npx ) then
767  do j=js1, je1
768  q2(j) = 0.5*(qin(npx-1,j) + qin(npx,j))
769  enddo
770  do j=js2, je1
771  qout(npx,j) = edge_e(j)*q2(j-1) + (1.-edge_e(j))*q2(j)
772  enddo
773  endif
774 
775  ! South Edges:
776  if ( js==1 ) then
777  do i=is1, ie1
778  q1(i) = 0.5*(qin(i,0) + qin(i,1))
779  enddo
780  do i=is2, ie1
781  qout(i,1) = edge_s(i)*q1(i-1) + (1.-edge_s(i))*q1(i)
782  enddo
783  endif
784 
785  ! North Edges:
786  if ( (je+1)==npy ) then
787  do i=is1, ie1
788  q1(i) = 0.5*(qin(i,npy-1) + qin(i,npy))
789  enddo
790  do i=is2, ie1
791  qout(i,npy) = edge_n(i)*q1(i-1) + (1.-edge_n(i))*q1(i)
792  enddo
793  endif
794 
795  end if
796 
797  else
798 
799  do j=js,je+1
800  do i=is,ie+1
801  qout(i,j) = 0.25*(qin(i-1,j-1)+qin(i,j-1)+qin(i-1,j)+qin(i,j))
802  enddo
803  enddo
804 
805  endif
806 
807 
808  if ( present(replace) ) then
809  if ( replace ) then
810  do j=js,je+1
811  do i=is,ie+1
812  qin(i,j) = qout(i,j)
813  enddo
814  enddo
815  endif
816  endif
817 
818  end subroutine a2b_ord2
819 
820  real function extrap_corner ( p0, p1, p2, q1, q2 )
821  real, intent(in ), dimension(2):: p0, p1, p2
822  real, intent(in ):: q1, q2
823  real:: x1, x2
824 
825  x1 = great_circle_dist( real(p1,kind=R_GRID), real(p0,kind=R_GRID) )
826  x2 = great_circle_dist( real(p2,kind=R_GRID), real(p0,kind=R_GRID) )
827 
828  extrap_corner = q1 + x1/(x2-x1) * (q1-q2)
829 
830  end function extrap_corner
831 
832 #ifdef TEST_VAND2
833  subroutine a2b_ord4(qin, qout, grid, agrid, npx, npy, is, ie, js, je, ng, replace)
834 ! use tp_core_mod, only: copy_corners
835  integer, intent(IN):: npx, npy, is, ie, js, je, ng
836  real, intent(INOUT):: qin(is-ng:ie+ng,js-ng:je+ng)
837  real, intent(INOUT):: qout(is-ng:ie+ng,js-ng:je+ng)
838  real, intent(in) :: grid(is-ng:ie+ng+1,js-ng:je+ng+1,2)
839  real, intent(in) :: agrid(is-ng:ie+ng,js-ng:je+ng,2)
840  logical, optional, intent(IN):: replace
841  real qx(is:ie+1,js-ng:je+ng)
842  real qy(is-ng:ie+ng,js:je+1)
843  real:: p0(2)
844  integer :: i, j
845 
846  real, pointer, dimension(:,:,:) :: grid, agrid
847  real, pointer, dimension(:,:) :: dxa, dya
848 
849  real, pointer, dimension(:) :: edge_w, edge_e, edge_s, edge_n
850 
851  edge_w => gridstruct%edge_w
852  edge_e => gridstruct%edge_e
853  edge_s => gridstruct%edge_s
854  edge_n => gridstruct%edge_n
855 
856  grid => gridstruct%grid
857  agrid => gridstruct%agrid
858  dxa => gridstruct%dxa
859  dya => gridstruct%dya
860 
861 
862  if (gridstruct%grid_type < 3) then
863 
864 !------------------------------------------
865 ! Copy fields to the phantom corner region:
866 !------------------------------------------
867 ! call copy_corners(qin, npx, npy, 1)
868 
869  do j=js,je+1
870  do i=is,ie+1
871 !SW:
872  if ( i==1 .and. j==1 ) goto 123
873  if ( i==2 .and. j==1 ) then
874  qin(0,-1) = qin(-1,2)
875  qin(0, 0) = qin(-1,1)
876  endif
877  if ( i==1 .and. j==2 ) then
878  qin(-1,0) = qin(2,-1)
879  qin( 0,0) = qin(1,-1)
880  endif
881  if ( i==2 .and. j==2 ) then
882  qin( 0,0) = qin(4,4)
883  endif
884 !SE:
885  if ( i==npx .and. j==1 ) goto 123
886  if ( i==npx-1 .and. j==1 ) then
887  qin(npx,-1) = qin(npx+1,2)
888  qin(npx, 0) = qin(npx+1,1)
889  endif
890  if ( i==npx-1 .and. j==2 ) then
891  qin(npx,0) = qin(npx-4,4)
892  endif
893  if ( i==npx .and. j==2 ) then
894  qin(npx+1,0) = qin(npx-2,-1)
895  qin(npx, 0) = qin(npx-1,-1)
896  endif
897 !NE:
898  if ( i==npx .and. j==npy ) goto 123
899  if ( i==npx-1 .and. j==npy-1 ) then
900  qin(npx,npy) = qin(npx-4,npy-4)
901  endif
902  if ( i==npx .and. j==npy-1 ) then
903  qin(npx+1,npy) = qin(npx-2,npy+1)
904  qin(npx, npy) = qin(npx-1,npy+1)
905  endif
906  if ( i==npx-1 .and. j==npy ) then
907  qin(npx,npy+1) = qin(npx+1,npy-2)
908  qin(npx,npy ) = qin(npx+1,npy-1)
909  endif
910 !NW:
911  if ( i==1 .and. j==npy ) goto 123
912  if ( i==1 .and. j==npy-1 ) then
913  qin(-1,npy) = qin(2,npy+1)
914  qin( 0,npy) = qin(1,npy+1)
915  endif
916  if ( i==2 .and. j==npy-1 ) then
917  qin(0,npy) = qin(4,npy-4)
918  endif
919  if ( i==2 .and. j==npy ) then
920  qin(0,npy+1) = qin(-1,npy-2)
921  qin(0,npy ) = qin(-1,npy-1)
922  endif
923 
924  qout(i,j) = van2(1, i,j)*qin(i-2,j-2) + van2(2, i,j)*qin(i-1,j-2) + &
925  van2(3, i,j)*qin(i ,j-2) + van2(4, i,j)*qin(i+1,j-2) + &
926  van2(5, i,j)*qin(i-2,j-1) + van2(6, i,j)*qin(i-1,j-1) + &
927  van2(7, i,j)*qin(i ,j-1) + van2(8, i,j)*qin(i+1,j-1) + &
928  van2(9, i,j)*qin(i-2,j ) + van2(10,i,j)*qin(i-1,j ) + &
929  van2(11,i,j)*qin(i ,j ) + van2(12,i,j)*qin(i+1,j ) + &
930  van2(13,i,j)*qin(i-2,j+1) + van2(14,i,j)*qin(i-1,j+1) + &
931  van2(15,i,j)*qin(i ,j+1) + van2(16,i,j)*qin(i+1,j+1)
932 123 continue
933  enddo
934  enddo
935 
936 ! 3-way extrapolation
937  if ( gridstruct%sw_corner ) then
938  p0(1:2) = grid(1,1,1:2)
939  qout(1,1) = (extrap_corner(p0, agrid(1,1,1:2), agrid( 2, 2,1:2), qin(1,1), qin( 2, 2)) + &
940  extrap_corner(p0, agrid(0,1,1:2), agrid(-1, 2,1:2), qin(0,1), qin(-1, 2)) + &
941  extrap_corner(p0, agrid(1,0,1:2), agrid( 2,-1,1:2), qin(1,0), qin( 2,-1)))*r3
942 
943  endif
944  if ( gridstruct%se_corner ) then
945  p0(1:2) = grid(npx,1,1:2)
946  qout(npx,1) = (extrap_corner(p0, agrid(npx-1,1,1:2), agrid(npx-2, 2,1:2), qin(npx-1,1), qin(npx-2, 2)) + &
947  extrap_corner(p0, agrid(npx-1,0,1:2), agrid(npx-2,-1,1:2), qin(npx-1,0), qin(npx-2,-1)) + &
948  extrap_corner(p0, agrid(npx ,1,1:2), agrid(npx+1, 2,1:2), qin(npx ,1), qin(npx+1, 2)))*r3
949  endif
950  if ( gridstruct%ne_corner ) then
951  p0(1:2) = grid(npx,npy,1:2)
952  qout(npx,npy) = (extrap_corner(p0, agrid(npx-1,npy-1,1:2), agrid(npx-2,npy-2,1:2), qin(npx-1,npy-1), qin(npx-2,npy-2)) + &
953  extrap_corner(p0, agrid(npx ,npy-1,1:2), agrid(npx+1,npy-2,1:2), qin(npx ,npy-1), qin(npx+1,npy-2)) + &
954  extrap_corner(p0, agrid(npx-1,npy ,1:2), agrid(npx-2,npy+1,1:2), qin(npx-1,npy ), qin(npx-2,npy+1)))*r3
955  endif
956  if ( gridstruct%nw_corner ) then
957  p0(1:2) = grid(1,npy,1:2)
958  qout(1,npy) = (extrap_corner(p0, agrid(1,npy-1,1:2), agrid( 2,npy-2,1:2), qin(1,npy-1), qin( 2,npy-2)) + &
959  extrap_corner(p0, agrid(0,npy-1,1:2), agrid(-1,npy-2,1:2), qin(0,npy-1), qin(-1,npy-2)) + &
960  extrap_corner(p0, agrid(1,npy, 1:2), agrid( 2,npy+1,1:2), qin(1,npy ), qin( 2,npy+1)))*r3
961  endif
962 
963  else ! grid_type>=3
964 
965 !------------------------
966 ! Doubly periodic domain:
967 !------------------------
968 ! X-sweep: PPM
969  do j=js-2,je+2
970  do i=is,ie+1
971  qx(i,j) = b1*(qin(i-1,j)+qin(i,j)) + b2*(qin(i-2,j)+qin(i+1,j))
972  enddo
973  enddo
974 ! Y-sweep: PPM
975  do j=js,je+1
976  do i=is-2,ie+2
977  qy(i,j) = b1*(qin(i,j-1)+qin(i,j)) + b2*(qin(i,j-2)+qin(i,j+1))
978  enddo
979  enddo
980 
981  do j=js,je+1
982  do i=is,ie+1
983  qout(i,j) = 0.5*( a1*(qx(i,j-1)+qx(i,j ) + qy(i-1,j)+qy(i, j)) + &
984  a2*(qx(i,j-2)+qx(i,j+1) + qy(i-2,j)+qy(i+1,j)) )
985  enddo
986  enddo
987 
988  endif
989 
990  if ( present(replace) ) then
991  if ( replace ) then
992  do j=js,je+1
993  do i=is,ie+1
994  qin(i,j) = qout(i,j)
995  enddo
996  enddo
997  endif
998  endif
999 
1000  end subroutine a2b_ord4
1001 #endif
1002 
1003 end module a2b_edge_mod
real, parameter b1
0.58333333
Definition: a2b_edge.F90:60
The type &#39;fv_grid_type&#39; is made up of grid-dependent information from fv_grid_tools and fv_grid_utils...
Definition: fv_arrays.F90:123
subroutine, public a2b_ord2(qin, qout, gridstruct, npx, npy, is, ie, js, je, ng, replace)
Definition: a2b_edge.F90:698
real, parameter b2
Definition: a2b_edge.F90:61
integer, parameter, public r_grid
Definition: fv_arrays.F90:34
The module &#39;a2b_edge&#39; performs FV-consistent interpolation of pressure to corners.
Definition: a2b_edge.F90:24
The module &#39;fv_arrays&#39; contains the &#39;fv_atmos_type&#39; and associated datatypes.
Definition: fv_arrays.F90:24
real function, public great_circle_dist(q1, q2, radius)
real, parameter r3
Definition: a2b_edge.F90:51
The module &#39;fv_grid_utils&#39; contains routines for setting up and computing grid-related quantities...
subroutine, public a2b_ord4(qin, qout, gridstruct, npx, npy, is, ie, js, je, ng, replace)
Definition: a2b_edge.F90:70
real, parameter a1
9/16
Definition: a2b_edge.F90:55
real, parameter a2
-1/16
Definition: a2b_edge.F90:56
real function extrap_corner(p0, p1, p2, q1, q2)
Definition: a2b_edge.F90:821