FV3 Bundle
fv_restart_adm.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU General Public License *
3 !* This file is a part of fvGFS. *
4 !* *
5 !* fvGFS is free software; you can redistribute it and/or modify it *
6 !* and are expected to follow the terms of the GNU General Public *
7 !* License as published by the Free Software Foundation; either *
8 !* version 2 of the License, or (at your option) any later version. *
9 !* *
10 !* fvGFS is distributed in the hope that it will be useful, but *
11 !* WITHOUT ANY WARRANTY; without even the implied warranty of *
12 !* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU *
13 !* General Public License for more details. *
14 !* *
15 !* For the full text of the GNU General Public License, *
16 !* write to: Free Software Foundation, Inc., *
17 !* 675 Mass Ave, Cambridge, MA 02139, USA. *
18 !* or see: http://www.gnu.org/licenses/gpl.html *
19 !***********************************************************************
21 
24 
25 implicit none
26 private
27 
28 public :: d2c_setup, d2a_setup
30 
31 CONTAINS
32 ! Differentiation of d2c_setup in reverse (adjoint) mode, forward sweep (with options split(a2b_edge_mod.a2b_ord4 a2b_edge_mo
33 !d.a2b_ord2 dyn_core_mod.dyn_core dyn_core_mod.pk3_halo dyn_core_mod.pln_halo dyn_core_mod.pe_halo dyn_core_mod.adv_pe dyn_core_m
34 !od.p_grad_c dyn_core_mod.nh_p_grad dyn_core_mod.split_p_grad dyn_core_mod.one_grad_p dyn_core_mod.grad1_p_update dyn_core_mod.mi
35 !x_dp dyn_core_mod.geopk dyn_core_mod.del2_cubed dyn_core_mod.Rayleigh_fast fv_dynamics_mod.fv_dynamics fv_dynamics_mod.Rayleigh_
36 !Super fv_dynamics_mod.Rayleigh_Friction fv_dynamics_mod.compute_aam fv_grid_utils_mod.cubed_to_latlon fv_grid_utils_mod.c2l_ord4
37 ! fv_grid_utils_mod.c2l_ord2 fv_mapz_mod.Lagrangian_to_Eulerian fv_mapz_mod.compute_total_energy fv_mapz_mod.pkez fv_mapz_mod.rem
38 !ap_z fv_mapz_mod.map_scalar fv_mapz_mod.map1_ppm fv_mapz_mod.mapn_tracer fv_mapz_mod.map1_q2 fv_mapz_mod.remap_2d fv
39 !_mapz_mod.scalar_profile fv_mapz_mod.cs_profile fv_mapz_mod.cs_limiters fv_mapz_mod.ppm_profile fv_mapz_mod.ppm_limiters f
40 !v_mapz_mod.steepz fv_mapz_mod.rst_remap fv_mapz_mod.mappm fv_mapz_mod.moist_cv fv_mapz_mod.moist_cp fv_mapz_mod.map1_cubic fv_re
41 !start_mod.d2c_setup fv_tracer2d_mod.tracer_2d_1L fv_tracer2d_mod.tracer_2d fv_tracer2d_mod.tracer_2d_nested fv_sg_mod.fv_subgrid
42 !_z main_mod.compute_pressures main_mod.run nh_core_mod.Riem_Solver3 nh_utils_mod.update_dz_c nh_utils_mod.update_dz_d nh_utils_m
43 !od.Riem_Solver_c nh_utils_mod.Riem_Solver3test nh_utils_mod.imp_diff_w nh_utils_mod.RIM_2D nh_utils_mod.SIM3_solver nh_utils_mod
44 !.SIM3p0_solver nh_utils_mod.SIM1_solver nh_utils_mod.SIM_solver nh_utils_mod.edge_scalar nh_utils_mod.edge_profile nh_utils_mod.
45 !nest_halo_nh sw_core_mod.c_sw sw_core_mod.d_sw sw_core_mod.divergence_corner sw_core_mod.divergence_corner_nest sw_core_mod.d2a
46 !2c_vect sw_core_mod.fill3_4corners sw_core_mod.fill2_4corners sw_core_mod.fill_4corners sw_core_mod.xtp_u sw_core_mod.ytp_v_f
47 !b sw_core_mod.compute_divergence_damping sw_core_mod.smag_corner tp_core_mod.mp_ghost_ew tp_core_mod.fv_tp_2d tp_core_m
48 !od.copy_corners tp_core_mod.xppm tp_core_mod.yppm tp_core_mod.deln_flux a2b_edge_mod.extrap_corner fv_grid_utils_
49 !mod.great_circle_dist sw_core_mod.edge_interpolate4)):
50 ! gradient of useful results: u v uc vc
51 ! with respect to varying inputs: u v uc vc
52  SUBROUTINE d2c_setup_fwd(u, v, ua, va, uc, vc, dord4, isd, ied, jsd, &
53 & jed, is, ie, js, je, npx, npy, grid_type, nested, se_corner, &
54 & sw_corner, ne_corner, nw_corner, rsin_u, rsin_v, cosa_s, rsin2)
55  IMPLICIT NONE
56  LOGICAL, INTENT(IN) :: dord4
57  INTEGER, INTENT(IN) :: isd, ied, jsd, jed, is, ie, js, je, npx, npy&
58 & , grid_type
59  REAL, INTENT(IN) :: u(isd:ied, jsd:jed+1)
60  REAL, INTENT(IN) :: v(isd:ied+1, jsd:jed)
61  REAL, DIMENSION(isd:ied, jsd:jed), INTENT(OUT) :: ua
62  REAL, DIMENSION(isd:ied, jsd:jed), INTENT(OUT) :: va
63  REAL, DIMENSION(isd:ied+1, jsd:jed) :: uc
64  REAL, DIMENSION(isd:ied, jsd:jed+1) :: vc
65  LOGICAL, INTENT(IN) :: nested, se_corner, sw_corner, ne_corner, &
66 & nw_corner
67  REAL, INTENT(IN) :: rsin_u(isd:ied+1, jsd:jed)
68  REAL, INTENT(IN) :: rsin_v(isd:ied, jsd:jed+1)
69  REAL, INTENT(IN) :: cosa_s(isd:ied, jsd:jed)
70  REAL, INTENT(IN) :: rsin2(isd:ied, jsd:jed)
71 ! Local
72  REAL, DIMENSION(isd:ied, jsd:jed) :: utmp, vtmp
73  REAL, PARAMETER :: t11=27./28., t12=-(13./28.), t13=3./7., t14=6./7.&
74 & , t15=3./28.
75  REAL, PARAMETER :: a1=0.5625
76  REAL, PARAMETER :: a2=-0.0625
77  REAL, PARAMETER :: c1=-(2./14.)
78  REAL, PARAMETER :: c2=11./14.
79  REAL, PARAMETER :: c3=5./14.
80  INTEGER :: npt, i, j, ifirst, ilast, id
81  INTRINSIC max
82  INTRINSIC min
83  INTEGER :: max1
84  INTEGER :: max2
85  INTEGER :: max3
86  INTEGER :: max4
87  INTEGER :: max5
88  INTEGER :: max6
89  INTEGER :: min1
90  INTEGER :: min2
91  INTEGER :: min3
92  INTEGER :: min4
93  INTEGER :: min5
94  INTEGER :: min6
95  INTEGER :: ad_from
96  INTEGER :: ad_from0
97 
98  utmp = 0.0
99  vtmp = 0.0
100  npt = 0
101  ifirst = 0
102  ilast = 0
103  id = 0
104  max1 = 0
105  max2 = 0
106  max3 = 0
107  max4 = 0
108  max5 = 0
109  max6 = 0
110  min1 = 0
111  min2 = 0
112  min3 = 0
113  min4 = 0
114  min5 = 0
115  min6 = 0
116  ad_from = 0
117  ad_from0 = 0
118 
119  IF (grid_type .LT. 3 .AND. (.NOT.nested)) THEN
120  npt = 4
121  ELSE
122  npt = -2
123  END IF
124  IF (nested) THEN
125  CALL pushcontrol(2,0)
126  ELSE
127 !----------
128 ! Interior:
129 !----------
130  IF (npt .LT. js - 1) THEN
131  max1 = js - 1
132  ELSE
133  max1 = npt
134  END IF
135  IF (npy - npt .GT. je + 1) THEN
136  min1 = je + 1
137  ELSE
138  min1 = npy - npt
139  END IF
140  DO j=max1,min1
141  IF (npt .LT. isd) THEN
142  max2 = isd
143  ELSE
144  max2 = npt
145  END IF
146  IF (npx - npt .GT. ied) THEN
147  min2 = ied
148  ELSE
149  min2 = npx - npt
150  END IF
151  ad_from = max2
152  i = min2 + 1
153  CALL pushinteger(i - 1)
154  CALL pushinteger(ad_from)
155  END DO
156  IF (npt .LT. jsd) THEN
157  max3 = jsd
158  ELSE
159  max3 = npt
160  END IF
161  IF (npy - npt .GT. jed) THEN
162  min3 = jed
163  ELSE
164  min3 = npy - npt
165  END IF
166  DO j=max3,min3
167  IF (npt .LT. is - 1) THEN
168  max4 = is - 1
169  ELSE
170  max4 = npt
171  END IF
172  IF (npx - npt .GT. ie + 1) THEN
173  min4 = ie + 1
174  ELSE
175  min4 = npx - npt
176  END IF
177  ad_from0 = max4
178  i = min4 + 1
179  CALL pushinteger(i - 1)
180  CALL pushinteger(ad_from0)
181  END DO
182 !----------
183 ! edges:
184 !----------
185  IF (grid_type .LT. 3) THEN
186  IF (js .EQ. 1 .OR. jsd .LT. npt) THEN
187  CALL pushcontrol(1,0)
188  ELSE
189  CALL pushcontrol(1,1)
190  END IF
191  IF (je + 1 .EQ. npy .OR. jed .GE. npy - npt) THEN
192  CALL pushcontrol(1,0)
193  ELSE
194  CALL pushcontrol(1,1)
195  END IF
196  IF (is .EQ. 1 .OR. isd .LT. npt) THEN
197  IF (npt .LT. jsd) THEN
198  max5 = jsd
199  ELSE
200  max5 = npt
201  END IF
202  IF (npy - npt .GT. jed) THEN
203  min5 = jed
204  ELSE
205  min5 = npy - npt
206  END IF
207  CALL pushcontrol(1,0)
208  ELSE
209  CALL pushcontrol(1,1)
210  END IF
211  IF (ie + 1 .EQ. npx .OR. ied .GE. npx - npt) THEN
212  IF (npt .LT. jsd) THEN
213  max6 = jsd
214  ELSE
215  max6 = npt
216  END IF
217  IF (npy - npt .GT. jed) THEN
218  min6 = jed
219  ELSE
220  min6 = npy - npt
221  END IF
222  CALL pushcontrol(2,3)
223  ELSE
224  CALL pushcontrol(2,2)
225  END IF
226  ELSE
227  CALL pushcontrol(2,1)
228  END IF
229  END IF
230 ! A -> C
231 !--------------
232 ! Fix the edges
233 !--------------
234 ! Xdir:
235  IF (sw_corner) THEN
236  CALL pushcontrol(1,0)
237  ELSE
238  CALL pushcontrol(1,1)
239  END IF
240  IF (se_corner) THEN
241  CALL pushcontrol(1,0)
242  ELSE
243  CALL pushcontrol(1,1)
244  END IF
245  IF (ne_corner) THEN
246  CALL pushcontrol(1,0)
247  ELSE
248  CALL pushcontrol(1,1)
249  END IF
250  IF (nw_corner) THEN
251  CALL pushcontrol(1,0)
252  ELSE
253  CALL pushcontrol(1,1)
254  END IF
255  IF (grid_type .LT. 3 .AND. (.NOT.nested)) THEN
256  IF (3 .LT. is - 1) THEN
257  ifirst = is - 1
258  ELSE
259  ifirst = 3
260  END IF
261  IF (npx - 2 .GT. ie + 2) THEN
262  CALL pushcontrol(1,1)
263  ilast = ie + 2
264  ELSE
265  CALL pushcontrol(1,1)
266  ilast = npx - 2
267  END IF
268  ELSE
269  CALL pushcontrol(1,0)
270  ifirst = is - 1
271  ilast = ie + 2
272  END IF
273  IF (grid_type .LT. 3) THEN
274 ! Xdir:
275  IF (is .EQ. 1 .AND. (.NOT.nested)) THEN
276  CALL pushcontrol(1,0)
277  ELSE
278  CALL pushcontrol(1,1)
279  END IF
280  IF (ie + 1 .EQ. npx .AND. (.NOT.nested)) THEN
281  CALL pushcontrol(2,0)
282  ELSE
283  CALL pushcontrol(2,1)
284  END IF
285  ELSE
286  CALL pushcontrol(2,2)
287  END IF
288 !------
289 ! Ydir:
290 !------
291  IF (sw_corner) THEN
292  CALL pushcontrol(1,0)
293  ELSE
294  CALL pushcontrol(1,1)
295  END IF
296  IF (nw_corner) THEN
297  CALL pushcontrol(1,0)
298  ELSE
299  CALL pushcontrol(1,1)
300  END IF
301  IF (se_corner) THEN
302  CALL pushcontrol(1,0)
303  ELSE
304  CALL pushcontrol(1,1)
305  END IF
306  IF (ne_corner) THEN
307  CALL pushcontrol(1,0)
308  ELSE
309  CALL pushcontrol(1,1)
310  END IF
311  IF (grid_type .LT. 3) THEN
312  DO j=js-1,je+2
313  IF (j .EQ. 1 .AND. (.NOT.nested)) THEN
314  CALL pushcontrol(3,4)
315  ELSE IF ((j .EQ. 0 .OR. j .EQ. npy - 1) .AND. (.NOT.nested)) &
316 & THEN
317  CALL pushcontrol(3,3)
318  ELSE IF ((j .EQ. 2 .OR. j .EQ. npy + 1) .AND. (.NOT.nested)) &
319 & THEN
320  CALL pushcontrol(3,2)
321  ELSE IF (j .EQ. npy .AND. (.NOT.nested)) THEN
322  CALL pushcontrol(3,1)
323  ELSE
324  CALL pushcontrol(3,0)
325  END IF
326  END DO
327  CALL pushinteger(npt)
328  CALL pushinteger(ifirst)
329  CALL pushinteger(min6)
330  CALL pushinteger(min5)
331  CALL pushinteger(min3)
332  CALL pushinteger(min1)
333  CALL pushinteger(id)
334  CALL pushinteger(ilast)
335  CALL pushinteger(max6)
336  CALL pushinteger(max5)
337  CALL pushinteger(max3)
338  CALL pushinteger(max1)
339  CALL pushcontrol(1,0)
340  ELSE
341  CALL pushinteger(npt)
342  CALL pushinteger(ifirst)
343  CALL pushinteger(min6)
344  CALL pushinteger(min5)
345  CALL pushinteger(min3)
346  CALL pushinteger(min1)
347  CALL pushinteger(id)
348  CALL pushinteger(ilast)
349  CALL pushinteger(max6)
350  CALL pushinteger(max5)
351  CALL pushinteger(max3)
352  CALL pushinteger(max1)
353  CALL pushcontrol(1,1)
354  END IF
355  END SUBROUTINE d2c_setup_fwd
356 ! Differentiation of d2c_setup in reverse (adjoint) mode, backward sweep (with options split(a2b_edge_mod.a2b_ord4 a2b_edge_m
357 !od.a2b_ord2 dyn_core_mod.dyn_core dyn_core_mod.pk3_halo dyn_core_mod.pln_halo dyn_core_mod.pe_halo dyn_core_mod.adv_pe dyn_core_
358 !mod.p_grad_c dyn_core_mod.nh_p_grad dyn_core_mod.split_p_grad dyn_core_mod.one_grad_p dyn_core_mod.grad1_p_update dyn_core_mod.m
359 !ix_dp dyn_core_mod.geopk dyn_core_mod.del2_cubed dyn_core_mod.Rayleigh_fast fv_dynamics_mod.fv_dynamics fv_dynamics_mod.Rayleigh
360 !_Super fv_dynamics_mod.Rayleigh_Friction fv_dynamics_mod.compute_aam fv_grid_utils_mod.cubed_to_latlon fv_grid_utils_mod.c2l_ord
361 !4 fv_grid_utils_mod.c2l_ord2 fv_mapz_mod.Lagrangian_to_Eulerian fv_mapz_mod.compute_total_energy fv_mapz_mod.pkez fv_mapz_mod.re
362 !map_z fv_mapz_mod.map_scalar fv_mapz_mod.map1_ppm fv_mapz_mod.mapn_tracer fv_mapz_mod.map1_q2 fv_mapz_mod.remap_2d f
363 !v_mapz_mod.scalar_profile fv_mapz_mod.cs_profile fv_mapz_mod.cs_limiters fv_mapz_mod.ppm_profile fv_mapz_mod.ppm_limiters
364 !fv_mapz_mod.steepz fv_mapz_mod.rst_remap fv_mapz_mod.mappm fv_mapz_mod.moist_cv fv_mapz_mod.moist_cp fv_mapz_mod.map1_cubic fv_r
365 !estart_mod.d2c_setup fv_tracer2d_mod.tracer_2d_1L fv_tracer2d_mod.tracer_2d fv_tracer2d_mod.tracer_2d_nested fv_sg_mod.fv_subgri
366 !d_z main_mod.compute_pressures main_mod.run nh_core_mod.Riem_Solver3 nh_utils_mod.update_dz_c nh_utils_mod.update_dz_d nh_utils_
367 !mod.Riem_Solver_c nh_utils_mod.Riem_Solver3test nh_utils_mod.imp_diff_w nh_utils_mod.RIM_2D nh_utils_mod.SIM3_solver nh_utils_mo
368 !d.SIM3p0_solver nh_utils_mod.SIM1_solver nh_utils_mod.SIM_solver nh_utils_mod.edge_scalar nh_utils_mod.edge_profile nh_utils_mod
369 !.nest_halo_nh sw_core_mod.c_sw sw_core_mod.d_sw sw_core_mod.divergence_corner sw_core_mod.divergence_corner_nest sw_core_mod.d2
370 !a2c_vect sw_core_mod.fill3_4corners sw_core_mod.fill2_4corners sw_core_mod.fill_4corners sw_core_mod.xtp_u sw_core_mod.ytp_v_
371 !fb sw_core_mod.compute_divergence_damping sw_core_mod.smag_corner tp_core_mod.mp_ghost_ew tp_core_mod.fv_tp_2d tp_core_
372 !mod.copy_corners tp_core_mod.xppm tp_core_mod.yppm tp_core_mod.deln_flux a2b_edge_mod.extrap_corner fv_grid_utils
373 !_mod.great_circle_dist sw_core_mod.edge_interpolate4)):
374 ! gradient of useful results: u v uc vc
375 ! with respect to varying inputs: u v uc vc
376  SUBROUTINE d2c_setup_bwd(u, u_ad, v, v_ad, ua, va, uc, uc_ad, vc, &
377 & vc_ad, dord4, isd, ied, jsd, jed, is, ie, js, je, npx, npy, &
378 & grid_type, nested, se_corner, sw_corner, ne_corner, nw_corner, &
379 & rsin_u, rsin_v, cosa_s, rsin2)
380  IMPLICIT NONE
381  LOGICAL, INTENT(IN) :: dord4
382  INTEGER, INTENT(IN) :: isd, ied, jsd, jed, is, ie, js, je, npx, npy&
383 & , grid_type
384  REAL, INTENT(IN) :: u(isd:ied, jsd:jed+1)
385  REAL :: u_ad(isd:ied, jsd:jed+1)
386  REAL, INTENT(IN) :: v(isd:ied+1, jsd:jed)
387  REAL :: v_ad(isd:ied+1, jsd:jed)
388  REAL, DIMENSION(isd:ied, jsd:jed), INTENT(OUT) :: ua
389  REAL, DIMENSION(isd:ied, jsd:jed), INTENT(OUT) :: va
390  REAL, DIMENSION(isd:ied+1, jsd:jed) :: uc
391  REAL, DIMENSION(isd:ied+1, jsd:jed) :: uc_ad
392  REAL, DIMENSION(isd:ied, jsd:jed+1) :: vc
393  REAL, DIMENSION(isd:ied, jsd:jed+1) :: vc_ad
394  LOGICAL, INTENT(IN) :: nested, se_corner, sw_corner, ne_corner, &
395 & nw_corner
396  REAL, INTENT(IN) :: rsin_u(isd:ied+1, jsd:jed)
397  REAL, INTENT(IN) :: rsin_v(isd:ied, jsd:jed+1)
398  REAL, INTENT(IN) :: cosa_s(isd:ied, jsd:jed)
399  REAL, INTENT(IN) :: rsin2(isd:ied, jsd:jed)
400  REAL, DIMENSION(isd:ied, jsd:jed) :: utmp, vtmp
401  REAL, DIMENSION(isd:ied, jsd:jed) :: utmp_ad, vtmp_ad
402  REAL, PARAMETER :: t11=27./28., t12=-(13./28.), t13=3./7., t14=6./7.&
403 & , t15=3./28.
404  REAL, PARAMETER :: a1=0.5625
405  REAL, PARAMETER :: a2=-0.0625
406  REAL, PARAMETER :: c1=-(2./14.)
407  REAL, PARAMETER :: c2=11./14.
408  REAL, PARAMETER :: c3=5./14.
409  INTEGER :: npt, i, j, ifirst, ilast, id
410  INTRINSIC max
411  INTRINSIC min
412  INTEGER :: max1
413  INTEGER :: max2
414  INTEGER :: max3
415  INTEGER :: max4
416  INTEGER :: max5
417  INTEGER :: max6
418  INTEGER :: min1
419  INTEGER :: min2
420  INTEGER :: min3
421  INTEGER :: min4
422  INTEGER :: min5
423  INTEGER :: min6
424  REAL :: temp_ad
425  REAL :: temp_ad0
426  REAL :: temp_ad1
427  REAL :: temp_ad2
428  INTEGER :: ad_from
429  INTEGER :: ad_to
430  INTEGER :: ad_from0
431  INTEGER :: ad_to0
432  INTEGER :: branch
433 
434  utmp = 0.0
435  vtmp = 0.0
436  npt = 0
437  ifirst = 0
438  ilast = 0
439  id = 0
440  max1 = 0
441  max2 = 0
442  max3 = 0
443  max4 = 0
444  max5 = 0
445  max6 = 0
446  min1 = 0
447  min2 = 0
448  min3 = 0
449  min4 = 0
450  min5 = 0
451  min6 = 0
452  ad_from = 0
453  ad_to = 0
454  ad_from0 = 0
455  ad_to0 = 0
456  branch = 0
457 
458  CALL popcontrol(1,branch)
459  IF (branch .EQ. 0) THEN
460  CALL popinteger(max1)
461  CALL popinteger(max3)
462  CALL popinteger(max5)
463  CALL popinteger(max6)
464  CALL popinteger(ilast)
465  CALL popinteger(id)
466  CALL popinteger(min1)
467  CALL popinteger(min3)
468  CALL popinteger(min5)
469  CALL popinteger(min6)
470  CALL popinteger(ifirst)
471  CALL popinteger(npt)
472  vtmp_ad = 0.0
473  DO j=je+2,js-1,-1
474  CALL popcontrol(3,branch)
475  IF (branch .LT. 2) THEN
476  IF (branch .EQ. 0) THEN
477  DO i=ie+1,is-1,-1
478  vtmp_ad(i, j-2) = vtmp_ad(i, j-2) + a2*vc_ad(i, j)
479  vtmp_ad(i, j+1) = vtmp_ad(i, j+1) + a2*vc_ad(i, j)
480  vtmp_ad(i, j-1) = vtmp_ad(i, j-1) + a1*vc_ad(i, j)
481  vtmp_ad(i, j) = vtmp_ad(i, j) + a1*vc_ad(i, j)
482  vc_ad(i, j) = 0.0
483  END DO
484  ELSE
485  DO i=ie+1,is-1,-1
486  temp_ad2 = rsin_v(i, npy)*vc_ad(i, npy)
487  vtmp_ad(i, npy-1) = vtmp_ad(i, npy-1) + t14*temp_ad2
488  vtmp_ad(i, npy) = vtmp_ad(i, npy) + t14*temp_ad2
489  vtmp_ad(i, npy-2) = vtmp_ad(i, npy-2) + t12*temp_ad2
490  vtmp_ad(i, npy+1) = vtmp_ad(i, npy+1) + t12*temp_ad2
491  vtmp_ad(i, npy-3) = vtmp_ad(i, npy-3) + t15*temp_ad2
492  vtmp_ad(i, npy+2) = vtmp_ad(i, npy+2) + t15*temp_ad2
493  vc_ad(i, npy) = 0.0
494  END DO
495  END IF
496  ELSE IF (branch .EQ. 2) THEN
497  DO i=ie+1,is-1,-1
498  vtmp_ad(i, j+1) = vtmp_ad(i, j+1) + c1*vc_ad(i, j)
499  vtmp_ad(i, j) = vtmp_ad(i, j) + c2*vc_ad(i, j)
500  vtmp_ad(i, j-1) = vtmp_ad(i, j-1) + c3*vc_ad(i, j)
501  vc_ad(i, j) = 0.0
502  END DO
503  ELSE IF (branch .EQ. 3) THEN
504  DO i=ie+1,is-1,-1
505  vtmp_ad(i, j-2) = vtmp_ad(i, j-2) + c1*vc_ad(i, j)
506  vtmp_ad(i, j-1) = vtmp_ad(i, j-1) + c2*vc_ad(i, j)
507  vtmp_ad(i, j) = vtmp_ad(i, j) + c3*vc_ad(i, j)
508  vc_ad(i, j) = 0.0
509  END DO
510  ELSE
511  DO i=ie+1,is-1,-1
512  temp_ad1 = rsin_v(i, 1)*vc_ad(i, 1)
513  vtmp_ad(i, 0) = vtmp_ad(i, 0) + t14*temp_ad1
514  vtmp_ad(i, 1) = vtmp_ad(i, 1) + t14*temp_ad1
515  vtmp_ad(i, -1) = vtmp_ad(i, -1) + t12*temp_ad1
516  vtmp_ad(i, 2) = vtmp_ad(i, 2) + t12*temp_ad1
517  vtmp_ad(i, -2) = vtmp_ad(i, -2) + t15*temp_ad1
518  vtmp_ad(i, 3) = vtmp_ad(i, 3) + t15*temp_ad1
519  vc_ad(i, 1) = 0.0
520  END DO
521  END IF
522  END DO
523  ELSE
524  CALL popinteger(max1)
525  CALL popinteger(max3)
526  CALL popinteger(max5)
527  CALL popinteger(max6)
528  CALL popinteger(ilast)
529  CALL popinteger(id)
530  CALL popinteger(min1)
531  CALL popinteger(min3)
532  CALL popinteger(min5)
533  CALL popinteger(min6)
534  CALL popinteger(ifirst)
535  CALL popinteger(npt)
536  vtmp_ad = 0.0
537  DO j=je+2,js-1,-1
538  DO i=ie+1,is-1,-1
539  vtmp_ad(i, j-2) = vtmp_ad(i, j-2) + a2*vc_ad(i, j)
540  vtmp_ad(i, j+1) = vtmp_ad(i, j+1) + a2*vc_ad(i, j)
541  vtmp_ad(i, j-1) = vtmp_ad(i, j-1) + a1*vc_ad(i, j)
542  vtmp_ad(i, j) = vtmp_ad(i, j) + a1*vc_ad(i, j)
543  vc_ad(i, j) = 0.0
544  END DO
545  END DO
546  END IF
547  CALL popcontrol(1,branch)
548  IF (branch .EQ. 0) THEN
549  utmp_ad = 0.0
550  DO j=2,0,-1
551  utmp_ad(ie-j, npy) = utmp_ad(ie-j, npy) - vtmp_ad(npx, npy+j)
552  vtmp_ad(npx, npy+j) = 0.0
553  END DO
554  ELSE
555  utmp_ad = 0.0
556  END IF
557  CALL popcontrol(1,branch)
558  IF (branch .EQ. 0) THEN
559  DO j=0,-2,-1
560  utmp_ad(ie+j, 0) = utmp_ad(ie+j, 0) + vtmp_ad(npx, j)
561  vtmp_ad(npx, j) = 0.0
562  END DO
563  END IF
564  CALL popcontrol(1,branch)
565  IF (branch .EQ. 0) THEN
566  DO j=2,0,-1
567  utmp_ad(j+1, npy) = utmp_ad(j+1, npy) + vtmp_ad(0, npy+j)
568  vtmp_ad(0, npy+j) = 0.0
569  END DO
570  END IF
571  CALL popcontrol(1,branch)
572  IF (branch .EQ. 0) THEN
573  DO j=0,-2,-1
574  utmp_ad(1-j, 0) = utmp_ad(1-j, 0) - vtmp_ad(0, j)
575  vtmp_ad(0, j) = 0.0
576  END DO
577  END IF
578  CALL popcontrol(2,branch)
579  IF (branch .EQ. 0) THEN
580  DO j=je+1,js-1,-1
581  utmp_ad(npx, j) = utmp_ad(npx, j) + c3*uc_ad(npx+1, j)
582  utmp_ad(npx+1, j) = utmp_ad(npx+1, j) + c2*uc_ad(npx+1, j)
583  utmp_ad(npx+2, j) = utmp_ad(npx+2, j) + c1*uc_ad(npx+1, j)
584  uc_ad(npx+1, j) = 0.0
585  temp_ad0 = rsin_u(npx, j)*uc_ad(npx, j)
586  utmp_ad(npx-1, j) = utmp_ad(npx-1, j) + t14*temp_ad0
587  utmp_ad(npx, j) = utmp_ad(npx, j) + t14*temp_ad0
588  utmp_ad(npx-2, j) = utmp_ad(npx-2, j) + t12*temp_ad0
589  utmp_ad(npx+1, j) = utmp_ad(npx+1, j) + t12*temp_ad0
590  utmp_ad(npx-3, j) = utmp_ad(npx-3, j) + t15*temp_ad0
591  utmp_ad(npx+2, j) = utmp_ad(npx+2, j) + t15*temp_ad0
592  uc_ad(npx, j) = 0.0
593  utmp_ad(npx-3, j) = utmp_ad(npx-3, j) + c1*uc_ad(npx-1, j)
594  utmp_ad(npx-2, j) = utmp_ad(npx-2, j) + c2*uc_ad(npx-1, j)
595  utmp_ad(npx-1, j) = utmp_ad(npx-1, j) + c3*uc_ad(npx-1, j)
596  uc_ad(npx-1, j) = 0.0
597  END DO
598  ELSE IF (branch .NE. 1) THEN
599  GOTO 100
600  END IF
601  CALL popcontrol(1,branch)
602  IF (branch .EQ. 0) THEN
603  DO j=je+1,js-1,-1
604  utmp_ad(3, j) = utmp_ad(3, j) + c1*uc_ad(2, j)
605  utmp_ad(2, j) = utmp_ad(2, j) + c2*uc_ad(2, j)
606  utmp_ad(1, j) = utmp_ad(1, j) + c3*uc_ad(2, j)
607  uc_ad(2, j) = 0.0
608  temp_ad = rsin_u(1, j)*uc_ad(1, j)
609  utmp_ad(0, j) = utmp_ad(0, j) + t14*temp_ad
610  utmp_ad(1, j) = utmp_ad(1, j) + t14*temp_ad
611  utmp_ad(-1, j) = utmp_ad(-1, j) + t12*temp_ad
612  utmp_ad(2, j) = utmp_ad(2, j) + t12*temp_ad
613  utmp_ad(-2, j) = utmp_ad(-2, j) + t15*temp_ad
614  utmp_ad(3, j) = utmp_ad(3, j) + t15*temp_ad
615  uc_ad(1, j) = 0.0
616  utmp_ad(-2, j) = utmp_ad(-2, j) + c1*uc_ad(0, j)
617  utmp_ad(-1, j) = utmp_ad(-1, j) + c2*uc_ad(0, j)
618  utmp_ad(0, j) = utmp_ad(0, j) + c3*uc_ad(0, j)
619  uc_ad(0, j) = 0.0
620  END DO
621  END IF
622  100 DO j=je+1,js-1,-1
623  DO i=ilast,ifirst,-1
624  utmp_ad(i-1, j) = utmp_ad(i-1, j) + a1*uc_ad(i, j)
625  utmp_ad(i, j) = utmp_ad(i, j) + a1*uc_ad(i, j)
626  utmp_ad(i-2, j) = utmp_ad(i-2, j) + a2*uc_ad(i, j)
627  utmp_ad(i+1, j) = utmp_ad(i+1, j) + a2*uc_ad(i, j)
628  uc_ad(i, j) = 0.0
629  END DO
630  END DO
631  CALL popcontrol(1,branch)
632  CALL popcontrol(1,branch)
633  IF (branch .EQ. 0) THEN
634  DO i=0,-2,-1
635  vtmp_ad(0, je+i) = vtmp_ad(0, je+i) + utmp_ad(i, npy)
636  utmp_ad(i, npy) = 0.0
637  END DO
638  END IF
639  CALL popcontrol(1,branch)
640  IF (branch .EQ. 0) THEN
641  DO i=2,0,-1
642  vtmp_ad(npx, je-i) = vtmp_ad(npx, je-i) - utmp_ad(npx+i, npy)
643  utmp_ad(npx+i, npy) = 0.0
644  END DO
645  END IF
646  CALL popcontrol(1,branch)
647  IF (branch .EQ. 0) THEN
648  DO i=2,0,-1
649  vtmp_ad(npx, i+1) = vtmp_ad(npx, i+1) + utmp_ad(npx+i, 0)
650  utmp_ad(npx+i, 0) = 0.0
651  END DO
652  END IF
653  CALL popcontrol(1,branch)
654  IF (branch .EQ. 0) THEN
655  DO i=0,-2,-1
656  vtmp_ad(0, 1-i) = vtmp_ad(0, 1-i) - utmp_ad(i, 0)
657  utmp_ad(i, 0) = 0.0
658  END DO
659  END IF
660  CALL popcontrol(2,branch)
661  IF (branch .LT. 2) THEN
662  IF (branch .EQ. 0) THEN
663  DO j=jed,jsd,-1
664  v_ad(ied, j) = v_ad(ied, j) + 0.5*vtmp_ad(ied, j)
665  v_ad(ied+1, j) = v_ad(ied+1, j) + 0.5*vtmp_ad(ied, j)
666  vtmp_ad(ied, j) = 0.0
667  v_ad(isd, j) = v_ad(isd, j) + 0.5*vtmp_ad(isd, j)
668  v_ad(isd+1, j) = v_ad(isd+1, j) + 0.5*vtmp_ad(isd, j)
669  vtmp_ad(isd, j) = 0.0
670  DO i=ied-1,isd+1,-1
671  v_ad(i-1, j) = v_ad(i-1, j) + a2*vtmp_ad(i, j)
672  v_ad(i+2, j) = v_ad(i+2, j) + a2*vtmp_ad(i, j)
673  v_ad(i, j) = v_ad(i, j) + a1*vtmp_ad(i, j)
674  v_ad(i+1, j) = v_ad(i+1, j) + a1*vtmp_ad(i, j)
675  vtmp_ad(i, j) = 0.0
676  END DO
677  END DO
678  DO i=ied,isd,-1
679  u_ad(i, jed) = u_ad(i, jed) + 0.5*utmp_ad(i, jed)
680  u_ad(i, jed+1) = u_ad(i, jed+1) + 0.5*utmp_ad(i, jed)
681  utmp_ad(i, jed) = 0.0
682  u_ad(i, jsd) = u_ad(i, jsd) + 0.5*utmp_ad(i, jsd)
683  u_ad(i, jsd+1) = u_ad(i, jsd+1) + 0.5*utmp_ad(i, jsd)
684  utmp_ad(i, jsd) = 0.0
685  END DO
686  DO j=jed-1,jsd+1,-1
687  DO i=ied,isd,-1
688  u_ad(i, j-1) = u_ad(i, j-1) + a2*utmp_ad(i, j)
689  u_ad(i, j+2) = u_ad(i, j+2) + a2*utmp_ad(i, j)
690  u_ad(i, j) = u_ad(i, j) + a1*utmp_ad(i, j)
691  u_ad(i, j+1) = u_ad(i, j+1) + a1*utmp_ad(i, j)
692  utmp_ad(i, j) = 0.0
693  END DO
694  END DO
695  GOTO 110
696  END IF
697  ELSE
698  IF (branch .NE. 2) THEN
699  DO j=min6,max6,-1
700  DO i=ied,npx-npt+1,-1
701  v_ad(i, j) = v_ad(i, j) + 0.5*vtmp_ad(i, j)
702  v_ad(i+1, j) = v_ad(i+1, j) + 0.5*vtmp_ad(i, j)
703  vtmp_ad(i, j) = 0.0
704  u_ad(i, j) = u_ad(i, j) + 0.5*utmp_ad(i, j)
705  u_ad(i, j+1) = u_ad(i, j+1) + 0.5*utmp_ad(i, j)
706  utmp_ad(i, j) = 0.0
707  END DO
708  END DO
709  END IF
710  CALL popcontrol(1,branch)
711  IF (branch .EQ. 0) THEN
712  DO j=min5,max5,-1
713  DO i=npt-1,isd,-1
714  v_ad(i, j) = v_ad(i, j) + 0.5*vtmp_ad(i, j)
715  v_ad(i+1, j) = v_ad(i+1, j) + 0.5*vtmp_ad(i, j)
716  vtmp_ad(i, j) = 0.0
717  u_ad(i, j) = u_ad(i, j) + 0.5*utmp_ad(i, j)
718  u_ad(i, j+1) = u_ad(i, j+1) + 0.5*utmp_ad(i, j)
719  utmp_ad(i, j) = 0.0
720  END DO
721  END DO
722  END IF
723  CALL popcontrol(1,branch)
724  IF (branch .EQ. 0) THEN
725  DO j=jed,npy-npt+1,-1
726  DO i=ied,isd,-1
727  v_ad(i, j) = v_ad(i, j) + 0.5*vtmp_ad(i, j)
728  v_ad(i+1, j) = v_ad(i+1, j) + 0.5*vtmp_ad(i, j)
729  vtmp_ad(i, j) = 0.0
730  u_ad(i, j) = u_ad(i, j) + 0.5*utmp_ad(i, j)
731  u_ad(i, j+1) = u_ad(i, j+1) + 0.5*utmp_ad(i, j)
732  utmp_ad(i, j) = 0.0
733  END DO
734  END DO
735  END IF
736  CALL popcontrol(1,branch)
737  IF (branch .EQ. 0) THEN
738  DO j=npt-1,jsd,-1
739  DO i=ied,isd,-1
740  v_ad(i, j) = v_ad(i, j) + 0.5*vtmp_ad(i, j)
741  v_ad(i+1, j) = v_ad(i+1, j) + 0.5*vtmp_ad(i, j)
742  vtmp_ad(i, j) = 0.0
743  u_ad(i, j) = u_ad(i, j) + 0.5*utmp_ad(i, j)
744  u_ad(i, j+1) = u_ad(i, j+1) + 0.5*utmp_ad(i, j)
745  utmp_ad(i, j) = 0.0
746  END DO
747  END DO
748  END IF
749  END IF
750  DO j=min3,max3,-1
751  CALL popinteger(ad_from0)
752  CALL popinteger(ad_to0)
753  DO i=ad_to0,ad_from0,-1
754  v_ad(i-1, j) = v_ad(i-1, j) + a2*vtmp_ad(i, j)
755  v_ad(i+2, j) = v_ad(i+2, j) + a2*vtmp_ad(i, j)
756  v_ad(i, j) = v_ad(i, j) + a1*vtmp_ad(i, j)
757  v_ad(i+1, j) = v_ad(i+1, j) + a1*vtmp_ad(i, j)
758  vtmp_ad(i, j) = 0.0
759  END DO
760  END DO
761  DO j=min1,max1,-1
762  CALL popinteger(ad_from)
763  CALL popinteger(ad_to)
764  DO i=ad_to,ad_from,-1
765  u_ad(i, j-1) = u_ad(i, j-1) + a2*utmp_ad(i, j)
766  u_ad(i, j+2) = u_ad(i, j+2) + a2*utmp_ad(i, j)
767  u_ad(i, j) = u_ad(i, j) + a1*utmp_ad(i, j)
768  u_ad(i, j+1) = u_ad(i, j+1) + a1*utmp_ad(i, j)
769  utmp_ad(i, j) = 0.0
770  END DO
771  END DO
772  110 CONTINUE
773  END SUBROUTINE d2c_setup_bwd
774  SUBROUTINE d2c_setup(u, v, ua, va, uc, vc, dord4, isd, ied, jsd, jed, &
775 & is, ie, js, je, npx, npy, grid_type, nested, se_corner, sw_corner, &
776 & ne_corner, nw_corner, rsin_u, rsin_v, cosa_s, rsin2)
777  IMPLICIT NONE
778  LOGICAL, INTENT(IN) :: dord4
779  INTEGER, INTENT(IN) :: isd, ied, jsd, jed, is, ie, js, je, npx, npy&
780 & , grid_type
781  REAL, INTENT(IN) :: u(isd:ied, jsd:jed+1)
782  REAL, INTENT(IN) :: v(isd:ied+1, jsd:jed)
783  REAL, DIMENSION(isd:ied, jsd:jed), INTENT(OUT) :: ua
784  REAL, DIMENSION(isd:ied, jsd:jed), INTENT(OUT) :: va
785  REAL, DIMENSION(isd:ied+1, jsd:jed), INTENT(OUT) :: uc
786  REAL, DIMENSION(isd:ied, jsd:jed+1), INTENT(OUT) :: vc
787  LOGICAL, INTENT(IN) :: nested, se_corner, sw_corner, ne_corner, &
788 & nw_corner
789  REAL, INTENT(IN) :: rsin_u(isd:ied+1, jsd:jed)
790  REAL, INTENT(IN) :: rsin_v(isd:ied, jsd:jed+1)
791  REAL, INTENT(IN) :: cosa_s(isd:ied, jsd:jed)
792  REAL, INTENT(IN) :: rsin2(isd:ied, jsd:jed)
793 ! Local
794  REAL, DIMENSION(isd:ied, jsd:jed) :: utmp, vtmp
795  REAL, PARAMETER :: t11=27./28., t12=-(13./28.), t13=3./7., t14=6./7.&
796 & , t15=3./28.
797  REAL, PARAMETER :: a1=0.5625
798  REAL, PARAMETER :: a2=-0.0625
799  REAL, PARAMETER :: c1=-(2./14.)
800  REAL, PARAMETER :: c2=11./14.
801  REAL, PARAMETER :: c3=5./14.
802  INTEGER :: npt, i, j, ifirst, ilast, id
803  INTRINSIC max
804  INTRINSIC min
805  INTEGER :: max1
806  INTEGER :: max2
807  INTEGER :: max3
808  INTEGER :: max4
809  INTEGER :: max5
810  INTEGER :: max6
811  INTEGER :: min1
812  INTEGER :: min2
813  INTEGER :: min3
814  INTEGER :: min4
815  INTEGER :: min5
816  INTEGER :: min6
817  IF (dord4) THEN
818  id = 1
819  ELSE
820  id = 0
821  END IF
822  IF (grid_type .LT. 3 .AND. (.NOT.nested)) THEN
823  npt = 4
824  ELSE
825  npt = -2
826  END IF
827  IF (nested) THEN
828  DO j=jsd+1,jed-1
829  DO i=isd,ied
830  utmp(i, j) = a2*(u(i, j-1)+u(i, j+2)) + a1*(u(i, j)+u(i, j+1))
831  END DO
832  END DO
833  DO i=isd,ied
834 !j = jsd
835  utmp(i, jsd) = 0.5*(u(i, jsd)+u(i, jsd+1))
836 !j = jed
837  utmp(i, jed) = 0.5*(u(i, jed)+u(i, jed+1))
838  END DO
839  DO j=jsd,jed
840  DO i=isd+1,ied-1
841  vtmp(i, j) = a2*(v(i-1, j)+v(i+2, j)) + a1*(v(i, j)+v(i+1, j))
842  END DO
843 !i = isd
844  vtmp(isd, j) = 0.5*(v(isd, j)+v(isd+1, j))
845 !i = ied
846  vtmp(ied, j) = 0.5*(v(ied, j)+v(ied+1, j))
847  END DO
848  DO j=jsd,jed
849  DO i=isd,ied
850  ua(i, j) = (utmp(i, j)-vtmp(i, j)*cosa_s(i, j))*rsin2(i, j)
851  va(i, j) = (vtmp(i, j)-utmp(i, j)*cosa_s(i, j))*rsin2(i, j)
852  END DO
853  END DO
854  ELSE
855 !----------
856 ! Interior:
857 !----------
858  utmp = 0.
859  vtmp = 0.
860  IF (npt .LT. js - 1) THEN
861  max1 = js - 1
862  ELSE
863  max1 = npt
864  END IF
865  IF (npy - npt .GT. je + 1) THEN
866  min1 = je + 1
867  ELSE
868  min1 = npy - npt
869  END IF
870  DO j=max1,min1
871  IF (npt .LT. isd) THEN
872  max2 = isd
873  ELSE
874  max2 = npt
875  END IF
876  IF (npx - npt .GT. ied) THEN
877  min2 = ied
878  ELSE
879  min2 = npx - npt
880  END IF
881  DO i=max2,min2
882  utmp(i, j) = a2*(u(i, j-1)+u(i, j+2)) + a1*(u(i, j)+u(i, j+1))
883  END DO
884  END DO
885  IF (npt .LT. jsd) THEN
886  max3 = jsd
887  ELSE
888  max3 = npt
889  END IF
890  IF (npy - npt .GT. jed) THEN
891  min3 = jed
892  ELSE
893  min3 = npy - npt
894  END IF
895  DO j=max3,min3
896  IF (npt .LT. is - 1) THEN
897  max4 = is - 1
898  ELSE
899  max4 = npt
900  END IF
901  IF (npx - npt .GT. ie + 1) THEN
902  min4 = ie + 1
903  ELSE
904  min4 = npx - npt
905  END IF
906  DO i=max4,min4
907  vtmp(i, j) = a2*(v(i-1, j)+v(i+2, j)) + a1*(v(i, j)+v(i+1, j))
908  END DO
909  END DO
910 !----------
911 ! edges:
912 !----------
913  IF (grid_type .LT. 3) THEN
914  IF (js .EQ. 1 .OR. jsd .LT. npt) THEN
915  DO j=jsd,npt-1
916  DO i=isd,ied
917  utmp(i, j) = 0.5*(u(i, j)+u(i, j+1))
918  vtmp(i, j) = 0.5*(v(i, j)+v(i+1, j))
919  END DO
920  END DO
921  END IF
922  IF (je + 1 .EQ. npy .OR. jed .GE. npy - npt) THEN
923  DO j=npy-npt+1,jed
924  DO i=isd,ied
925  utmp(i, j) = 0.5*(u(i, j)+u(i, j+1))
926  vtmp(i, j) = 0.5*(v(i, j)+v(i+1, j))
927  END DO
928  END DO
929  END IF
930  IF (is .EQ. 1 .OR. isd .LT. npt) THEN
931  IF (npt .LT. jsd) THEN
932  max5 = jsd
933  ELSE
934  max5 = npt
935  END IF
936  IF (npy - npt .GT. jed) THEN
937  min5 = jed
938  ELSE
939  min5 = npy - npt
940  END IF
941  DO j=max5,min5
942  DO i=isd,npt-1
943  utmp(i, j) = 0.5*(u(i, j)+u(i, j+1))
944  vtmp(i, j) = 0.5*(v(i, j)+v(i+1, j))
945  END DO
946  END DO
947  END IF
948  IF (ie + 1 .EQ. npx .OR. ied .GE. npx - npt) THEN
949  IF (npt .LT. jsd) THEN
950  max6 = jsd
951  ELSE
952  max6 = npt
953  END IF
954  IF (npy - npt .GT. jed) THEN
955  min6 = jed
956  ELSE
957  min6 = npy - npt
958  END IF
959  DO j=max6,min6
960  DO i=npx-npt+1,ied
961  utmp(i, j) = 0.5*(u(i, j)+u(i, j+1))
962  vtmp(i, j) = 0.5*(v(i, j)+v(i+1, j))
963  END DO
964  END DO
965  END IF
966  END IF
967  DO j=js-1-id,je+1+id
968  DO i=is-1-id,ie+1+id
969  ua(i, j) = (utmp(i, j)-vtmp(i, j)*cosa_s(i, j))*rsin2(i, j)
970  va(i, j) = (vtmp(i, j)-utmp(i, j)*cosa_s(i, j))*rsin2(i, j)
971  END DO
972  END DO
973  END IF
974 ! A -> C
975 !--------------
976 ! Fix the edges
977 !--------------
978 ! Xdir:
979  IF (sw_corner) THEN
980  DO i=-2,0
981  utmp(i, 0) = -vtmp(0, 1-i)
982  END DO
983  END IF
984  IF (se_corner) THEN
985  DO i=0,2
986  utmp(npx+i, 0) = vtmp(npx, i+1)
987  END DO
988  END IF
989  IF (ne_corner) THEN
990  DO i=0,2
991  utmp(npx+i, npy) = -vtmp(npx, je-i)
992  END DO
993  END IF
994  IF (nw_corner) THEN
995  DO i=-2,0
996  utmp(i, npy) = vtmp(0, je+i)
997  END DO
998  END IF
999  IF (grid_type .LT. 3 .AND. (.NOT.nested)) THEN
1000  IF (3 .LT. is - 1) THEN
1001  ifirst = is - 1
1002  ELSE
1003  ifirst = 3
1004  END IF
1005  IF (npx - 2 .GT. ie + 2) THEN
1006  ilast = ie + 2
1007  ELSE
1008  ilast = npx - 2
1009  END IF
1010  ELSE
1011  ifirst = is - 1
1012  ilast = ie + 2
1013  END IF
1014 !---------------------------------------------
1015 ! 4th order interpolation for interior points:
1016 !---------------------------------------------
1017  DO j=js-1,je+1
1018  DO i=ifirst,ilast
1019  uc(i, j) = a1*(utmp(i-1, j)+utmp(i, j)) + a2*(utmp(i-2, j)+utmp(&
1020 & i+1, j))
1021  END DO
1022  END DO
1023  IF (grid_type .LT. 3) THEN
1024 ! Xdir:
1025  IF (is .EQ. 1 .AND. (.NOT.nested)) THEN
1026  DO j=js-1,je+1
1027  uc(0, j) = c1*utmp(-2, j) + c2*utmp(-1, j) + c3*utmp(0, j)
1028  uc(1, j) = (t14*(utmp(0, j)+utmp(1, j))+t12*(utmp(-1, j)+utmp(&
1029 & 2, j))+t15*(utmp(-2, j)+utmp(3, j)))*rsin_u(1, j)
1030  uc(2, j) = c1*utmp(3, j) + c2*utmp(2, j) + c3*utmp(1, j)
1031  END DO
1032  END IF
1033  IF (ie + 1 .EQ. npx .AND. (.NOT.nested)) THEN
1034  DO j=js-1,je+1
1035  uc(npx-1, j) = c1*utmp(npx-3, j) + c2*utmp(npx-2, j) + c3*utmp&
1036 & (npx-1, j)
1037  uc(npx, j) = (t14*(utmp(npx-1, j)+utmp(npx, j))+t12*(utmp(npx-&
1038 & 2, j)+utmp(npx+1, j))+t15*(utmp(npx-3, j)+utmp(npx+2, j)))*&
1039 & rsin_u(npx, j)
1040  uc(npx+1, j) = c3*utmp(npx, j) + c2*utmp(npx+1, j) + c1*utmp(&
1041 & npx+2, j)
1042  END DO
1043  END IF
1044  END IF
1045 !------
1046 ! Ydir:
1047 !------
1048  IF (sw_corner) THEN
1049  DO j=-2,0
1050  vtmp(0, j) = -utmp(1-j, 0)
1051  END DO
1052  END IF
1053  IF (nw_corner) THEN
1054  DO j=0,2
1055  vtmp(0, npy+j) = utmp(j+1, npy)
1056  END DO
1057  END IF
1058  IF (se_corner) THEN
1059  DO j=-2,0
1060  vtmp(npx, j) = utmp(ie+j, 0)
1061  END DO
1062  END IF
1063  IF (ne_corner) THEN
1064  DO j=0,2
1065  vtmp(npx, npy+j) = -utmp(ie-j, npy)
1066  END DO
1067  END IF
1068  IF (grid_type .LT. 3) THEN
1069  DO j=js-1,je+2
1070  IF (j .EQ. 1 .AND. (.NOT.nested)) THEN
1071  DO i=is-1,ie+1
1072  vc(i, 1) = (t14*(vtmp(i, 0)+vtmp(i, 1))+t12*(vtmp(i, -1)+&
1073 & vtmp(i, 2))+t15*(vtmp(i, -2)+vtmp(i, 3)))*rsin_v(i, 1)
1074  END DO
1075  ELSE IF ((j .EQ. 0 .OR. j .EQ. npy - 1) .AND. (.NOT.nested)) &
1076 & THEN
1077  DO i=is-1,ie+1
1078  vc(i, j) = c1*vtmp(i, j-2) + c2*vtmp(i, j-1) + c3*vtmp(i, j)
1079  END DO
1080  ELSE IF ((j .EQ. 2 .OR. j .EQ. npy + 1) .AND. (.NOT.nested)) &
1081 & THEN
1082  DO i=is-1,ie+1
1083  vc(i, j) = c1*vtmp(i, j+1) + c2*vtmp(i, j) + c3*vtmp(i, j-1)
1084  END DO
1085  ELSE IF (j .EQ. npy .AND. (.NOT.nested)) THEN
1086  DO i=is-1,ie+1
1087  vc(i, npy) = (t14*(vtmp(i, npy-1)+vtmp(i, npy))+t12*(vtmp(i&
1088 & , npy-2)+vtmp(i, npy+1))+t15*(vtmp(i, npy-3)+vtmp(i, npy+2&
1089 & )))*rsin_v(i, npy)
1090  END DO
1091  ELSE
1092 ! 4th order interpolation for interior points:
1093  DO i=is-1,ie+1
1094  vc(i, j) = a2*(vtmp(i, j-2)+vtmp(i, j+1)) + a1*(vtmp(i, j-1)&
1095 & +vtmp(i, j))
1096  END DO
1097  END IF
1098  END DO
1099  ELSE
1100 ! 4th order interpolation:
1101  DO j=js-1,je+2
1102  DO i=is-1,ie+1
1103  vc(i, j) = a2*(vtmp(i, j-2)+vtmp(i, j+1)) + a1*(vtmp(i, j-1)+&
1104 & vtmp(i, j))
1105  END DO
1106  END DO
1107  END IF
1108  END SUBROUTINE d2c_setup
1109  SUBROUTINE d2a_setup(u, v, ua, va, dord4, isd, ied, jsd, jed, is, ie, &
1110 & js, je, npx, npy, grid_type, nested, cosa_s, rsin2)
1111  IMPLICIT NONE
1112  LOGICAL, INTENT(IN) :: dord4
1113  INTEGER, INTENT(IN) :: isd, ied, jsd, jed, is, ie, js, je, npx, npy&
1114 & , grid_type
1115  REAL, INTENT(IN) :: u(isd:ied, jsd:jed+1)
1116  REAL, INTENT(IN) :: v(isd:ied+1, jsd:jed)
1117  REAL, DIMENSION(isd:ied, jsd:jed), INTENT(OUT) :: ua
1118  REAL, DIMENSION(isd:ied, jsd:jed), INTENT(OUT) :: va
1119  REAL, INTENT(IN) :: cosa_s(isd:ied, jsd:jed)
1120  REAL, INTENT(IN) :: rsin2(isd:ied, jsd:jed)
1121  LOGICAL, INTENT(IN) :: nested
1122 ! Local
1123  REAL, DIMENSION(isd:ied, jsd:jed) :: utmp, vtmp
1124  REAL, PARAMETER :: t11=27./28., t12=-(13./28.), t13=3./7., t14=6./7.&
1125 & , t15=3./28.
1126  REAL, PARAMETER :: a1=0.5625
1127  REAL, PARAMETER :: a2=-0.0625
1128  REAL, PARAMETER :: c1=-(2./14.)
1129  REAL, PARAMETER :: c2=11./14.
1130  REAL, PARAMETER :: c3=5./14.
1131  INTEGER :: npt, i, j, ifirst, ilast, id
1132  INTRINSIC max
1133  INTRINSIC min
1134  INTEGER :: max1
1135  INTEGER :: max2
1136  INTEGER :: max3
1137  INTEGER :: max4
1138  INTEGER :: max5
1139  INTEGER :: max6
1140  INTEGER :: min1
1141  INTEGER :: min2
1142  INTEGER :: min3
1143  INTEGER :: min4
1144  INTEGER :: min5
1145  INTEGER :: min6
1146  IF (dord4) THEN
1147  id = 1
1148  ELSE
1149  id = 0
1150  END IF
1151  IF (grid_type .LT. 3 .AND. (.NOT.nested)) THEN
1152  npt = 4
1153  ELSE
1154  npt = -2
1155  END IF
1156  IF (nested) THEN
1157  DO j=jsd+1,jed-1
1158  DO i=isd,ied
1159  utmp(i, j) = a2*(u(i, j-1)+u(i, j+2)) + a1*(u(i, j)+u(i, j+1))
1160  END DO
1161  END DO
1162  DO i=isd,ied
1163 !j = jsd
1164  utmp(i, jsd) = 0.5*(u(i, jsd)+u(i, jsd+1))
1165 !j = jed
1166  utmp(i, jed) = 0.5*(u(i, jed)+u(i, jed+1))
1167  END DO
1168  DO j=jsd,jed
1169  DO i=isd+1,ied-1
1170  vtmp(i, j) = a2*(v(i-1, j)+v(i+2, j)) + a1*(v(i, j)+v(i+1, j))
1171  END DO
1172 !i = isd
1173  vtmp(isd, j) = 0.5*(v(isd, j)+v(isd+1, j))
1174 !i = ied
1175  vtmp(ied, j) = 0.5*(v(ied, j)+v(ied+1, j))
1176  END DO
1177  ELSE
1178  IF (npt .LT. js - 1) THEN
1179  max1 = js - 1
1180  ELSE
1181  max1 = npt
1182  END IF
1183  IF (npy - npt .GT. je + 1) THEN
1184  min1 = je + 1
1185  ELSE
1186  min1 = npy - npt
1187  END IF
1188 !----------
1189 ! Interior:
1190 !----------
1191  DO j=max1,min1
1192  IF (npt .LT. isd) THEN
1193  max2 = isd
1194  ELSE
1195  max2 = npt
1196  END IF
1197  IF (npx - npt .GT. ied) THEN
1198  min2 = ied
1199  ELSE
1200  min2 = npx - npt
1201  END IF
1202  DO i=max2,min2
1203  utmp(i, j) = a2*(u(i, j-1)+u(i, j+2)) + a1*(u(i, j)+u(i, j+1))
1204  END DO
1205  END DO
1206  IF (npt .LT. jsd) THEN
1207  max3 = jsd
1208  ELSE
1209  max3 = npt
1210  END IF
1211  IF (npy - npt .GT. jed) THEN
1212  min3 = jed
1213  ELSE
1214  min3 = npy - npt
1215  END IF
1216  DO j=max3,min3
1217  IF (npt .LT. is - 1) THEN
1218  max4 = is - 1
1219  ELSE
1220  max4 = npt
1221  END IF
1222  IF (npx - npt .GT. ie + 1) THEN
1223  min4 = ie + 1
1224  ELSE
1225  min4 = npx - npt
1226  END IF
1227  DO i=max4,min4
1228  vtmp(i, j) = a2*(v(i-1, j)+v(i+2, j)) + a1*(v(i, j)+v(i+1, j))
1229  END DO
1230  END DO
1231 !----------
1232 ! edges:
1233 !----------
1234  IF (grid_type .LT. 3) THEN
1235  IF (js .EQ. 1 .OR. jsd .LT. npt) THEN
1236  DO j=jsd,npt-1
1237  DO i=isd,ied
1238  utmp(i, j) = 0.5*(u(i, j)+u(i, j+1))
1239  vtmp(i, j) = 0.5*(v(i, j)+v(i+1, j))
1240  END DO
1241  END DO
1242  END IF
1243  IF (je + 1 .EQ. npy .OR. jed .GE. npy - npt) THEN
1244  DO j=npy-npt+1,jed
1245  DO i=isd,ied
1246  utmp(i, j) = 0.5*(u(i, j)+u(i, j+1))
1247  vtmp(i, j) = 0.5*(v(i, j)+v(i+1, j))
1248  END DO
1249  END DO
1250  END IF
1251  IF (is .EQ. 1 .OR. isd .LT. npt) THEN
1252  IF (npt .LT. jsd) THEN
1253  max5 = jsd
1254  ELSE
1255  max5 = npt
1256  END IF
1257  IF (npy - npt .GT. jed) THEN
1258  min5 = jed
1259  ELSE
1260  min5 = npy - npt
1261  END IF
1262  DO j=max5,min5
1263  DO i=isd,npt-1
1264  utmp(i, j) = 0.5*(u(i, j)+u(i, j+1))
1265  vtmp(i, j) = 0.5*(v(i, j)+v(i+1, j))
1266  END DO
1267  END DO
1268  END IF
1269  IF (ie + 1 .EQ. npx .OR. ied .GE. npx - npt) THEN
1270  IF (npt .LT. jsd) THEN
1271  max6 = jsd
1272  ELSE
1273  max6 = npt
1274  END IF
1275  IF (npy - npt .GT. jed) THEN
1276  min6 = jed
1277  ELSE
1278  min6 = npy - npt
1279  END IF
1280  DO j=max6,min6
1281  DO i=npx-npt+1,ied
1282  utmp(i, j) = 0.5*(u(i, j)+u(i, j+1))
1283  vtmp(i, j) = 0.5*(v(i, j)+v(i+1, j))
1284  END DO
1285  END DO
1286  END IF
1287  END IF
1288  END IF
1289  DO j=js-1-id,je+1+id
1290  DO i=is-1-id,ie+1+id
1291  ua(i, j) = (utmp(i, j)-vtmp(i, j)*cosa_s(i, j))*rsin2(i, j)
1292  va(i, j) = (vtmp(i, j)-utmp(i, j)*cosa_s(i, j))*rsin2(i, j)
1293  END DO
1294  END DO
1295  END SUBROUTINE d2a_setup
1296 end module fv_restart_adm_mod
subroutine, public d2c_setup(u, v, ua, va, uc, vc, dord4, isd, ied, jsd, jed, is, ie, js, je, npx, npy, grid_type, nested, se_corner, sw_corner, ne_corner, nw_corner, rsin_u, rsin_v, cosa_s, rsin2)
subroutine, public pushcontrol(ctype, field)
subroutine, public d2c_setup_fwd(u, v, ua, va, uc, vc, dord4, isd, ied, jsd, jed, is, ie, js, je, npx, npy, grid_type, nested, se_corner, sw_corner, ne_corner, nw_corner, rsin_u, rsin_v, cosa_s, rsin2)
subroutine, public d2c_setup_bwd(u, u_ad, v, v_ad, ua, va, uc, uc_ad, vc, vc_ad, dord4, isd, ied, jsd, jed, is, ie, js, je, npx, npy, grid_type, nested, se_corner, sw_corner, ne_corner, nw_corner, rsin_u, rsin_v, cosa_s, rsin2)
#define max(a, b)
Definition: mosaic_util.h:33
#define min(a, b)
Definition: mosaic_util.h:32
subroutine, public popcontrol(ctype, field)
Derived type containing the data.
subroutine, public d2a_setup(u, v, ua, va, dord4, isd, ied, jsd, jed, is, ie, js, je, npx, npy, grid_type, nested, cosa_s, rsin2)