SCALE-RM
mod_cpl_vars.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
20  use scale_debug
22  use scale_tracer
24  !-----------------------------------------------------------------------------
25  implicit none
26  private
27  !-----------------------------------------------------------------------------
28  !
29  !++ Public procedure
30  !
31  public :: cpl_vars_setup
32  public :: cpl_vars_finalize
33 
34  public :: cpl_putatm
35  public :: cpl_putocn
36  public :: cpl_putlnd
37  public :: cpl_puturb
38  public :: cpl_getsfc_atm
39  public :: cpl_getatm_ocn
40  public :: cpl_getatm_lnd
41  public :: cpl_getatm_urb
42 
43  !-----------------------------------------------------------------------------
44  !
45  !++ Public parameters & variables
46  !
47 
48  ! Input from ocean model
49  real(rp), public, allocatable :: ocn_sfc_temp (:,:) ! ocean surface skin temperature [K]
50  real(rp), public, allocatable :: ocn_sfc_albedo(:,:,:,:) ! ocean surface albedo (direct/diffuse,IR/near-IR/VIS) (0-1)
51  real(rp), public, allocatable :: ocn_sfc_z0m (:,:) ! ocean surface roughness length for momemtum [m]
52  real(rp), public, allocatable :: ocn_sfc_z0h (:,:) ! ocean surface roughness length for heat [m]
53  real(rp), public, allocatable :: ocn_sfc_z0e (:,:) ! ocean surface roughness length for vapor [m]
54  real(rp), public, allocatable :: ocn_sflx_mw (:,:) ! ocean surface w-momentum flux [kg/m/s2]
55  real(rp), public, allocatable :: ocn_sflx_mu (:,:) ! ocean surface u-momentum flux [kg/m/s2]
56  real(rp), public, allocatable :: ocn_sflx_mv (:,:) ! ocean surface v-momentum flux [kg/m/s2]
57  real(rp), public, allocatable :: ocn_sflx_sh (:,:) ! ocean surface sensible heat flux [J/m/s2]
58  real(rp), public, allocatable :: ocn_sflx_lh (:,:) ! ocean surface latent heat flux [J/m2/s]
59  real(rp), public, allocatable :: ocn_sflx_gh (:,:) ! ocean surface water heat flux [J/m2/s]
60  real(rp), public, allocatable :: ocn_sflx_qtrc (:,:,:) ! ocean surface tracer flux [kg/m2/s]
61  real(rp), public, allocatable :: ocn_sflx_engi (:,:) ! ocean surface internal energy flux [J/m2/s]
62  real(rp), public, allocatable :: ocn_u10 (:,:) ! ocean velocity u at 10m [m/s]
63  real(rp), public, allocatable :: ocn_v10 (:,:) ! ocean velocity v at 10m [m/s]
64  real(rp), public, allocatable :: ocn_t2 (:,:) ! ocean temperature at 2m [K]
65  real(rp), public, allocatable :: ocn_q2 (:,:) ! ocean water vapor at 2m [kg/kg]
66 
67  ! Input from land model
68  real(rp), public, allocatable :: lnd_sfc_temp (:,:) ! land surface skin temperature [K]
69  real(rp), public, allocatable :: lnd_sfc_albedo(:,:,:,:) ! land surface albedo (direct/diffuse,IR/near-IR/VIS) (0-1)
70  real(rp), public, allocatable :: lnd_sfc_z0m (:,:) ! land surface roughness length for momemtum [m]
71  real(rp), public, allocatable :: lnd_sfc_z0h (:,:) ! land surface roughness length for heat [m]
72  real(rp), public, allocatable :: lnd_sfc_z0e (:,:) ! land surface roughness length for vapor [m]
73  real(rp), public, allocatable :: lnd_sflx_mw (:,:) ! land surface w-momentum flux [kg/m/s2]
74  real(rp), public, allocatable :: lnd_sflx_mu (:,:) ! land surface u-momentum flux [kg/m/s2]
75  real(rp), public, allocatable :: lnd_sflx_mv (:,:) ! land surface v-momentum flux [kg/m/s2]
76  real(rp), public, allocatable :: lnd_sflx_sh (:,:) ! land surface sensible heat flux [J/m2/s]
77  real(rp), public, allocatable :: lnd_sflx_lh (:,:) ! land surface latent heat flux [J/m2/s]
78  real(rp), public, allocatable :: lnd_sflx_gh (:,:) ! land surface ground heat flux [J/m2/s]
79  real(rp), public, allocatable :: lnd_sflx_qtrc (:,:,:) ! land surface tracer flux [kg/m2/s]
80  real(rp), public, allocatable :: lnd_sflx_engi (:,:) ! land surface internal energy flux [J/m2/s]
81  real(rp), public, allocatable :: lnd_u10 (:,:) ! land velocity u at 10m [m/s]
82  real(rp), public, allocatable :: lnd_v10 (:,:) ! land velocity v at 10m [m/s]
83  real(rp), public, allocatable :: lnd_t2 (:,:) ! land temperature at 2m [K]
84  real(rp), public, allocatable :: lnd_q2 (:,:) ! land water vapor at 2m [kg/kg]
85 
86  ! Input from urban model
87  real(rp), public, allocatable :: urb_sfc_temp (:,:) ! urban surface skin temperature [K]
88  real(rp), public, allocatable :: urb_sfc_albedo(:,:,:,:) ! urban surface albedo (direct/diffuse,IR/near-IR/VIS) (0-1)
89  real(rp), public, allocatable :: urb_sfc_z0m (:,:) ! urban surface roughness length for momemtum [m]
90  real(rp), public, allocatable :: urb_sfc_z0h (:,:) ! urban surface roughness length for heat [m]
91  real(rp), public, allocatable :: urb_sfc_z0e (:,:) ! urban surface roughness length for vapor [m]
92  real(rp), public, allocatable :: urb_sflx_mw (:,:) ! urban surface w-momentum flux [kg/m/s2]
93  real(rp), public, allocatable :: urb_sflx_mu (:,:) ! urban surface u-momentum flux [kg/m/s2]
94  real(rp), public, allocatable :: urb_sflx_mv (:,:) ! urban surface v-momentum flux [kg/m/s2]
95  real(rp), public, allocatable :: urb_sflx_sh (:,:) ! urban surface sensible heat flux [J/m2/s] (urban grid average)
96  real(rp), public, allocatable :: urb_sflx_lh (:,:) ! urban surface latent heat flux [J/m2/s] (urban grid average)
97  real(rp), public, allocatable :: urb_sflx_shex (:,:) ! extra (anthropogenic) sensible heat flux from urban [J/m2/s] (grid average value)
98  real(rp), public, allocatable :: urb_sflx_lhex (:,:) ! extra (anthropogenic) latent heat flux from urban [J/m2/s] (grid average value)
99  real(rp), public, allocatable :: urb_sflx_qvex (:,:) ! extra (anthropogenic) water vapour flux from urban [kg/kg/m2/s] (grid average value)
100  real(rp), public, allocatable :: urb_sflx_gh (:,:) ! urban surface ground heat flux [J/m2/s]
101  real(rp), public, allocatable :: urb_sflx_qtrc (:,:,:) ! urban surface tracer flux [kg/m2/s]
102  real(rp), public, allocatable :: urb_sflx_engi (:,:) ! urban surface internal energy flux [J/m2/s]
103  real(rp), public, allocatable :: urb_u10 (:,:) ! urban velocity u at 10m [m/s]
104  real(rp), public, allocatable :: urb_v10 (:,:) ! urban velocity v at 10m [m/s]
105  real(rp), public, allocatable :: urb_t2 (:,:) ! urban temperature at 2m [K]
106  real(rp), public, allocatable :: urb_q2 (:,:) ! urban water vapor at 2m [kg/kg]
107 
108  ! Output to ocean model
109  real(rp), public, allocatable :: ocn_atm_temp (:,:) ! temperature at the lowermost atmosphere layer [K]
110  real(rp), public, allocatable :: ocn_atm_pres (:,:) ! pressure at the lowermost atmosphere layer [Pa]
111  real(rp), public, allocatable :: ocn_atm_w (:,:) ! velocity w at the lowermost atmosphere layer [m/s]
112  real(rp), public, allocatable :: ocn_atm_u (:,:) ! velocity u at the lowermost atmosphere layer [m/s]
113  real(rp), public, allocatable :: ocn_atm_v (:,:) ! velocity v at the lowermost atmosphere layer [m/s]
114  real(rp), public, allocatable :: ocn_atm_dens (:,:) ! density at the lowermost atmosphere layer [kg/m3]
115  real(rp), public, allocatable :: ocn_atm_qv (:,:) ! water vapor at the lowermost atmosphere layer [kg/kg]
116  real(rp), public, allocatable :: ocn_atm_pbl (:,:) ! the top of atmospheric mixing layer [m]
117  real(rp), public, allocatable :: ocn_atm_sfc_dens (:,:) ! surface density [kg/m3]
118  real(rp), public, allocatable :: ocn_atm_sfc_pres (:,:) ! surface pressure [Pa]
119  real(rp), public, allocatable :: ocn_atm_sflx_rad_dn(:,:,:,:) ! downward radiation flux (direct/diffuse,IR/near-IR/VIS) [J/m2/s]
120  real(rp), public, allocatable :: ocn_atm_cossza (:,:) ! cos(solar zenith angle) (0-1)
121  real(rp), public, allocatable :: ocn_atm_sflx_water (:,:) ! precipitation flux [kg/m2/s]
122  real(rp), public, allocatable :: ocn_atm_sflx_engi (:,:) ! internal energy flux [J/m2/s]
123 
124  ! Output to land model
125  real(rp), public, allocatable :: lnd_atm_temp (:,:) ! temperature at the lowermost atmosphere layer [K]
126  real(rp), public, allocatable :: lnd_atm_pres (:,:) ! pressure at the lowermost atmosphere layer [Pa]
127  real(rp), public, allocatable :: lnd_atm_w (:,:) ! velocity w at the lowermost atmosphere layer [m/s]
128  real(rp), public, allocatable :: lnd_atm_u (:,:) ! velocity u at the lowermost atmosphere layer [m/s]
129  real(rp), public, allocatable :: lnd_atm_v (:,:) ! velocity v at the lowermost atmosphere layer [m/s]
130  real(rp), public, allocatable :: lnd_atm_dens (:,:) ! density at the lowermost atmosphere layer [kg/m3]
131  real(rp), public, allocatable :: lnd_atm_qv (:,:) ! water vapor at the lowermost atmosphere layer [kg/kg]
132  real(rp), public, allocatable :: lnd_atm_pbl (:,:) ! the top of atmospheric mixing layer [m]
133  real(rp), public, allocatable :: lnd_atm_sfc_dens (:,:) ! surface density [kg/m3]
134  real(rp), public, allocatable :: lnd_atm_sfc_pres (:,:) ! surface pressure [Pa]
135  real(rp), public, allocatable :: lnd_atm_sflx_rad_dn(:,:,:,:) ! downward radiation flux (direct/diffuse,IR/near-IR/VIS) [J/m2/s]
136  real(rp), public, allocatable :: lnd_atm_cossza (:,:) ! cos(solar zenith angle) (0-1)
137  real(rp), public, allocatable :: lnd_atm_sflx_water (:,:) ! precipitation flux [kg/m2/s]
138  real(rp), public, allocatable :: lnd_atm_sflx_engi (:,:) ! internal energy flux [J/m2/s]
139 
140  ! Output to urban model
141  real(rp), public, allocatable :: urb_atm_temp (:,:) ! temperature at the lowermost atmosphere layer [K]
142  real(rp), public, allocatable :: urb_atm_pres (:,:) ! pressure at the lowermost atmosphere layer [Pa]
143  real(rp), public, allocatable :: urb_atm_w (:,:) ! velocity w at the lowermost atmosphere layer [m/s]
144  real(rp), public, allocatable :: urb_atm_u (:,:) ! velocity u at the lowermost atmosphere layer [m/s]
145  real(rp), public, allocatable :: urb_atm_v (:,:) ! velocity v at the lowermost atmosphere layer [m/s]
146  real(rp), public, allocatable :: urb_atm_dens (:,:) ! density at the lowermost atmosphere layer [kg/m3]
147  real(rp), public, allocatable :: urb_atm_qv (:,:) ! water vapor at the lowermost atmosphere layer [kg/kg]
148  real(rp), public, allocatable :: urb_atm_pbl (:,:) ! the top of atmospheric mixing layer [m]
149  real(rp), public, allocatable :: urb_atm_sfc_dens (:,:) ! surface density [kg/m3]
150  real(rp), public, allocatable :: urb_atm_sfc_pres (:,:) ! surface pressure [Pa]
151  real(rp), public, allocatable :: urb_atm_sflx_rad_dn(:,:,:,:) ! downward radiation flux (direct/diffuse,IR/near-IR/VIS) [J/m2/s]
152  real(rp), public, allocatable :: urb_atm_cossza (:,:) ! cos(solar zenith angle) (0-1)
153  real(rp), public, allocatable :: urb_atm_sflx_water (:,:) ! precipitation flux [kg/m2/s]
154  real(rp), public, allocatable :: urb_atm_sflx_engi (:,:) ! internal energy flux [J/m2/s]
155 
156  ! counter
157  real(rp), public :: cnt_putatm_ocn ! put counter for atmos to ocean
158  real(rp), public :: cnt_putatm_lnd ! put counter for atmos to land
159  real(rp), public :: cnt_putatm_urb ! put counter for atmos to urban
160  real(rp), public :: cnt_putocn ! put counter for ocean
161  real(rp), public :: cnt_putlnd ! put counter for land
162  real(rp), public :: cnt_puturb ! put counter for urban
163 
164  !-----------------------------------------------------------------------------
165  !
166  !++ Private procedure
167  !
168  !-----------------------------------------------------------------------------
169  !
170  !++ Private parameters & variables
171  !
172  !-----------------------------------------------------------------------------
173 contains
174  !-----------------------------------------------------------------------------
176  subroutine cpl_vars_setup
177  use scale_const, only: &
178  undef => const_undef
179  use scale_prc, only: &
180  prc_abort
181  use scale_landuse, only: &
186  use scale_atmos_hydrometeor, only: &
188  use mod_ocean_admin, only: &
189  ocean_do
190  use mod_land_admin, only: &
191  land_do
192  use mod_urban_admin, only: &
193  urban_do
194  use mod_lake_admin, only: &
195  lake_do
196  implicit none
197 
198  real(rp) :: checkfact
199  !---------------------------------------------------------------------------
200 
201  log_newline
202  log_info("CPL_vars_setup",*) 'Setup'
203  log_info("CPL_vars_setup",*) 'No namelists.'
204 
205  ! Check consistency of OCEAN_do and LANDUSE_fact_ocean
206  checkfact = maxval( landuse_fact_ocean(:,:) )
207  if ( .NOT. ocean_do .AND. checkfact > 0.0_rp ) then
208  log_error("CPL_vars_setup",*) 'Ocean fraction exists, but ocean component has not been called. Please check this inconsistency. STOP.', checkfact
209  call prc_abort
210  endif
211 
212  ! Check consistency of LAND_do and LANDUSE_fact_land
213  checkfact = maxval( landuse_fact_land(:,:) )
214  if ( .NOT. land_do .AND. checkfact > 0.0_rp ) then
215  log_error("CPL_vars_setup",*) 'Land fraction exists, but land component has not been called. Please check this inconsistency. STOP.', checkfact
216  call prc_abort
217  endif
218 
219  ! Check consistency of URBAN_do and LANDUSE_fact_urban
220  checkfact = maxval( landuse_fact_urban(:,:) )
221  if ( .NOT. urban_do .AND. checkfact > 0.0_rp ) then
222  log_error("CPL_vars_setup",*) 'Urban fraction exists, but urban component has not been called. Please check this inconsistency. STOP.', checkfact
223  call prc_abort
224  endif
225 
226  ! Check consistency of LAKE_do and LANDUSE_fact_lake
227  checkfact = maxval( landuse_fact_lake(:,:) )
228  if ( .NOT. lake_do .AND. checkfact > 0.0_rp ) then
229  log_error("CPL_vars_setup",*) 'Lake fraction exists, but lake component has not been called. Please check this inconsistency. STOP.', checkfact
230  call prc_abort
231  endif
232 
233 
234  allocate( ocn_sfc_temp(ia,ja) )
235  allocate( ocn_sfc_albedo(ia,ja,n_rad_dir,n_rad_rgn) )
236  allocate( ocn_sfc_z0m(ia,ja) )
237  allocate( ocn_sfc_z0h(ia,ja) )
238  allocate( ocn_sfc_z0e(ia,ja) )
239  allocate( ocn_sflx_mu(ia,ja) )
240  allocate( ocn_sflx_mv(ia,ja) )
241  allocate( ocn_sflx_mw(ia,ja) )
242  allocate( ocn_sflx_sh(ia,ja) )
243  allocate( ocn_sflx_lh(ia,ja) )
244  allocate( ocn_sflx_gh(ia,ja) )
245  allocate( ocn_sflx_qtrc(ia,ja,max(qa,1)) )
246  allocate( ocn_sflx_engi(ia,ja) )
247  allocate( ocn_u10(ia,ja) )
248  allocate( ocn_v10(ia,ja) )
249  allocate( ocn_t2(ia,ja) )
250  allocate( ocn_q2(ia,ja) )
251  ocn_sfc_temp(:,:) = undef
252  ocn_sfc_albedo(:,:,:,:) = undef
253  ocn_sfc_z0m(:,:) = undef
254  ocn_sfc_z0h(:,:) = undef
255  ocn_sfc_z0e(:,:) = undef
256  ocn_sflx_mu(:,:) = undef
257  ocn_sflx_mv(:,:) = undef
258  ocn_sflx_mw(:,:) = undef
259  ocn_sflx_sh(:,:) = undef
260  ocn_sflx_lh(:,:) = undef
261  ocn_sflx_gh(:,:) = undef
262  ocn_sflx_qtrc(:,:,:) = undef
263  ocn_sflx_engi(:,:) = undef
264  ocn_u10(:,:) = undef
265  ocn_v10(:,:) = undef
266  ocn_t2(:,:) = undef
267  ocn_q2(:,:) = undef
268  !$acc enter data create(OCN_SFC_TEMP,OCN_SFC_albedo,OCN_SFC_Z0M,OCN_SFC_Z0H,OCN_SFC_Z0E,OCN_SFLX_MU,OCN_SFLX_MV,OCN_SFLX_MW,OCN_SFLX_SH,OCN_SFLX_LH,OCN_SFLX_GH,OCN_SFLX_QTRC,OCN_SFLX_ENGI,OCN_U10,OCN_V10,OCN_T2,OCN_Q2)
269 
270  allocate( lnd_sfc_temp(ia,ja) )
271  allocate( lnd_sfc_albedo(ia,ja,n_rad_dir,n_rad_rgn) )
272  allocate( lnd_sfc_z0m(ia,ja) )
273  allocate( lnd_sfc_z0h(ia,ja) )
274  allocate( lnd_sfc_z0e(ia,ja) )
275  allocate( lnd_sflx_mu(ia,ja) )
276  allocate( lnd_sflx_mv(ia,ja) )
277  allocate( lnd_sflx_mw(ia,ja) )
278  allocate( lnd_sflx_sh(ia,ja) )
279  allocate( lnd_sflx_lh(ia,ja) )
280  allocate( lnd_sflx_gh(ia,ja) )
281  allocate( lnd_sflx_qtrc(ia,ja,max(qa,1)) )
282  allocate( lnd_sflx_engi(ia,ja) )
283  allocate( lnd_u10(ia,ja) )
284  allocate( lnd_v10(ia,ja) )
285  allocate( lnd_t2(ia,ja) )
286  allocate( lnd_q2(ia,ja) )
287  lnd_sfc_temp(:,:) = undef
288  lnd_sfc_albedo(:,:,:,:) = undef
289  lnd_sfc_z0m(:,:) = undef
290  lnd_sfc_z0h(:,:) = undef
291  lnd_sfc_z0e(:,:) = undef
292  lnd_sflx_mu(:,:) = undef
293  lnd_sflx_mv(:,:) = undef
294  lnd_sflx_mw(:,:) = undef
295  lnd_sflx_sh(:,:) = undef
296  lnd_sflx_lh(:,:) = undef
297  lnd_sflx_gh(:,:) = undef
298  lnd_sflx_qtrc(:,:,:) = undef
299  lnd_sflx_engi(:,:) = undef
300  lnd_u10(:,:) = undef
301  lnd_v10(:,:) = undef
302  lnd_t2(:,:) = undef
303  lnd_q2(:,:) = undef
304  !$acc enter data create(LND_SFC_TEMP,LND_SFC_albedo,LND_SFC_Z0M,LND_SFC_Z0H,LND_SFC_Z0E,LND_SFLX_MU,LND_SFLX_MV,LND_SFLX_MW,LND_SFLX_SH,LND_SFLX_LH,LND_SFLX_GH,LND_SFLX_QTRC,LND_SFLX_ENGI,LND_U10,LND_V10,LND_T2,LND_Q2)
305 
306  allocate( urb_sfc_temp(ia,ja) )
307  allocate( urb_sfc_albedo(ia,ja,n_rad_dir,n_rad_rgn) )
308  allocate( urb_sfc_z0m(ia,ja) )
309  allocate( urb_sfc_z0h(ia,ja) )
310  allocate( urb_sfc_z0e(ia,ja) )
311  allocate( urb_sflx_mu(ia,ja) )
312  allocate( urb_sflx_mv(ia,ja) )
313  allocate( urb_sflx_mw(ia,ja) )
314  allocate( urb_sflx_sh(ia,ja) )
315  allocate( urb_sflx_lh(ia,ja) )
316  allocate( urb_sflx_shex(ia,ja) )
317  allocate( urb_sflx_lhex(ia,ja) )
318  allocate( urb_sflx_qvex(ia,ja) )
319  allocate( urb_sflx_gh(ia,ja) )
320  allocate( urb_sflx_qtrc(ia,ja,max(qa,1)) )
321  allocate( urb_sflx_engi(ia,ja) )
322  allocate( urb_u10(ia,ja) )
323  allocate( urb_v10(ia,ja) )
324  allocate( urb_t2(ia,ja) )
325  allocate( urb_q2(ia,ja) )
326  urb_sfc_temp(:,:) = undef
327  urb_sfc_albedo(:,:,:,:) = undef
328  urb_sfc_z0m(:,:) = undef
329  urb_sfc_z0h(:,:) = undef
330  urb_sfc_z0e(:,:) = undef
331  urb_sflx_mu(:,:) = undef
332  urb_sflx_mv(:,:) = undef
333  urb_sflx_mw(:,:) = undef
334  urb_sflx_sh(:,:) = undef
335  urb_sflx_lh(:,:) = undef
336  urb_sflx_shex(:,:) = undef
337  urb_sflx_lhex(:,:) = undef
338  urb_sflx_qvex(:,:) = undef
339  urb_sflx_gh(:,:) = undef
340  urb_sflx_qtrc(:,:,:) = undef
341  urb_sflx_engi(:,:) = undef
342  urb_u10(:,:) = undef
343  urb_v10(:,:) = undef
344  urb_t2(:,:) = undef
345  urb_q2(:,:) = undef
346  !$acc enter data create(URB_SFC_TEMP,URB_SFC_albedo,URB_SFC_Z0M,URB_SFC_Z0H,URB_SFC_Z0E,URB_SFLX_MU,URB_SFLX_MV,URB_SFLX_MW,URB_SFLX_SH,URB_SFLX_LH,URB_SFLX_SHEX,URB_SFLX_LHEX,URB_SFLX_QVEX,URB_SFLX_GH,URB_SFLX_QTRC,URB_SFLX_ENGI,URB_U10,URB_V10,URB_T2,URB_Q2)
347 
348  allocate( ocn_atm_temp(ia,ja) )
349  allocate( ocn_atm_pres(ia,ja) )
350  allocate( ocn_atm_w(ia,ja) )
351  allocate( ocn_atm_u(ia,ja) )
352  allocate( ocn_atm_v(ia,ja) )
353  allocate( ocn_atm_dens(ia,ja) )
354  allocate( ocn_atm_qv(ia,ja) )
355  allocate( ocn_atm_pbl(ia,ja) )
356  allocate( ocn_atm_sfc_dens(ia,ja) )
357  allocate( ocn_atm_sfc_pres(ia,ja) )
359  allocate( ocn_atm_cossza(ia,ja) )
360  allocate( ocn_atm_sflx_water(ia,ja) )
361  allocate( ocn_atm_sflx_engi(ia,ja) )
362  ocn_atm_temp(:,:) = undef
363  ocn_atm_pres(:,:) = undef
364  ocn_atm_w(:,:) = undef
365  ocn_atm_u(:,:) = undef
366  ocn_atm_v(:,:) = undef
367  ocn_atm_dens(:,:) = undef
368  ocn_atm_qv(:,:) = undef
369  ocn_atm_pbl(:,:) = undef
370  ocn_atm_sfc_dens(:,:) = undef
371  ocn_atm_sfc_pres(:,:) = undef
372  ocn_atm_sflx_rad_dn(:,:,:,:) = undef
373  ocn_atm_cossza(:,:) = undef
374  ocn_atm_sflx_water(:,:) = undef
375  ocn_atm_sflx_engi(:,:) = undef
376  !$acc enter data create(OCN_ATM_TEMP,OCN_ATM_PRES,OCN_ATM_W,OCN_ATM_U,OCN_ATM_V,OCN_ATM_DENS,OCN_ATM_QV,OCN_ATM_PBL,OCN_ATM_SFC_DENS,OCN_ATM_SFC_PRES,OCN_ATM_SFLX_rad_dn,OCN_ATM_cosSZA,OCN_ATM_SFLX_water,OCN_ATM_SFLX_ENGI)
377 
378  allocate( lnd_atm_temp(ia,ja) )
379  allocate( lnd_atm_pres(ia,ja) )
380  allocate( lnd_atm_w(ia,ja) )
381  allocate( lnd_atm_u(ia,ja) )
382  allocate( lnd_atm_v(ia,ja) )
383  allocate( lnd_atm_dens(ia,ja) )
384  allocate( lnd_atm_qv(ia,ja) )
385  allocate( lnd_atm_pbl(ia,ja) )
386  allocate( lnd_atm_sfc_dens(ia,ja) )
387  allocate( lnd_atm_sfc_pres(ia,ja) )
389  allocate( lnd_atm_cossza(ia,ja) )
390  allocate( lnd_atm_sflx_water(ia,ja) )
391  allocate( lnd_atm_sflx_engi(ia,ja) )
392  lnd_atm_temp(:,:) = undef
393  lnd_atm_pres(:,:) = undef
394  lnd_atm_w(:,:) = undef
395  lnd_atm_u(:,:) = undef
396  lnd_atm_v(:,:) = undef
397  lnd_atm_dens(:,:) = undef
398  lnd_atm_qv(:,:) = undef
399  lnd_atm_pbl(:,:) = undef
400  lnd_atm_sfc_dens(:,:) = undef
401  lnd_atm_sfc_pres(:,:) = undef
402  lnd_atm_sflx_rad_dn(:,:,:,:) = undef
403  lnd_atm_cossza(:,:) = undef
404  lnd_atm_sflx_water(:,:) = undef
405  lnd_atm_sflx_engi(:,:) = undef
406  !$acc enter data create(LND_ATM_TEMP,LND_ATM_PRES,LND_ATM_W,LND_ATM_U,LND_ATM_V,LND_ATM_DENS,LND_ATM_QV,LND_ATM_PBL,LND_ATM_SFC_DENS,LND_ATM_SFC_PRES,LND_ATM_SFLX_rad_dn,LND_ATM_cosSZA,LND_ATM_SFLX_water,LND_ATM_SFLX_ENGI)
407 
408  allocate( urb_atm_temp(ia,ja) )
409  allocate( urb_atm_pres(ia,ja) )
410  allocate( urb_atm_w(ia,ja) )
411  allocate( urb_atm_u(ia,ja) )
412  allocate( urb_atm_v(ia,ja) )
413  allocate( urb_atm_dens(ia,ja) )
414  allocate( urb_atm_qv(ia,ja) )
415  allocate( urb_atm_pbl(ia,ja) )
416  allocate( urb_atm_sfc_dens(ia,ja) )
417  allocate( urb_atm_sfc_pres(ia,ja) )
419  allocate( urb_atm_cossza(ia,ja) )
420  allocate( urb_atm_sflx_water(ia,ja) )
421  allocate( urb_atm_sflx_engi(ia,ja) )
422  urb_atm_temp(:,:) = undef
423  urb_atm_pres(:,:) = undef
424  urb_atm_w(:,:) = undef
425  urb_atm_u(:,:) = undef
426  urb_atm_v(:,:) = undef
427  urb_atm_dens(:,:) = undef
428  urb_atm_qv(:,:) = undef
429  urb_atm_pbl(:,:) = undef
430  urb_atm_sfc_dens(:,:) = undef
431  urb_atm_sfc_pres(:,:) = undef
432  urb_atm_sflx_rad_dn(:,:,:,:) = undef
433  urb_atm_cossza(:,:) = undef
434  urb_atm_sflx_water(:,:) = undef
435  urb_atm_sflx_engi(:,:) = undef
436  !$acc enter data create(URB_ATM_TEMP,URB_ATM_PRES,URB_ATM_W,URB_ATM_U,URB_ATM_V,URB_ATM_DENS,URB_ATM_QV,URB_ATM_PBL,URB_ATM_SFC_DENS,URB_ATM_SFC_PRES,URB_ATM_SFLX_rad_dn,URB_ATM_cosSZA,URB_ATM_SFLX_water,URB_ATM_SFLX_ENGI)
437 
438  ! counter intialize
439  cnt_putatm_ocn = 0.0_rp
440  cnt_putatm_lnd = 0.0_rp
441  cnt_putatm_urb = 0.0_rp
442  cnt_putocn = 0.0_rp
443  cnt_putlnd = 0.0_rp
444  cnt_puturb = 0.0_rp
445 
446  if ( atmos_hydrometeor_dry ) then
447  ocn_atm_qv = 0.0_rp
448  lnd_atm_qv = 0.0_rp
449  urb_atm_qv = 0.0_rp
450  endif
451 
452  return
453  end subroutine cpl_vars_setup
454 
455  !-----------------------------------------------------------------------------
457  subroutine cpl_vars_finalize
458  implicit none
459  !---------------------------------------------------------------------------
460 
461  log_newline
462  log_info("CPL_vars_finalize",*) 'Finalize'
463 
464  !$acc exit data delete(OCN_SFC_TEMP,OCN_SFC_albedo,OCN_SFC_Z0M,OCN_SFC_Z0H,OCN_SFC_Z0E,OCN_SFLX_MU,OCN_SFLX_MV,OCN_SFLX_MW,OCN_SFLX_SH,OCN_SFLX_LH,OCN_SFLX_GH,OCN_SFLX_QTRC,OCN_SFLX_ENGI,OCN_U10,OCN_V10,OCN_T2,OCN_Q2)
465  deallocate( ocn_sfc_temp )
466  deallocate( ocn_sfc_albedo )
467  deallocate( ocn_sfc_z0m )
468  deallocate( ocn_sfc_z0h )
469  deallocate( ocn_sfc_z0e )
470  deallocate( ocn_sflx_mu )
471  deallocate( ocn_sflx_mv )
472  deallocate( ocn_sflx_mw )
473  deallocate( ocn_sflx_sh )
474  deallocate( ocn_sflx_lh )
475  deallocate( ocn_sflx_gh )
476  deallocate( ocn_sflx_qtrc )
477  deallocate( ocn_sflx_engi )
478  deallocate( ocn_u10 )
479  deallocate( ocn_v10 )
480  deallocate( ocn_t2 )
481  deallocate( ocn_q2 )
482 
483  !$acc exit data delete(LND_SFC_TEMP,LND_SFC_albedo,LND_SFC_Z0M,LND_SFC_Z0H,LND_SFC_Z0E,LND_SFLX_MU,LND_SFLX_MV,LND_SFLX_MW,LND_SFLX_SH,LND_SFLX_LH,LND_SFLX_GH,LND_SFLX_QTRC,LND_SFLX_ENGI,LND_U10,LND_V10,LND_T2,LND_Q2)
484  deallocate( lnd_sfc_temp )
485  deallocate( lnd_sfc_albedo )
486  deallocate( lnd_sfc_z0m )
487  deallocate( lnd_sfc_z0h )
488  deallocate( lnd_sfc_z0e )
489  deallocate( lnd_sflx_mu )
490  deallocate( lnd_sflx_mv )
491  deallocate( lnd_sflx_mw )
492  deallocate( lnd_sflx_sh )
493  deallocate( lnd_sflx_lh )
494  deallocate( lnd_sflx_gh )
495  deallocate( lnd_sflx_qtrc )
496  deallocate( lnd_sflx_engi )
497  deallocate( lnd_u10 )
498  deallocate( lnd_v10 )
499  deallocate( lnd_t2 )
500  deallocate( lnd_q2 )
501 
502  !$acc exit data delete(URB_SFC_TEMP,URB_SFC_albedo,URB_SFC_Z0M,URB_SFC_Z0H,URB_SFC_Z0E,URB_SFLX_MU,URB_SFLX_MV,URB_SFLX_MW,URB_SFLX_SH,URB_SFLX_LH,URB_SFLX_SHEX,URB_SFLX_LHEX,URB_SFLX_QVEX,URB_SFLX_GH,URB_SFLX_QTRC,URB_SFLX_ENGI,URB_U10,URB_V10,URB_T2,URB_Q2)
503  deallocate( urb_sfc_temp )
504  deallocate( urb_sfc_albedo )
505  deallocate( urb_sfc_z0m )
506  deallocate( urb_sfc_z0h )
507  deallocate( urb_sfc_z0e )
508  deallocate( urb_sflx_mu )
509  deallocate( urb_sflx_mv )
510  deallocate( urb_sflx_mw )
511  deallocate( urb_sflx_sh )
512  deallocate( urb_sflx_lh )
513  deallocate( urb_sflx_shex )
514  deallocate( urb_sflx_lhex )
515  deallocate( urb_sflx_qvex )
516  deallocate( urb_sflx_gh )
517  deallocate( urb_sflx_qtrc )
518  deallocate( urb_sflx_engi )
519  deallocate( urb_u10 )
520  deallocate( urb_v10 )
521  deallocate( urb_t2 )
522  deallocate( urb_q2 )
523 
524  !$acc exit data delete(OCN_ATM_TEMP,OCN_ATM_PRES,OCN_ATM_W,OCN_ATM_U,OCN_ATM_V,OCN_ATM_DENS,OCN_ATM_QV,OCN_ATM_PBL,OCN_ATM_SFC_DENS,OCN_ATM_SFC_PRES,OCN_ATM_SFLX_rad_dn,OCN_ATM_cosSZA,OCN_ATM_SFLX_water,OCN_ATM_SFLX_ENGI)
525  deallocate( ocn_atm_temp )
526  deallocate( ocn_atm_pres )
527  deallocate( ocn_atm_w )
528  deallocate( ocn_atm_u )
529  deallocate( ocn_atm_v )
530  deallocate( ocn_atm_dens )
531  deallocate( ocn_atm_qv )
532  deallocate( ocn_atm_pbl )
533  deallocate( ocn_atm_sfc_dens )
534  deallocate( ocn_atm_sfc_pres )
535  deallocate( ocn_atm_sflx_rad_dn )
536  deallocate( ocn_atm_cossza )
537  deallocate( ocn_atm_sflx_water )
538  deallocate( ocn_atm_sflx_engi )
539 
540  !$acc exit data delete(LND_ATM_TEMP,LND_ATM_PRES,LND_ATM_W,LND_ATM_U,LND_ATM_V,LND_ATM_DENS,LND_ATM_QV,LND_ATM_PBL,LND_ATM_SFC_DENS,LND_ATM_SFC_PRES,LND_ATM_SFLX_rad_dn,LND_ATM_cosSZA,LND_ATM_SFLX_water,LND_ATM_SFLX_ENGI)
541  deallocate( lnd_atm_temp )
542  deallocate( lnd_atm_pres )
543  deallocate( lnd_atm_w )
544  deallocate( lnd_atm_u )
545  deallocate( lnd_atm_v )
546  deallocate( lnd_atm_dens )
547  deallocate( lnd_atm_qv )
548  deallocate( lnd_atm_pbl )
549  deallocate( lnd_atm_sfc_dens )
550  deallocate( lnd_atm_sfc_pres )
551  deallocate( lnd_atm_sflx_rad_dn )
552  deallocate( lnd_atm_cossza )
553  deallocate( lnd_atm_sflx_water )
554  deallocate( lnd_atm_sflx_engi )
555 
556  !$acc exit data delete(URB_ATM_TEMP,URB_ATM_PRES,URB_ATM_W,URB_ATM_U,URB_ATM_V,URB_ATM_DENS,URB_ATM_QV,URB_ATM_PBL,URB_ATM_SFC_DENS,URB_ATM_SFC_PRES,URB_ATM_SFLX_rad_dn,URB_ATM_cosSZA,URB_ATM_SFLX_water,URB_ATM_SFLX_ENGI)
557  deallocate( urb_atm_temp )
558  deallocate( urb_atm_pres )
559  deallocate( urb_atm_w )
560  deallocate( urb_atm_u )
561  deallocate( urb_atm_v )
562  deallocate( urb_atm_dens )
563  deallocate( urb_atm_qv )
564  deallocate( urb_atm_pbl )
565  deallocate( urb_atm_sfc_dens )
566  deallocate( urb_atm_sfc_pres )
567  deallocate( urb_atm_sflx_rad_dn )
568  deallocate( urb_atm_cossza )
569  deallocate( urb_atm_sflx_water )
570  deallocate( urb_atm_sflx_engi )
571 
572  return
573  end subroutine cpl_vars_finalize
574 
575  !-----------------------------------------------------------------------------
576  subroutine cpl_putatm( &
577  TEMP, &
578  PRES, &
579  W, &
580  U, &
581  V, &
582  DENS, &
583  QV, &
584  PBL, &
585  SFC_DENS, &
586  SFC_PRES, &
587  SFLX_rad_dn, &
588  cosSZA, &
589  SFLX_water, &
590  SFLX_ENGI, &
591  countup )
592  use scale_atmos_hydrometeor, only: &
594  cv_water, &
595  cv_ice, &
596  lhf
597  implicit none
598 
599  real(rp), intent(in) :: temp (ia,ja)
600  real(rp), intent(in) :: pres (ia,ja)
601  real(rp), intent(in) :: w (ia,ja)
602  real(rp), intent(in) :: u (ia,ja)
603  real(rp), intent(in) :: v (ia,ja)
604  real(rp), intent(in) :: dens (ia,ja)
605  real(rp), intent(in) :: qv (ia,ja)
606  real(rp), intent(in) :: pbl (ia,ja)
607  real(rp), intent(in) :: sfc_dens (ia,ja)
608  real(rp), intent(in) :: sfc_pres (ia,ja)
609  real(rp), intent(in) :: sflx_rad_dn(ia,ja,n_rad_dir,n_rad_rgn)
610  real(rp), intent(in) :: cossza (ia,ja)
611  real(rp), intent(in) :: sflx_water (ia,ja)
612  real(rp), intent(in) :: sflx_engi (ia,ja)
613  logical, intent(in) :: countup
614 
615  integer :: i, j, idir, irgn
616  !---------------------------------------------------------------------------
617 
618  !$omp parallel do default(none) private(i,j,idir,irgn) OMP_SCHEDULE_ &
619  !$omp shared(JS,JE,IS,IE,OCN_ATM_TEMP,OCN_ATM_PRES,OCN_ATM_W,OCN_ATM_U) &
620  !$omp shared(OCN_ATM_V,OCN_ATM_DENS,CNT_putATM_OCN,TEMP,PRES,W,U,V,DENS) &
621  !$omp shared(ATMOS_HYDROMETEOR_dry,OCN_ATM_QV,OCN_ATM_PBL,OCN_ATM_SFC_DENS,OCN_ATM_SFC_PRES,OCN_ATM_SFLX_rad_dn) &
622  !$omp shared(OCN_ATM_cosSZA,OCN_ATM_SFLX_water,OCN_ATM_SFLX_ENGI,QV,PBL,SFC_DENS,SFC_PRES) &
623  !$omp shared(SFLX_rad_dn,cosSZA,SFLX_water,SFLX_ENGI) &
624  !$omp shared(LND_ATM_TEMP,LND_ATM_PRES,LND_ATM_W,LND_ATM_U,LND_ATM_V,LND_ATM_DENS) &
625  !$omp shared(LND_ATM_QV,LND_ATM_PBL,LND_ATM_SFC_DENS,LND_ATM_SFC_PRES,LND_ATM_SFLX_rad_dn,LND_ATM_cosSZA) &
626  !$omp shared(LND_ATM_SFLX_water,LND_ATM_SFLX_ENGI) &
627  !$omp shared(URB_ATM_TEMP,URB_ATM_PRES,URB_ATM_W,URB_ATM_U,URB_ATM_V,URB_ATM_DENS) &
628  !$omp shared(URB_ATM_QV,URB_ATM_PBL,URB_ATM_SFC_DENS,URB_ATM_SFC_PRES,URB_ATM_SFLX_rad_dn,URB_ATM_cosSZA) &
629  !$omp shared(URB_ATM_SFLX_water,URB_ATM_SFLX_ENGI,CNT_putATM_URB,CNT_putATM_LND)
630  !$acc kernels
631  do j = js, je
632  do i = is, ie
633  ! for ocean
634  ocn_atm_temp(i,j) = ocn_atm_temp(i,j) * cnt_putatm_ocn + temp(i,j)
635  ocn_atm_pres(i,j) = ocn_atm_pres(i,j) * cnt_putatm_ocn + pres(i,j)
636  ocn_atm_w(i,j) = ocn_atm_w(i,j) * cnt_putatm_ocn + w(i,j)
637  ocn_atm_u(i,j) = ocn_atm_u(i,j) * cnt_putatm_ocn + u(i,j)
638  ocn_atm_v(i,j) = ocn_atm_v(i,j) * cnt_putatm_ocn + v(i,j)
639  ocn_atm_dens(i,j) = ocn_atm_dens(i,j) * cnt_putatm_ocn + dens(i,j)
640  if ( .NOT. atmos_hydrometeor_dry ) &
641  ocn_atm_qv(i,j) = ocn_atm_qv(i,j) * cnt_putatm_ocn + qv(i,j)
642  ocn_atm_pbl(i,j) = ocn_atm_pbl(i,j) * cnt_putatm_ocn + pbl(i,j)
643  ocn_atm_sfc_dens(i,j) = ocn_atm_sfc_dens(i,j) * cnt_putatm_ocn + sfc_dens(i,j)
644  ocn_atm_sfc_pres(i,j) = ocn_atm_sfc_pres(i,j) * cnt_putatm_ocn + sfc_pres(i,j)
645  ocn_atm_cossza(i,j) = ocn_atm_cossza(i,j) * cnt_putatm_ocn + cossza(i,j)
646  ocn_atm_sflx_water(i,j) = ocn_atm_sflx_water(i,j) * cnt_putatm_ocn + sflx_water(i,j)
647  ocn_atm_sflx_engi(i,j) = ocn_atm_sflx_engi(i,j) * cnt_putatm_ocn + sflx_engi(i,j)
648  do irgn = i_r_ir, i_r_vis
649  do idir = i_r_direct, i_r_diffuse
650  ocn_atm_sflx_rad_dn(i,j,idir,irgn) = ocn_atm_sflx_rad_dn(i,j,idir,irgn) * cnt_putatm_ocn + sflx_rad_dn(i,j,idir,irgn)
651  enddo
652  enddo
653 
654  ! for land
655  lnd_atm_temp(i,j) = lnd_atm_temp(i,j) * cnt_putatm_lnd + temp(i,j)
656  lnd_atm_pres(i,j) = lnd_atm_pres(i,j) * cnt_putatm_lnd + pres(i,j)
657  lnd_atm_w(i,j) = lnd_atm_w(i,j) * cnt_putatm_lnd + w(i,j)
658  lnd_atm_u(i,j) = lnd_atm_u(i,j) * cnt_putatm_lnd + u(i,j)
659  lnd_atm_v(i,j) = lnd_atm_v(i,j) * cnt_putatm_lnd + v(i,j)
660  lnd_atm_dens(i,j) = lnd_atm_dens(i,j) * cnt_putatm_lnd + dens(i,j)
661  if ( .NOT. atmos_hydrometeor_dry ) &
662  lnd_atm_qv(i,j) = lnd_atm_qv(i,j) * cnt_putatm_lnd + qv(i,j)
663  lnd_atm_pbl(i,j) = lnd_atm_pbl(i,j) * cnt_putatm_lnd + pbl(i,j)
664  lnd_atm_sfc_dens(i,j) = lnd_atm_sfc_dens(i,j) * cnt_putatm_lnd + sfc_dens(i,j)
665  lnd_atm_sfc_pres(i,j) = lnd_atm_sfc_pres(i,j) * cnt_putatm_lnd + sfc_pres(i,j)
666  lnd_atm_cossza(i,j) = lnd_atm_cossza(i,j) * cnt_putatm_lnd + cossza(i,j)
667  lnd_atm_sflx_water(i,j) = lnd_atm_sflx_water(i,j) * cnt_putatm_lnd + sflx_water(i,j)
668  lnd_atm_sflx_engi(i,j) = lnd_atm_sflx_engi(i,j) * cnt_putatm_lnd + sflx_engi(i,j)
669  do irgn = i_r_ir, i_r_vis
670  do idir = i_r_direct, i_r_diffuse
671  lnd_atm_sflx_rad_dn(i,j,idir,irgn) = lnd_atm_sflx_rad_dn(i,j,idir,irgn) * cnt_putatm_lnd + sflx_rad_dn(i,j,idir,irgn)
672  enddo
673  enddo
674 
675  ! for urban
676  urb_atm_temp(i,j) = urb_atm_temp(i,j) * cnt_putatm_urb + temp(i,j)
677  urb_atm_pres(i,j) = urb_atm_pres(i,j) * cnt_putatm_urb + pres(i,j)
678  urb_atm_w(i,j) = urb_atm_w(i,j) * cnt_putatm_urb + w(i,j)
679  urb_atm_u(i,j) = urb_atm_u(i,j) * cnt_putatm_urb + u(i,j)
680  urb_atm_v(i,j) = urb_atm_v(i,j) * cnt_putatm_urb + v(i,j)
681  urb_atm_dens(i,j) = urb_atm_dens(i,j) * cnt_putatm_urb + dens(i,j)
682  if ( .NOT. atmos_hydrometeor_dry ) &
683  urb_atm_qv(i,j) = urb_atm_qv(i,j) * cnt_putatm_urb + qv(i,j)
684  urb_atm_pbl(i,j) = urb_atm_pbl(i,j) * cnt_putatm_urb + pbl(i,j)
685  urb_atm_sfc_dens(i,j) = urb_atm_sfc_dens(i,j) * cnt_putatm_urb + sfc_dens(i,j)
686  urb_atm_sfc_pres(i,j) = urb_atm_sfc_pres(i,j) * cnt_putatm_urb + sfc_pres(i,j)
687  urb_atm_cossza(i,j) = urb_atm_cossza(i,j) * cnt_putatm_urb + cossza(i,j)
688  urb_atm_sflx_water(i,j) = urb_atm_sflx_water(i,j) * cnt_putatm_urb + sflx_water(i,j)
689  urb_atm_sflx_engi(i,j) = urb_atm_sflx_engi(i,j) * cnt_putatm_urb + sflx_engi(i,j)
690  do irgn = i_r_ir, i_r_vis
691  do idir = i_r_direct, i_r_diffuse
692  urb_atm_sflx_rad_dn(i,j,idir,irgn) = urb_atm_sflx_rad_dn(i,j,idir,irgn) * cnt_putatm_urb + sflx_rad_dn(i,j,idir,irgn)
693  enddo
694  enddo
695 
696 
697  ! for ocean
698  ocn_atm_temp(i,j) = ocn_atm_temp(i,j) / ( cnt_putatm_ocn + 1.0_rp )
699  ocn_atm_pres(i,j) = ocn_atm_pres(i,j) / ( cnt_putatm_ocn + 1.0_rp )
700  ocn_atm_w(i,j) = ocn_atm_w(i,j) / ( cnt_putatm_ocn + 1.0_rp )
701  ocn_atm_u(i,j) = ocn_atm_u(i,j) / ( cnt_putatm_ocn + 1.0_rp )
702  ocn_atm_v(i,j) = ocn_atm_v(i,j) / ( cnt_putatm_ocn + 1.0_rp )
703  ocn_atm_dens(i,j) = ocn_atm_dens(i,j) / ( cnt_putatm_ocn + 1.0_rp )
704  ocn_atm_qv(i,j) = ocn_atm_qv(i,j) / ( cnt_putatm_ocn + 1.0_rp )
705  ocn_atm_pbl(i,j) = ocn_atm_pbl(i,j) / ( cnt_putatm_ocn + 1.0_rp )
706  ocn_atm_sfc_dens(i,j) = ocn_atm_sfc_dens(i,j) / ( cnt_putatm_ocn + 1.0_rp )
707  ocn_atm_sfc_pres(i,j) = ocn_atm_sfc_pres(i,j) / ( cnt_putatm_ocn + 1.0_rp )
708  ocn_atm_cossza(i,j) = ocn_atm_cossza(i,j) / ( cnt_putatm_ocn + 1.0_rp )
709  ocn_atm_sflx_water(i,j) = ocn_atm_sflx_water(i,j) / ( cnt_putatm_ocn + 1.0_rp )
710  ocn_atm_sflx_engi(i,j) = ocn_atm_sflx_engi(i,j) / ( cnt_putatm_ocn + 1.0_rp )
711  do irgn = i_r_ir, i_r_vis
712  do idir = i_r_direct, i_r_diffuse
713  ocn_atm_sflx_rad_dn(i,j,idir,irgn) = ocn_atm_sflx_rad_dn(i,j,idir,irgn) / ( cnt_putatm_ocn + 1.0_rp )
714  enddo
715  enddo
716 
717  ! for land
718  lnd_atm_temp(i,j) = lnd_atm_temp(i,j) / ( cnt_putatm_lnd + 1.0_rp )
719  lnd_atm_pres(i,j) = lnd_atm_pres(i,j) / ( cnt_putatm_lnd + 1.0_rp )
720  lnd_atm_w(i,j) = lnd_atm_w(i,j) / ( cnt_putatm_lnd + 1.0_rp )
721  lnd_atm_u(i,j) = lnd_atm_u(i,j) / ( cnt_putatm_lnd + 1.0_rp )
722  lnd_atm_v(i,j) = lnd_atm_v(i,j) / ( cnt_putatm_lnd + 1.0_rp )
723  lnd_atm_dens(i,j) = lnd_atm_dens(i,j) / ( cnt_putatm_lnd + 1.0_rp )
724  lnd_atm_qv(i,j) = lnd_atm_qv(i,j) / ( cnt_putatm_lnd + 1.0_rp )
725  lnd_atm_pbl(i,j) = lnd_atm_pbl(i,j) / ( cnt_putatm_lnd + 1.0_rp )
726  lnd_atm_sfc_dens(i,j) = lnd_atm_sfc_dens(i,j) / ( cnt_putatm_lnd + 1.0_rp )
727  lnd_atm_sfc_pres(i,j) = lnd_atm_sfc_pres(i,j) / ( cnt_putatm_lnd + 1.0_rp )
728  lnd_atm_cossza(i,j) = lnd_atm_cossza(i,j) / ( cnt_putatm_lnd + 1.0_rp )
729  lnd_atm_sflx_water(i,j) = lnd_atm_sflx_water(i,j) / ( cnt_putatm_lnd + 1.0_rp )
730  lnd_atm_sflx_engi(i,j) = lnd_atm_sflx_engi(i,j) / ( cnt_putatm_lnd + 1.0_rp )
731  do irgn = i_r_ir, i_r_vis
732  do idir = i_r_direct, i_r_diffuse
733  lnd_atm_sflx_rad_dn(i,j,idir,irgn) = lnd_atm_sflx_rad_dn(i,j,idir,irgn) / ( cnt_putatm_lnd + 1.0_rp )
734  enddo
735  enddo
736 
737  ! for urban
738  urb_atm_temp(i,j) = urb_atm_temp(i,j) / ( cnt_putatm_urb + 1.0_rp )
739  urb_atm_pres(i,j) = urb_atm_pres(i,j) / ( cnt_putatm_urb + 1.0_rp )
740  urb_atm_w(i,j) = urb_atm_w(i,j) / ( cnt_putatm_urb + 1.0_rp )
741  urb_atm_u(i,j) = urb_atm_u(i,j) / ( cnt_putatm_urb + 1.0_rp )
742  urb_atm_v(i,j) = urb_atm_v(i,j) / ( cnt_putatm_urb + 1.0_rp )
743  urb_atm_dens(i,j) = urb_atm_dens(i,j) / ( cnt_putatm_urb + 1.0_rp )
744  urb_atm_qv(i,j) = urb_atm_qv(i,j) / ( cnt_putatm_urb + 1.0_rp )
745  urb_atm_pbl(i,j) = urb_atm_pbl(i,j) / ( cnt_putatm_urb + 1.0_rp )
746  urb_atm_sfc_dens(i,j) = urb_atm_sfc_dens(i,j) / ( cnt_putatm_urb + 1.0_rp )
747  urb_atm_sfc_pres(i,j) = urb_atm_sfc_pres(i,j) / ( cnt_putatm_urb + 1.0_rp )
748  urb_atm_cossza(i,j) = urb_atm_cossza(i,j) / ( cnt_putatm_urb + 1.0_rp )
749  urb_atm_sflx_water(i,j) = urb_atm_sflx_water(i,j) / ( cnt_putatm_urb + 1.0_rp )
750  urb_atm_sflx_engi(i,j) = urb_atm_sflx_engi(i,j) / ( cnt_putatm_urb + 1.0_rp )
751  do irgn = i_r_ir, i_r_vis
752  do idir = i_r_direct, i_r_diffuse
753  urb_atm_sflx_rad_dn(i,j,idir,irgn) = urb_atm_sflx_rad_dn(i,j,idir,irgn) / ( cnt_putatm_urb + 1.0_rp )
754  enddo
755  enddo
756  enddo
757  enddo
758  !$acc end kernels
759 
760  if( countup ) then
761  cnt_putatm_ocn = cnt_putatm_ocn + 1.0_rp
762  cnt_putatm_lnd = cnt_putatm_lnd + 1.0_rp
763  cnt_putatm_urb = cnt_putatm_urb + 1.0_rp
764  endif
765 
766  return
767  end subroutine cpl_putatm
768 
769  !-----------------------------------------------------------------------------
770  subroutine cpl_putocn( &
771  SFC_TEMP, &
772  SFC_albedo, &
773  SFC_Z0M, &
774  SFC_Z0H, &
775  SFC_Z0E, &
776  SFLX_MW, &
777  SFLX_MU, &
778  SFLX_MV, &
779  SFLX_SH, &
780  SFLX_LH, &
781  SFLX_GH, &
782  SFLX_QTRC, &
783  U10, &
784  V10, &
785  T2, &
786  Q2, &
787  mask, &
788  countup )
789  implicit none
790 
791  real(rp), intent(in) :: sfc_temp (ia,ja)
792  real(rp), intent(in) :: sfc_albedo(ia,ja,n_rad_dir,n_rad_rgn)
793  real(rp), intent(in) :: sfc_z0m (ia,ja)
794  real(rp), intent(in) :: sfc_z0h (ia,ja)
795  real(rp), intent(in) :: sfc_z0e (ia,ja)
796  real(rp), intent(in) :: sflx_mw (ia,ja)
797  real(rp), intent(in) :: sflx_mu (ia,ja)
798  real(rp), intent(in) :: sflx_mv (ia,ja)
799  real(rp), intent(in) :: sflx_sh (ia,ja)
800  real(rp), intent(in) :: sflx_lh (ia,ja)
801  real(rp), intent(in) :: sflx_gh (ia,ja)
802  real(rp), intent(in) :: sflx_qtrc (ia,ja,qa)
803  real(rp), intent(in) :: u10 (ia,ja)
804  real(rp), intent(in) :: v10 (ia,ja)
805  real(rp), intent(in) :: t2 (ia,ja)
806  real(rp), intent(in) :: q2 (ia,ja)
807  logical, intent(in) :: mask (ia,ja)
808  logical, intent(in) :: countup
809 
810  integer :: i, j, iq, idir, irgn
811  !---------------------------------------------------------------------------
812 
813  !$omp parallel do default(none) private(i,j,idir,irgn) OMP_SCHEDULE_ &
814  !$omp shared(JS,JE,IS,IE,QA) &
815  !$omp shared(OCN_SFC_TEMP,OCN_SFC_albedo,OCN_SFC_Z0M,OCN_SFC_Z0H,OCN_SFC_Z0E) &
816  !$omp shared(OCN_SFLX_MW,OCN_SFLX_MU,OCN_SFLX_MV,OCN_SFLX_SH,OCN_SFLX_LH,OCN_SFLX_GH,OCN_SFLX_QTRC,OCN_SFLX_ENGI,OCN_U10,OCN_V10,OCN_T2,OCN_Q2) &
817  !$omp shared(SFC_TEMP,SFC_albedo,SFC_Z0M,SFC_Z0H,SFC_Z0E,SFLX_MW,SFLX_MU,SFLX_MV,SFLX_SH,SFLX_LH,SFLX_GH,SFLX_QTRC,U10,V10,T2,Q2) &
818  !$omp shared(mask) &
819  !$omp shared(CNT_putOCN) &
820  !$omp shared(TRACER_CV,TRACER_ENGI0)
821  !$acc kernels
822  do j = js, je
823  do i = is, ie
824  if ( mask(i,j) ) then
825  ocn_sfc_temp(i,j) = ocn_sfc_temp(i,j) * cnt_putocn + sfc_temp(i,j)
826  ocn_sfc_z0m(i,j) = ocn_sfc_z0m(i,j) * cnt_putocn + sfc_z0m(i,j)
827  ocn_sfc_z0h(i,j) = ocn_sfc_z0h(i,j) * cnt_putocn + sfc_z0h(i,j)
828  ocn_sfc_z0e(i,j) = ocn_sfc_z0e(i,j) * cnt_putocn + sfc_z0e(i,j)
829  ocn_sflx_mw(i,j) = ocn_sflx_mw(i,j) * cnt_putocn + sflx_mw(i,j)
830  ocn_sflx_mu(i,j) = ocn_sflx_mu(i,j) * cnt_putocn + sflx_mu(i,j)
831  ocn_sflx_mv(i,j) = ocn_sflx_mv(i,j) * cnt_putocn + sflx_mv(i,j)
832  ocn_sflx_sh(i,j) = ocn_sflx_sh(i,j) * cnt_putocn + sflx_sh(i,j)
833  ocn_sflx_lh(i,j) = ocn_sflx_lh(i,j) * cnt_putocn + sflx_lh(i,j)
834  ocn_sflx_gh(i,j) = ocn_sflx_gh(i,j) * cnt_putocn + sflx_gh(i,j)
836  do iq = 1, qa
837  ocn_sflx_qtrc(i,j,iq) = ocn_sflx_qtrc(i,j,iq) * cnt_putocn + sflx_qtrc(i,j,iq)
838  ocn_sflx_engi(i,j) = ocn_sflx_engi(i,j) + sflx_qtrc(i,j,iq) * ( tracer_cv(iq) * sfc_temp(i,j) + tracer_engi0(iq) )
839  end do
840  ocn_u10(i,j) = ocn_u10(i,j) * cnt_putocn + u10(i,j)
841  ocn_v10(i,j) = ocn_v10(i,j) * cnt_putocn + v10(i,j)
842  ocn_t2(i,j) = ocn_t2(i,j) * cnt_putocn + t2(i,j)
843  ocn_q2(i,j) = ocn_q2(i,j) * cnt_putocn + q2(i,j)
844  do irgn = i_r_ir, i_r_vis
845  do idir = i_r_direct, i_r_diffuse
846  ocn_sfc_albedo(i,j,idir,irgn) = ocn_sfc_albedo(i,j,idir,irgn) * cnt_putocn + sfc_albedo(i,j,idir,irgn)
847  enddo
848  enddo
849 
850  ocn_sfc_temp(i,j) = ocn_sfc_temp(i,j) / ( cnt_putocn + 1.0_rp )
851  ocn_sfc_z0m(i,j) = ocn_sfc_z0m(i,j) / ( cnt_putocn + 1.0_rp )
852  ocn_sfc_z0h(i,j) = ocn_sfc_z0h(i,j) / ( cnt_putocn + 1.0_rp )
853  ocn_sfc_z0e(i,j) = ocn_sfc_z0e(i,j) / ( cnt_putocn + 1.0_rp )
854  ocn_sflx_mw(i,j) = ocn_sflx_mw(i,j) / ( cnt_putocn + 1.0_rp )
855  ocn_sflx_mu(i,j) = ocn_sflx_mu(i,j) / ( cnt_putocn + 1.0_rp )
856  ocn_sflx_mv(i,j) = ocn_sflx_mv(i,j) / ( cnt_putocn + 1.0_rp )
857  ocn_sflx_sh(i,j) = ocn_sflx_sh(i,j) / ( cnt_putocn + 1.0_rp )
858  ocn_sflx_lh(i,j) = ocn_sflx_lh(i,j) / ( cnt_putocn + 1.0_rp )
859  ocn_sflx_gh(i,j) = ocn_sflx_gh(i,j) / ( cnt_putocn + 1.0_rp )
860  ocn_sflx_qtrc(i,j,:) = ocn_sflx_qtrc(i,j,:) / ( cnt_putocn + 1.0_rp )
861  ocn_sflx_engi(i,j) = ocn_sflx_engi(i,j) / ( cnt_putocn + 1.0_rp )
862  ocn_u10(i,j) = ocn_u10(i,j) / ( cnt_putocn + 1.0_rp )
863  ocn_v10(i,j) = ocn_v10(i,j) / ( cnt_putocn + 1.0_rp )
864  ocn_t2(i,j) = ocn_t2(i,j) / ( cnt_putocn + 1.0_rp )
865  ocn_q2(i,j) = ocn_q2(i,j) / ( cnt_putocn + 1.0_rp )
866  do irgn = i_r_ir, i_r_vis
867  do idir = i_r_direct, i_r_diffuse
868  ocn_sfc_albedo(i,j,idir,irgn) = ocn_sfc_albedo(i,j,idir,irgn) / ( cnt_putocn + 1.0_rp )
869  enddo
870  enddo
871  end if
872  enddo
873  enddo
874  !$acc end kernels
875 
876  if( countup ) then
877  cnt_putocn = cnt_putocn + 1.0_rp
878  endif
879 
880  return
881  end subroutine cpl_putocn
882 
883  !-----------------------------------------------------------------------------
884  subroutine cpl_putlnd( &
885  SFC_TEMP, &
886  SFC_albedo, &
887  SFC_Z0M, &
888  SFC_Z0H, &
889  SFC_Z0E, &
890  SFLX_MW, &
891  SFLX_MU, &
892  SFLX_MV, &
893  SFLX_SH, &
894  SFLX_LH, &
895  SFLX_GH, &
896  SFLX_QTRC, &
897  U10, &
898  V10, &
899  T2, &
900  Q2, &
901  mask, &
902  countup )
903  implicit none
904 
905  real(rp), intent(in) :: sfc_temp (ia,ja)
906  real(rp), intent(in) :: sfc_albedo(ia,ja,n_rad_dir,n_rad_rgn)
907  real(rp), intent(in) :: sfc_z0m (ia,ja)
908  real(rp), intent(in) :: sfc_z0h (ia,ja)
909  real(rp), intent(in) :: sfc_z0e (ia,ja)
910  real(rp), intent(in) :: sflx_mw (ia,ja)
911  real(rp), intent(in) :: sflx_mu (ia,ja)
912  real(rp), intent(in) :: sflx_mv (ia,ja)
913  real(rp), intent(in) :: sflx_sh (ia,ja)
914  real(rp), intent(in) :: sflx_lh (ia,ja)
915  real(rp), intent(in) :: sflx_gh (ia,ja)
916  real(rp), intent(in) :: sflx_qtrc (ia,ja,qa)
917  real(rp), intent(in) :: u10 (ia,ja)
918  real(rp), intent(in) :: v10 (ia,ja)
919  real(rp), intent(in) :: t2 (ia,ja)
920  real(rp), intent(in) :: q2 (ia,ja)
921  logical, intent(in) :: mask (ia,ja)
922  logical, intent(in) :: countup
923 
924  integer :: i, j, iq, idir, irgn
925  !---------------------------------------------------------------------------
926 
927  !$omp parallel do default(none) &
928  !$omp shared(JS,JE,IS,IE,QA) &
929  !$omp shared(LND_SFC_TEMP,LND_SFC_albedo,LND_SFC_Z0M,LND_SFC_Z0H,LND_SFC_Z0E) &
930  !$omp shared(LND_SFLX_MW,LND_SFLX_MU,LND_SFLX_MV,LND_SFLX_SH,LND_SFLX_LH,LND_SFLX_GH,LND_SFLX_QTRC,LND_SFLX_ENGI) &
931  !$omp shared(LND_U10,LND_V10,LND_T2,LND_Q2,CNT_putLND,SFC_TEMP,SFC_albedo,SFC_Z0M,SFC_Z0H) &
932  !$omp shared(SFC_Z0E,SFLX_MW,SFLX_MU,SFLX_MV,SFLX_SH,SFLX_LH,SFLX_GH,SFLX_QTRC,U10,V10,T2,Q2) &
933  !$omp shared(mask) &
934  !$omp shared(TRACER_CV,TRACER_ENGI0) &
935  !$omp private(i,j,idir,irgn) OMP_SCHEDULE_
936  !$acc kernels
937  do j = js, je
938  do i = is, ie
939  if ( mask(i,j) ) then
940  lnd_sfc_temp(i,j) = lnd_sfc_temp(i,j) * cnt_putlnd + sfc_temp(i,j)
941  lnd_sfc_z0m(i,j) = lnd_sfc_z0m(i,j) * cnt_putlnd + sfc_z0m(i,j)
942  lnd_sfc_z0h(i,j) = lnd_sfc_z0h(i,j) * cnt_putlnd + sfc_z0h(i,j)
943  lnd_sfc_z0e(i,j) = lnd_sfc_z0e(i,j) * cnt_putlnd + sfc_z0e(i,j)
944  lnd_sflx_mw(i,j) = lnd_sflx_mw(i,j) * cnt_putlnd + sflx_mw(i,j)
945  lnd_sflx_mu(i,j) = lnd_sflx_mu(i,j) * cnt_putlnd + sflx_mu(i,j)
946  lnd_sflx_mv(i,j) = lnd_sflx_mv(i,j) * cnt_putlnd + sflx_mv(i,j)
947  lnd_sflx_sh(i,j) = lnd_sflx_sh(i,j) * cnt_putlnd + sflx_sh(i,j)
948  lnd_sflx_lh(i,j) = lnd_sflx_lh(i,j) * cnt_putlnd + sflx_lh(i,j)
949  lnd_sflx_gh(i,j) = lnd_sflx_gh(i,j) * cnt_putlnd + sflx_gh(i,j)
951  do iq = 1, qa
952  lnd_sflx_qtrc(i,j,iq) = lnd_sflx_qtrc(i,j,iq) * cnt_putlnd + sflx_qtrc(i,j,iq)
953  lnd_sflx_engi(i,j) = lnd_sflx_engi(i,j) + sflx_qtrc(i,j,iq) * ( tracer_cv(iq) * sfc_temp(i,j) + tracer_engi0(iq) )
954  end do
955  lnd_u10(i,j) = lnd_u10(i,j) * cnt_putlnd + u10(i,j)
956  lnd_v10(i,j) = lnd_v10(i,j) * cnt_putlnd + v10(i,j)
957  lnd_t2(i,j) = lnd_t2(i,j) * cnt_putlnd + t2(i,j)
958  lnd_q2(i,j) = lnd_q2(i,j) * cnt_putlnd + q2(i,j)
959  do irgn = i_r_ir, i_r_vis
960  do idir = i_r_direct, i_r_diffuse
961  lnd_sfc_albedo(i,j,idir,irgn) = lnd_sfc_albedo(i,j,idir,irgn) * cnt_putlnd + sfc_albedo(i,j,idir,irgn)
962  enddo
963  enddo
964 
965  lnd_sfc_temp(i,j) = lnd_sfc_temp(i,j) / ( cnt_putlnd + 1.0_rp )
966  lnd_sfc_z0m(i,j) = lnd_sfc_z0m(i,j) / ( cnt_putlnd + 1.0_rp )
967  lnd_sfc_z0h(i,j) = lnd_sfc_z0h(i,j) / ( cnt_putlnd + 1.0_rp )
968  lnd_sfc_z0e(i,j) = lnd_sfc_z0e(i,j) / ( cnt_putlnd + 1.0_rp )
969  lnd_sflx_mw(i,j) = lnd_sflx_mw(i,j) / ( cnt_putlnd + 1.0_rp )
970  lnd_sflx_mu(i,j) = lnd_sflx_mu(i,j) / ( cnt_putlnd + 1.0_rp )
971  lnd_sflx_mv(i,j) = lnd_sflx_mv(i,j) / ( cnt_putlnd + 1.0_rp )
972  lnd_sflx_sh(i,j) = lnd_sflx_sh(i,j) / ( cnt_putlnd + 1.0_rp )
973  lnd_sflx_lh(i,j) = lnd_sflx_lh(i,j) / ( cnt_putlnd + 1.0_rp )
974  lnd_sflx_gh(i,j) = lnd_sflx_gh(i,j) / ( cnt_putlnd + 1.0_rp )
975  lnd_sflx_qtrc(i,j,:) = lnd_sflx_qtrc(i,j,:) / ( cnt_putlnd + 1.0_rp )
976  lnd_sflx_engi(i,j) = lnd_sflx_engi(i,j) / ( cnt_putlnd + 1.0_rp )
977  lnd_u10(i,j) = lnd_u10(i,j) / ( cnt_putlnd + 1.0_rp )
978  lnd_v10(i,j) = lnd_v10(i,j) / ( cnt_putlnd + 1.0_rp )
979  lnd_t2(i,j) = lnd_t2(i,j) / ( cnt_putlnd + 1.0_rp )
980  lnd_q2(i,j) = lnd_q2(i,j) / ( cnt_putlnd + 1.0_rp )
981  do irgn = i_r_ir, i_r_vis
982  do idir = i_r_direct, i_r_diffuse
983  lnd_sfc_albedo(i,j,idir,irgn) = lnd_sfc_albedo(i,j,idir,irgn) / ( cnt_putlnd + 1.0_rp )
984  enddo
985  enddo
986  end if
987  enddo
988  enddo
989  !$acc end kernels
990 
991  if( countup ) then
992  cnt_putlnd = cnt_putlnd + 1.0_rp
993  endif
994 
995  return
996  end subroutine cpl_putlnd
997 
998  !-----------------------------------------------------------------------------
999  subroutine cpl_puturb( &
1000  SFC_TEMP, &
1001  SFC_albedo, &
1002  SFC_Z0M, &
1003  SFC_Z0H, &
1004  SFC_Z0E, &
1005  SFLX_MW, &
1006  SFLX_MU, &
1007  SFLX_MV, &
1008  SFLX_SH, &
1009  SFLX_LH, &
1010  SFLX_SHEX, &
1011  SFLX_LHEX, &
1012  SFLX_QVEX, &
1013  SFLX_GH, &
1014  SFLX_QTRC, &
1015  U10, &
1016  V10, &
1017  T2, &
1018  Q2, &
1019  mask, &
1020  countup )
1021  implicit none
1022 
1023  real(rp), intent(in) :: sfc_temp (ia,ja)
1024  real(rp), intent(in) :: sfc_albedo(ia,ja,n_rad_dir,n_rad_rgn)
1025  real(rp), intent(in) :: sfc_z0m (ia,ja)
1026  real(rp), intent(in) :: sfc_z0h (ia,ja)
1027  real(rp), intent(in) :: sfc_z0e (ia,ja)
1028  real(rp), intent(in) :: sflx_mw (ia,ja)
1029  real(rp), intent(in) :: sflx_mu (ia,ja)
1030  real(rp), intent(in) :: sflx_mv (ia,ja)
1031  real(rp), intent(in) :: sflx_sh (ia,ja)
1032  real(rp), intent(in) :: sflx_lh (ia,ja)
1033  real(rp), intent(in) :: sflx_shex (ia,ja)
1034  real(rp), intent(in) :: sflx_lhex (ia,ja)
1035  real(rp), intent(in) :: sflx_qvex (ia,ja)
1036  real(rp), intent(in) :: sflx_gh (ia,ja)
1037  real(rp), intent(in) :: sflx_qtrc (ia,ja,qa)
1038  real(rp), intent(in) :: u10 (ia,ja)
1039  real(rp), intent(in) :: v10 (ia,ja)
1040  real(rp), intent(in) :: t2 (ia,ja)
1041  real(rp), intent(in) :: q2 (ia,ja)
1042  logical, intent(in) :: mask (ia,ja)
1043  logical, intent(in) :: countup
1044 
1045  integer :: i, j, iq, idir, irgn
1046  !---------------------------------------------------------------------------
1047 
1048  !$omp parallel do default(none) OMP_SCHEDULE_ &
1049  !$omp shared(JS,JE,IS,IE,QA, &
1050  !$omp URB_SFC_TEMP,URB_SFC_albedo,URB_SFC_Z0M,URB_SFC_Z0H,URB_SFC_Z0E, &
1051  !$omp URB_SFLX_MW,URB_SFLX_MU,URB_SFLX_MV,URB_SFLX_SH,URB_SFLX_LH,URB_SFLX_SHEX,URB_SFLX_LHEX,URB_SFLX_QVEX,URB_SFLX_GH, &
1052  !$omp URB_SFLX_QTRC,URB_SFLX_ENGI,URB_U10,URB_V10,URB_T2,URB_Q2,CNT_putURB, &
1053  !$omp SFC_TEMP,SFC_albedo,SFC_Z0M,SFC_Z0H,SFC_Z0E,SFLX_MW,SFLX_MU,SFLX_MV,SFLX_SH,SFLX_LH,SFLX_SHEX,SFLX_LHEX,SFLX_QVEX,SFLX_GH,SFLX_QTRC,U10,V10,T2,Q2, &
1054  !$omp mask, &
1055  !$omp TRACER_CV,TRACER_ENGI0)
1056  !$acc kernels
1057  do j = js, je
1058  do i = is, ie
1059  if ( mask(i,j) ) then
1060  urb_sfc_temp(i,j) = urb_sfc_temp(i,j) * cnt_puturb + sfc_temp(i,j)
1061  urb_sfc_z0m(i,j) = urb_sfc_z0m(i,j) * cnt_puturb + sfc_z0m(i,j)
1062  urb_sfc_z0h(i,j) = urb_sfc_z0h(i,j) * cnt_puturb + sfc_z0h(i,j)
1063  urb_sfc_z0e(i,j) = urb_sfc_z0e(i,j) * cnt_puturb + sfc_z0e(i,j)
1064  urb_sflx_mw(i,j) = urb_sflx_mw(i,j) * cnt_puturb + sflx_mw(i,j)
1065  urb_sflx_mu(i,j) = urb_sflx_mu(i,j) * cnt_puturb + sflx_mu(i,j)
1066  urb_sflx_mv(i,j) = urb_sflx_mv(i,j) * cnt_puturb + sflx_mv(i,j)
1067  urb_sflx_sh(i,j) = urb_sflx_sh(i,j) * cnt_puturb + sflx_sh(i,j)
1068  urb_sflx_lh(i,j) = urb_sflx_lh(i,j) * cnt_puturb + sflx_lh(i,j)
1069  urb_sflx_shex(i,j) = urb_sflx_shex(i,j) * cnt_puturb + sflx_shex(i,j)
1070  urb_sflx_lhex(i,j) = urb_sflx_lhex(i,j) * cnt_puturb + sflx_lhex(i,j)
1071  urb_sflx_qvex(i,j) = urb_sflx_qvex(i,j) * cnt_puturb + sflx_qvex(i,j)
1072  urb_sflx_gh(i,j) = urb_sflx_gh(i,j) * cnt_puturb + sflx_gh(i,j)
1073  urb_sflx_engi(i,j) = urb_sflx_engi(i,j) * cnt_puturb
1074  do iq = 1, qa
1075  urb_sflx_qtrc(i,j,iq) = urb_sflx_qtrc(i,j,iq) * cnt_puturb + sflx_qtrc(i,j,iq)
1076  urb_sflx_engi(i,j) = urb_sflx_engi(i,j) + sflx_qtrc(i,j,iq) * ( tracer_cv(iq) * sfc_temp(i,j) + tracer_engi0(iq) )
1077  end do
1078  urb_u10(i,j) = urb_u10(i,j) * cnt_puturb + u10(i,j)
1079  urb_v10(i,j) = urb_v10(i,j) * cnt_puturb + v10(i,j)
1080  urb_t2(i,j) = urb_t2(i,j) * cnt_puturb + t2(i,j)
1081  urb_q2(i,j) = urb_q2(i,j) * cnt_puturb + q2(i,j)
1082  do irgn = i_r_ir, i_r_vis
1083  do idir = i_r_direct, i_r_diffuse
1084  urb_sfc_albedo(i,j,idir,irgn) = urb_sfc_albedo(i,j,idir,irgn) * cnt_puturb + sfc_albedo(i,j,idir,irgn)
1085  enddo
1086  enddo
1087 
1088  urb_sfc_temp(i,j) = urb_sfc_temp(i,j) / ( cnt_puturb + 1.0_rp )
1089  urb_sfc_z0m(i,j) = urb_sfc_z0m(i,j) / ( cnt_puturb + 1.0_rp )
1090  urb_sfc_z0h(i,j) = urb_sfc_z0h(i,j) / ( cnt_puturb + 1.0_rp )
1091  urb_sfc_z0e(i,j) = urb_sfc_z0e(i,j) / ( cnt_puturb + 1.0_rp )
1092  urb_sflx_mw(i,j) = urb_sflx_mw(i,j) / ( cnt_puturb + 1.0_rp )
1093  urb_sflx_mu(i,j) = urb_sflx_mu(i,j) / ( cnt_puturb + 1.0_rp )
1094  urb_sflx_mv(i,j) = urb_sflx_mv(i,j) / ( cnt_puturb + 1.0_rp )
1095  urb_sflx_sh(i,j) = urb_sflx_sh(i,j) / ( cnt_puturb + 1.0_rp )
1096  urb_sflx_lh(i,j) = urb_sflx_lh(i,j) / ( cnt_puturb + 1.0_rp )
1097  urb_sflx_shex(i,j) = urb_sflx_shex(i,j) / ( cnt_puturb + 1.0_rp )
1098  urb_sflx_lhex(i,j) = urb_sflx_lhex(i,j) / ( cnt_puturb + 1.0_rp )
1099  urb_sflx_qvex(i,j) = urb_sflx_qvex(i,j) / ( cnt_puturb + 1.0_rp )
1100  urb_sflx_gh(i,j) = urb_sflx_gh(i,j) / ( cnt_puturb + 1.0_rp )
1101  urb_sflx_qtrc(i,j,:) = urb_sflx_qtrc(i,j,:) / ( cnt_puturb + 1.0_rp )
1102  urb_sflx_engi(i,j) = urb_sflx_engi(i,j) / ( cnt_puturb + 1.0_rp )
1103  urb_u10(i,j) = urb_u10(i,j) / ( cnt_puturb + 1.0_rp )
1104  urb_v10(i,j) = urb_v10(i,j) / ( cnt_puturb + 1.0_rp )
1105  urb_t2(i,j) = urb_t2(i,j) / ( cnt_puturb + 1.0_rp )
1106  urb_q2(i,j) = urb_q2(i,j) / ( cnt_puturb + 1.0_rp )
1107  do irgn = i_r_ir, i_r_vis
1108  do idir = i_r_direct, i_r_diffuse
1109  urb_sfc_albedo(i,j,idir,irgn) = urb_sfc_albedo(i,j,idir,irgn) / ( cnt_puturb + 1.0_rp )
1110  enddo
1111  enddo
1112  end if
1113  enddo
1114  enddo
1115  !$acc end kernels
1116 
1117  if( countup ) then
1118  cnt_puturb = cnt_puturb + 1.0_rp
1119  endif
1120 
1121  return
1122  end subroutine cpl_puturb
1123 
1124  !-----------------------------------------------------------------------------
1125  subroutine cpl_getsfc_atm( &
1126  SFC_TEMP, &
1127  SFC_albedo, &
1128  SFC_Z0M, &
1129  SFC_Z0H, &
1130  SFC_Z0E, &
1131  SFLX_MW, &
1132  SFLX_MU, &
1133  SFLX_MV, &
1134  SFLX_SH, &
1135  SFLX_LH, &
1136  SFLX_SHEX, &
1137  SFLX_LHEX, &
1138  SFLX_QVEX, &
1139  SFLX_GH, &
1140  SFLX_QTRC, &
1141  SFLX_ENGI, &
1142  U10, &
1143  V10, &
1144  T2, &
1145  Q2 )
1146  use scale_landuse, only: &
1147  fact_ocean => landuse_fact_ocean, &
1148  fact_land => landuse_fact_land, &
1149  fact_urban => landuse_fact_urban
1150  use scale_const, only: &
1151  eps => const_eps
1152  use scale_atmos_hydrometeor, only: &
1153  i_qv
1154  implicit none
1155 
1156  real(rp), intent(out) :: sfc_temp (ia,ja)
1157  real(rp), intent(out) :: sfc_albedo(ia,ja,n_rad_dir,n_rad_rgn)
1158  real(rp), intent(out) :: sfc_z0m (ia,ja)
1159  real(rp), intent(out) :: sfc_z0h (ia,ja)
1160  real(rp), intent(out) :: sfc_z0e (ia,ja)
1161  real(rp), intent(out) :: sflx_mw (ia,ja)
1162  real(rp), intent(out) :: sflx_mu (ia,ja)
1163  real(rp), intent(out) :: sflx_mv (ia,ja)
1164  real(rp), intent(out) :: sflx_sh (ia,ja)
1165  real(rp), intent(out) :: sflx_lh (ia,ja)
1166  real(rp), intent(out) :: sflx_shex (ia,ja)
1167  real(rp), intent(out) :: sflx_lhex (ia,ja)
1168  real(rp), intent(out) :: sflx_qvex (ia,ja)
1169  real(rp), intent(out) :: sflx_gh (ia,ja)
1170  real(rp), intent(out) :: sflx_qtrc (ia,ja,qa)
1171  real(rp), intent(out) :: sflx_engi (ia,ja)
1172  real(rp), intent(out) :: u10 (ia,ja)
1173  real(rp), intent(out) :: v10 (ia,ja)
1174  real(rp), intent(out) :: t2 (ia,ja)
1175  real(rp), intent(out) :: q2 (ia,ja)
1176 
1177  integer :: i, j, idir, irgn, iq
1178  !---------------------------------------------------------------------------
1179 
1180  !$omp parallel do default(none) &
1181  !$omp shared(JS,JE,IS,IE,QA,I_QV,EPS,SFLX_QTRC,SFLX_ENGI,SFC_TEMP,SFC_albedo,SFC_Z0M,SFC_Z0H,SFC_Z0E) &
1182  !$omp shared(SFLX_MW,SFLX_MU,SFLX_MV,SFLX_SH,SFLX_LH,SFLX_SHEX,SFLX_LHEX,SFLX_QVEX,SFLX_GH,U10,V10,T2,Q2) &
1183  !$omp shared(fact_ocean,fact_land,fact_urban,OCN_SFC_TEMP,LND_SFC_TEMP,URB_SFC_TEMP,OCN_SFC_albedo) &
1184  !$omp shared(LND_SFC_albedo,URB_SFC_albedo,OCN_SFC_Z0M,LND_SFC_Z0M,URB_SFC_Z0M) &
1185  !$omp shared(OCN_SFC_Z0H,LND_SFC_Z0H,URB_SFC_Z0H,OCN_SFC_Z0E,LND_SFC_Z0E,URB_SFC_Z0E) &
1186  !$omp shared(OCN_SFLX_MW,LND_SFLX_MW,URB_SFLX_MW,OCN_SFLX_MU,LND_SFLX_MU,URB_SFLX_MU) &
1187  !$omp shared(OCN_SFLX_MV,LND_SFLX_MV,URB_SFLX_MV) &
1188  !$omp shared(OCN_SFLX_SH,LND_SFLX_SH,URB_SFLX_SH,URB_SFLX_SHEX,OCN_SFLX_LH,LND_SFLX_LH,URB_SFLX_LH,URB_SFLX_LHEX,URB_SFLX_QVEX) &
1189  !$omp shared(OCN_SFLX_GH,LND_SFLX_GH,URB_SFLX_GH,OCN_SFLX_QTRC,LND_SFLX_QTRC,URB_SFLX_QTRC,OCN_SFLX_ENGI,LND_SFLX_ENGI,URB_SFLX_ENGI) &
1190  !$omp shared(OCN_U10,LND_U10,URB_U10,OCN_V10,LND_V10,URB_V10,OCN_T2,LND_T2,URB_T2,OCN_Q2,LND_Q2,URB_Q2) &
1191  !$omp private(i,j,iq) OMP_SCHEDULE_
1192  !$acc kernels
1193  do j = js, je
1194  do i = is, ie
1195  sfc_temp(i,j) = fact_ocean(i,j) * ocn_sfc_temp(i,j) &
1196  + fact_land(i,j) * lnd_sfc_temp(i,j) &
1197  + fact_urban(i,j) * urb_sfc_temp(i,j)
1198 
1199  do irgn = i_r_ir, i_r_vis
1200  do idir = i_r_direct, i_r_diffuse
1201  sfc_albedo(i,j,idir,irgn) = fact_ocean(i,j) * ocn_sfc_albedo(i,j,idir,irgn) &
1202  + fact_land(i,j) * lnd_sfc_albedo(i,j,idir,irgn) &
1203  + fact_urban(i,j) * urb_sfc_albedo(i,j,idir,irgn)
1204  enddo
1205  enddo
1206 
1207  sfc_z0m(i,j) = fact_ocean(i,j) * ocn_sfc_z0m(i,j) &
1208  + fact_land(i,j) * lnd_sfc_z0m(i,j) &
1209  + fact_urban(i,j) * urb_sfc_z0m(i,j)
1210 
1211  sfc_z0h(i,j) = fact_ocean(i,j) * ocn_sfc_z0h(i,j) &
1212  + fact_land(i,j) * lnd_sfc_z0h(i,j) &
1213  + fact_urban(i,j) * urb_sfc_z0h(i,j)
1214 
1215  sfc_z0e(i,j) = fact_ocean(i,j) * ocn_sfc_z0e(i,j) &
1216  + fact_land(i,j) * lnd_sfc_z0e(i,j) &
1217  + fact_urban(i,j) * urb_sfc_z0e(i,j)
1218 
1219  sflx_mw(i,j) = fact_ocean(i,j) * ocn_sflx_mw(i,j) &
1220  + fact_land(i,j) * lnd_sflx_mw(i,j) &
1221  + fact_urban(i,j) * urb_sflx_mw(i,j)
1222 
1223  sflx_mu(i,j) = fact_ocean(i,j) * ocn_sflx_mu(i,j) &
1224  + fact_land(i,j) * lnd_sflx_mu(i,j) &
1225  + fact_urban(i,j) * urb_sflx_mu(i,j)
1226 
1227  sflx_mv(i,j) = fact_ocean(i,j) * ocn_sflx_mv(i,j) &
1228  + fact_land(i,j) * lnd_sflx_mv(i,j) &
1229  + fact_urban(i,j) * urb_sflx_mv(i,j)
1230 
1231  sflx_shex(i,j) = real(ceiling(fact_urban(i,j)-eps), kind=rp) * urb_sflx_shex(i,j) ! only a grid cell where urban exists
1232 
1233  sflx_lhex(i,j) = real(ceiling(fact_urban(i,j)-eps), kind=rp) * urb_sflx_lhex(i,j) ! only a grid cell where urban exists
1234 
1235  sflx_qvex(i,j) = real(ceiling(fact_urban(i,j)-eps), kind=rp) * urb_sflx_qvex(i,j) ! only a grid cell where urban exists
1236 
1237  sflx_sh(i,j) = fact_ocean(i,j) * ocn_sflx_sh(i,j) &
1238  + fact_land(i,j) * lnd_sflx_sh(i,j) &
1239  + fact_urban(i,j) * urb_sflx_sh(i,j) &
1240  + sflx_shex(i,j) ! grid average
1241 
1242  sflx_lh(i,j) = fact_ocean(i,j) * ocn_sflx_lh(i,j) &
1243  + fact_land(i,j) * lnd_sflx_lh(i,j) &
1244  + fact_urban(i,j) * urb_sflx_lh(i,j) &
1245  + sflx_lhex(i,j) ! grid average
1246 
1247  sflx_gh(i,j) = fact_ocean(i,j) * ocn_sflx_gh(i,j) &
1248  + fact_land(i,j) * lnd_sflx_gh(i,j) &
1249  + fact_urban(i,j) * urb_sflx_gh(i,j)
1250 
1251  do iq = 1, qa
1252  sflx_qtrc(i,j,iq) = fact_ocean(i,j) * ocn_sflx_qtrc(i,j,iq) &
1253  + fact_land(i,j) * lnd_sflx_qtrc(i,j,iq) &
1254  + fact_urban(i,j) * urb_sflx_qtrc(i,j,iq)
1255  enddo
1256  sflx_qtrc(i,j,i_qv) = sflx_qtrc(i,j,i_qv) + sflx_qvex(i,j)
1257 
1258  sflx_engi(i,j) = fact_ocean(i,j) * ocn_sflx_engi(i,j) &
1259  + fact_land(i,j) * lnd_sflx_engi(i,j) &
1260  + fact_urban(i,j) * urb_sflx_engi(i,j)
1261 
1262  u10(i,j) = fact_ocean(i,j) * ocn_u10(i,j) &
1263  + fact_land(i,j) * lnd_u10(i,j) &
1264  + fact_urban(i,j) * urb_u10(i,j)
1265 
1266  v10(i,j) = fact_ocean(i,j) * ocn_v10(i,j) &
1267  + fact_land(i,j) * lnd_v10(i,j) &
1268  + fact_urban(i,j) * urb_v10(i,j)
1269 
1270  t2(i,j) = fact_ocean(i,j) * ocn_t2(i,j) &
1271  + fact_land(i,j) * lnd_t2(i,j) &
1272  + fact_urban(i,j) * urb_t2(i,j)
1273 
1274  q2(i,j) = fact_ocean(i,j) * ocn_q2(i,j) &
1275  + fact_land(i,j) * lnd_q2(i,j) &
1276  + fact_urban(i,j) * urb_q2(i,j)
1277  enddo
1278  enddo
1279  !$acc end kernels
1280 
1281  cnt_putocn = 0.0_rp
1282  cnt_putlnd = 0.0_rp
1283  cnt_puturb = 0.0_rp
1284 
1285  return
1286  end subroutine cpl_getsfc_atm
1287 
1288  !-----------------------------------------------------------------------------
1289  subroutine cpl_getatm_ocn( &
1290  TEMP, &
1291  PRES, &
1292  W, &
1293  U, &
1294  V, &
1295  DENS, &
1296  QV, &
1297  PBL, &
1298  SFC_DENS, &
1299  SFC_PRES, &
1300  SFLX_rad_dn, &
1301  cosSZA, &
1302  SFLX_water, &
1303  SFLX_ENGI )
1304  implicit none
1305 
1306  real(rp), intent(out) :: temp (ia,ja)
1307  real(rp), intent(out) :: pres (ia,ja)
1308  real(rp), intent(out) :: w (ia,ja)
1309  real(rp), intent(out) :: u (ia,ja)
1310  real(rp), intent(out) :: v (ia,ja)
1311  real(rp), intent(out) :: dens (ia,ja)
1312  real(rp), intent(out) :: qv (ia,ja)
1313  real(rp), intent(out) :: pbl (ia,ja)
1314  real(rp), intent(out) :: sfc_dens (ia,ja)
1315  real(rp), intent(out) :: sfc_pres (ia,ja)
1316  real(rp), intent(out) :: sflx_rad_dn(ia,ja,n_rad_dir,n_rad_rgn)
1317  real(rp), intent(out) :: cossza (ia,ja)
1318  real(rp), intent(out) :: sflx_water (ia,ja)
1319  real(rp), intent(out) :: sflx_engi (ia,ja)
1320 
1321  integer :: i, j, idir, irgn
1322  !---------------------------------------------------------------------------
1323 
1324 !OCL XFILL
1325  !$omp parallel do default(none) private(i,j,idir,irgn) OMP_SCHEDULE_ &
1326  !$omp shared(JS,JE,IS,IE,TEMP,PRES,W,U,V,DENS,QV,PBL,SFC_DENS,SFC_PRES,SFLX_rad_dn,cosSZA,SFLX_water,SFLX_ENGI) &
1327  !$omp shared(OCN_ATM_TEMP,OCN_ATM_PRES,OCN_ATM_W,OCN_ATM_U,OCN_ATM_V,OCN_ATM_DENS,OCN_ATM_QV) &
1328  !$omp shared(OCN_ATM_PBL,OCN_ATM_SFC_DENS,OCN_ATM_SFC_PRES,OCN_ATM_SFLX_rad_dn,OCN_ATM_cosSZA,OCN_ATM_SFLX_water,OCN_ATM_SFLX_ENGI)
1329  !$acc kernels
1330  do j = js, je
1331  do i = is, ie
1332  temp(i,j) = ocn_atm_temp(i,j)
1333  pres(i,j) = ocn_atm_pres(i,j)
1334  w(i,j) = ocn_atm_w(i,j)
1335  u(i,j) = ocn_atm_u(i,j)
1336  v(i,j) = ocn_atm_v(i,j)
1337  dens(i,j) = ocn_atm_dens(i,j)
1338  qv(i,j) = ocn_atm_qv(i,j)
1339  pbl(i,j) = ocn_atm_pbl(i,j)
1340  sfc_dens(i,j) = ocn_atm_sfc_dens(i,j)
1341  sfc_pres(i,j) = ocn_atm_sfc_pres(i,j)
1342  cossza(i,j) = ocn_atm_cossza(i,j)
1343  sflx_water(i,j) = ocn_atm_sflx_water(i,j)
1344  sflx_engi(i,j) = ocn_atm_sflx_engi(i,j)
1345  do irgn = i_r_ir, i_r_vis
1346  do idir = i_r_direct, i_r_diffuse
1347  sflx_rad_dn(i,j,idir,irgn) = ocn_atm_sflx_rad_dn(i,j,idir,irgn)
1348  enddo
1349  enddo
1350  enddo
1351  enddo
1352  !$acc end kernels
1353 
1354  cnt_putatm_ocn = 0.0_rp
1355 
1356  return
1357  end subroutine cpl_getatm_ocn
1358 
1359  !-----------------------------------------------------------------------------
1360  subroutine cpl_getatm_lnd( &
1361  TEMP, &
1362  PRES, &
1363  W, &
1364  U, &
1365  V, &
1366  DENS, &
1367  QV, &
1368  PBL, &
1369  SFC_DENS, &
1370  SFC_PRES, &
1371  SFLX_rad_dn, &
1372  cosSZA, &
1373  SFLX_water, &
1374  SFLX_ENGI )
1375  implicit none
1376 
1377  real(rp), intent(out) :: temp (ia,ja)
1378  real(rp), intent(out) :: pres (ia,ja)
1379  real(rp), intent(out) :: w (ia,ja)
1380  real(rp), intent(out) :: u (ia,ja)
1381  real(rp), intent(out) :: v (ia,ja)
1382  real(rp), intent(out) :: dens (ia,ja)
1383  real(rp), intent(out) :: qv (ia,ja)
1384  real(rp), intent(out) :: pbl (ia,ja)
1385  real(rp), intent(out) :: sfc_dens (ia,ja)
1386  real(rp), intent(out) :: sfc_pres (ia,ja)
1387  real(rp), intent(out) :: sflx_rad_dn(ia,ja,n_rad_dir,n_rad_rgn)
1388  real(rp), intent(out) :: cossza (ia,ja)
1389  real(rp), intent(out) :: sflx_water (ia,ja)
1390  real(rp), intent(out) :: sflx_engi (ia,ja)
1391 
1392  integer :: i, j, idir, irgn
1393  !---------------------------------------------------------------------------
1394 
1395 !OCL XFILL
1396  !$omp parallel do default(none) private(i,j,idir,irgn) OMP_SCHEDULE_ &
1397  !$omp shared(JS,JE,IS,IE,TEMP,PRES,W,U,V,DENS,QV,PBL,SFC_DENS,SFC_PRES,SFLX_rad_dn,cosSZA,SFLX_water,SFLX_ENGI) &
1398  !$omp shared(LND_ATM_TEMP,LND_ATM_PRES,LND_ATM_W,LND_ATM_U,LND_ATM_V,LND_ATM_DENS,LND_ATM_QV) &
1399  !$omp shared(LND_ATM_PBL,LND_ATM_SFC_DENS,LND_ATM_SFC_PRES,LND_ATM_SFLX_rad_dn,LND_ATM_cosSZA,LND_ATM_SFLX_water,LND_ATM_SFLX_ENGI)
1400  !$acc kernels
1401  do j = js, je
1402  do i = is, ie
1403  temp(i,j) = lnd_atm_temp(i,j)
1404  pres(i,j) = lnd_atm_pres(i,j)
1405  w(i,j) = lnd_atm_w(i,j)
1406  u(i,j) = lnd_atm_u(i,j)
1407  v(i,j) = lnd_atm_v(i,j)
1408  dens(i,j) = lnd_atm_dens(i,j)
1409  qv(i,j) = lnd_atm_qv(i,j)
1410  pbl(i,j) = lnd_atm_pbl(i,j)
1411  sfc_dens(i,j) = lnd_atm_sfc_dens(i,j)
1412  sfc_pres(i,j) = lnd_atm_sfc_pres(i,j)
1413  cossza(i,j) = lnd_atm_cossza(i,j)
1414  sflx_water(i,j) = lnd_atm_sflx_water(i,j)
1415  sflx_engi(i,j) = lnd_atm_sflx_engi(i,j)
1416  do irgn = i_r_ir, i_r_vis
1417  do idir = i_r_direct, i_r_diffuse
1418  sflx_rad_dn(i,j,idir,irgn) = lnd_atm_sflx_rad_dn(i,j,idir,irgn)
1419  enddo
1420  enddo
1421  enddo
1422  enddo
1423  !$acc end kernels
1424 
1425  cnt_putatm_lnd = 0.0_rp
1426 
1427  return
1428  end subroutine cpl_getatm_lnd
1429 
1430  !-----------------------------------------------------------------------------
1431  subroutine cpl_getatm_urb( &
1432  TEMP, &
1433  PRES, &
1434  W, &
1435  U, &
1436  V, &
1437  DENS, &
1438  QV, &
1439  PBL, &
1440  SFC_DENS, &
1441  SFC_PRES, &
1442  SFLX_rad_dn, &
1443  cosSZA, &
1444  SFLX_water, &
1445  SFLX_ENGI )
1446  implicit none
1447 
1448  real(rp), intent(out) :: temp (ia,ja)
1449  real(rp), intent(out) :: pres (ia,ja)
1450  real(rp), intent(out) :: w (ia,ja)
1451  real(rp), intent(out) :: u (ia,ja)
1452  real(rp), intent(out) :: v (ia,ja)
1453  real(rp), intent(out) :: dens (ia,ja)
1454  real(rp), intent(out) :: qv (ia,ja)
1455  real(rp), intent(out) :: pbl (ia,ja)
1456  real(rp), intent(out) :: sfc_dens (ia,ja)
1457  real(rp), intent(out) :: sfc_pres (ia,ja)
1458  real(rp), intent(out) :: sflx_rad_dn(ia,ja,n_rad_dir,n_rad_rgn)
1459  real(rp), intent(out) :: cossza (ia,ja)
1460  real(rp), intent(out) :: sflx_water (ia,ja)
1461  real(rp), intent(out) :: sflx_engi (ia,ja)
1462 
1463  integer :: i, j, idir, irgn
1464  !---------------------------------------------------------------------------
1465 
1466 !OCL XFILL
1467  !$omp parallel do default(none) private(i,j,idir,irgn) OMP_SCHEDULE_ &
1468  !$omp shared(JS,JE,IS,IE,TEMP,PRES,W,U,V,DENS,QV,PBL,SFC_DENS,SFC_PRES,SFLX_rad_dn,cosSZA,SFLX_water,SFLX_ENGI) &
1469  !$omp shared(URB_ATM_TEMP,URB_ATM_PRES,URB_ATM_W,URB_ATM_U,URB_ATM_V,URB_ATM_DENS,URB_ATM_QV) &
1470  !$omp shared(URB_ATM_PBL,URB_ATM_SFC_DENS,URB_ATM_SFC_PRES,URB_ATM_SFLX_rad_dn,URB_ATM_cosSZA,URB_ATM_SFLX_water,URB_ATM_SFLX_ENGI)
1471  !$acc kernels
1472  do j = js, je
1473  do i = is, ie
1474  temp(i,j) = urb_atm_temp(i,j)
1475  pres(i,j) = urb_atm_pres(i,j)
1476  w(i,j) = urb_atm_w(i,j)
1477  u(i,j) = urb_atm_u(i,j)
1478  v(i,j) = urb_atm_v(i,j)
1479  dens(i,j) = urb_atm_dens(i,j)
1480  qv(i,j) = urb_atm_qv(i,j)
1481  pbl(i,j) = urb_atm_pbl(i,j)
1482  sfc_dens(i,j) = urb_atm_sfc_dens(i,j)
1483  sfc_pres(i,j) = urb_atm_sfc_pres(i,j)
1484  cossza(i,j) = urb_atm_cossza(i,j)
1485  sflx_water(i,j) = urb_atm_sflx_water(i,j)
1486  sflx_engi(i,j) = urb_atm_sflx_engi(i,j)
1487  do irgn = i_r_ir, i_r_vis
1488  do idir = i_r_direct, i_r_diffuse
1489  sflx_rad_dn(i,j,idir,irgn) = urb_atm_sflx_rad_dn(i,j,idir,irgn)
1490  enddo
1491  enddo
1492  enddo
1493  enddo
1494  !$acc end kernels
1495 
1496  cnt_putatm_urb = 0.0_rp
1497 
1498  return
1499  end subroutine cpl_getatm_urb
1500 
1501 end module mod_cpl_vars
module COUPLER Variables
real(rp), dimension(:,:), allocatable, public ocn_sfc_z0e
real(rp), dimension(:,:), allocatable, public urb_u10
real(rp), dimension(:,:), allocatable, public ocn_atm_qv
real(rp), dimension(:,:), allocatable, public lnd_sflx_mw
real(rp), dimension(:,:), allocatable, public lnd_sflx_gh
real(rp), dimension(:,:), allocatable, public ocn_atm_pbl
real(rp), dimension(:,:), allocatable, public ocn_sflx_sh
real(rp), dimension(:,:), allocatable, public lnd_sfc_temp
subroutine, public cpl_puturb(SFC_TEMP, SFC_albedo, SFC_Z0M, SFC_Z0H, SFC_Z0E, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_SHEX, SFLX_LHEX, SFLX_QVEX, SFLX_GH, SFLX_QTRC, U10, V10, T2, Q2, mask, countup)
real(rp), dimension(:,:), allocatable, public lnd_sflx_sh
real(rp), dimension(:,:), allocatable, public lnd_q2
real(rp), dimension(:,:), allocatable, public ocn_sflx_lh
real(rp), public cnt_putocn
real(rp), public cnt_putatm_urb
real(rp), dimension(:,:), allocatable, public urb_atm_sflx_water
real(rp), dimension(:,:), allocatable, public urb_atm_pres
real(rp), dimension(:,:), allocatable, public urb_sflx_engi
real(rp), dimension(:,:), allocatable, public ocn_v10
real(rp), dimension(:,:), allocatable, public ocn_atm_w
real(rp), public cnt_puturb
real(rp), dimension(:,:), allocatable, public ocn_t2
real(rp), dimension(:,:), allocatable, public lnd_atm_w
real(rp), dimension(:,:,:), allocatable, public lnd_sflx_qtrc
real(rp), dimension(:,:,:,:), allocatable, public ocn_sfc_albedo
real(rp), public cnt_putlnd
real(rp), dimension(:,:), allocatable, public ocn_sfc_z0h
real(rp), dimension(:,:), allocatable, public lnd_sflx_mv
real(rp), dimension(:,:), allocatable, public ocn_atm_sflx_water
real(rp), dimension(:,:), allocatable, public lnd_atm_sflx_water
real(rp), dimension(:,:), allocatable, public lnd_atm_v
real(rp), dimension(:,:), allocatable, public urb_atm_dens
real(rp), dimension(:,:), allocatable, public urb_v10
real(rp), dimension(:,:), allocatable, public ocn_sflx_mv
real(rp), dimension(:,:), allocatable, public urb_sfc_z0e
real(rp), dimension(:,:), allocatable, public ocn_sflx_engi
real(rp), dimension(:,:), allocatable, public urb_sfc_temp
real(rp), dimension(:,:), allocatable, public urb_sflx_mv
real(rp), dimension(:,:,:), allocatable, public urb_sflx_qtrc
real(rp), dimension(:,:), allocatable, public lnd_atm_dens
real(rp), dimension(:,:), allocatable, public ocn_atm_sfc_dens
real(rp), public cnt_putatm_ocn
real(rp), dimension(:,:), allocatable, public urb_t2
real(rp), dimension(:,:), allocatable, public lnd_sfc_z0h
real(rp), dimension(:,:), allocatable, public urb_atm_sfc_pres
real(rp), dimension(:,:), allocatable, public ocn_sflx_mu
real(rp), dimension(:,:,:,:), allocatable, public ocn_atm_sflx_rad_dn
real(rp), dimension(:,:), allocatable, public urb_sflx_shex
real(rp), dimension(:,:), allocatable, public ocn_atm_v
real(rp), dimension(:,:,:), allocatable, public ocn_sflx_qtrc
subroutine, public cpl_putatm(TEMP, PRES, W, U, V, DENS, QV, PBL, SFC_DENS, SFC_PRES, SFLX_rad_dn, cosSZA, SFLX_water, SFLX_ENGI, countup)
subroutine, public cpl_vars_setup
Setup.
real(rp), dimension(:,:), allocatable, public urb_atm_v
real(rp), dimension(:,:), allocatable, public lnd_sflx_lh
real(rp), dimension(:,:), allocatable, public urb_sflx_gh
real(rp), dimension(:,:,:,:), allocatable, public lnd_atm_sflx_rad_dn
subroutine, public cpl_getatm_ocn(TEMP, PRES, W, U, V, DENS, QV, PBL, SFC_DENS, SFC_PRES, SFLX_rad_dn, cosSZA, SFLX_water, SFLX_ENGI)
real(rp), dimension(:,:), allocatable, public lnd_v10
real(rp), dimension(:,:), allocatable, public urb_sflx_lhex
real(rp), dimension(:,:), allocatable, public urb_atm_temp
real(rp), dimension(:,:), allocatable, public lnd_sflx_engi
real(rp), dimension(:,:), allocatable, public urb_sflx_mu
real(rp), dimension(:,:), allocatable, public urb_sflx_lh
subroutine, public cpl_getsfc_atm(SFC_TEMP, SFC_albedo, SFC_Z0M, SFC_Z0H, SFC_Z0E, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_SHEX, SFLX_LHEX, SFLX_QVEX, SFLX_GH, SFLX_QTRC, SFLX_ENGI, U10, V10, T2, Q2)
real(rp), dimension(:,:), allocatable, public lnd_sfc_z0e
real(rp), dimension(:,:), allocatable, public urb_sflx_mw
real(rp), dimension(:,:), allocatable, public ocn_atm_pres
real(rp), dimension(:,:), allocatable, public ocn_atm_sfc_pres
real(rp), dimension(:,:), allocatable, public ocn_sfc_z0m
real(rp), dimension(:,:), allocatable, public ocn_atm_u
real(rp), dimension(:,:), allocatable, public urb_atm_cossza
real(rp), dimension(:,:), allocatable, public lnd_atm_qv
subroutine, public cpl_putocn(SFC_TEMP, SFC_albedo, SFC_Z0M, SFC_Z0H, SFC_Z0E, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_GH, SFLX_QTRC, U10, V10, T2, Q2, mask, countup)
subroutine, public cpl_putlnd(SFC_TEMP, SFC_albedo, SFC_Z0M, SFC_Z0H, SFC_Z0E, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_GH, SFLX_QTRC, U10, V10, T2, Q2, mask, countup)
real(rp), dimension(:,:), allocatable, public lnd_atm_pbl
real(rp), dimension(:,:), allocatable, public lnd_atm_sfc_dens
real(rp), dimension(:,:), allocatable, public urb_q2
subroutine, public cpl_vars_finalize
Finalize.
real(rp), dimension(:,:), allocatable, public ocn_sflx_mw
subroutine, public cpl_getatm_urb(TEMP, PRES, W, U, V, DENS, QV, PBL, SFC_DENS, SFC_PRES, SFLX_rad_dn, cosSZA, SFLX_water, SFLX_ENGI)
real(rp), dimension(:,:), allocatable, public lnd_atm_pres
real(rp), dimension(:,:), allocatable, public ocn_atm_temp
real(rp), dimension(:,:), allocatable, public lnd_atm_sfc_pres
real(rp), dimension(:,:), allocatable, public lnd_t2
real(rp), dimension(:,:), allocatable, public lnd_atm_sflx_engi
real(rp), dimension(:,:), allocatable, public urb_atm_w
real(rp), dimension(:,:,:,:), allocatable, public urb_sfc_albedo
real(rp), dimension(:,:), allocatable, public urb_sflx_sh
real(rp), dimension(:,:), allocatable, public lnd_atm_u
real(rp), dimension(:,:,:,:), allocatable, public lnd_sfc_albedo
real(rp), dimension(:,:), allocatable, public urb_atm_qv
real(rp), dimension(:,:), allocatable, public urb_sfc_z0m
subroutine, public cpl_getatm_lnd(TEMP, PRES, W, U, V, DENS, QV, PBL, SFC_DENS, SFC_PRES, SFLX_rad_dn, cosSZA, SFLX_water, SFLX_ENGI)
real(rp), dimension(:,:), allocatable, public ocn_q2
real(rp), dimension(:,:), allocatable, public urb_sflx_qvex
real(rp), dimension(:,:), allocatable, public lnd_atm_temp
real(rp), dimension(:,:), allocatable, public lnd_atm_cossza
real(rp), dimension(:,:), allocatable, public ocn_atm_sflx_engi
real(rp), public cnt_putatm_lnd
real(rp), dimension(:,:), allocatable, public urb_atm_sfc_dens
real(rp), dimension(:,:), allocatable, public lnd_sflx_mu
real(rp), dimension(:,:), allocatable, public ocn_atm_cossza
real(rp), dimension(:,:,:,:), allocatable, public urb_atm_sflx_rad_dn
real(rp), dimension(:,:), allocatable, public ocn_sfc_temp
real(rp), dimension(:,:), allocatable, public urb_atm_pbl
real(rp), dimension(:,:), allocatable, public lnd_sfc_z0m
real(rp), dimension(:,:), allocatable, public urb_atm_sflx_engi
real(rp), dimension(:,:), allocatable, public urb_atm_u
real(rp), dimension(:,:), allocatable, public ocn_u10
real(rp), dimension(:,:), allocatable, public lnd_u10
real(rp), dimension(:,:), allocatable, public ocn_atm_dens
real(rp), dimension(:,:), allocatable, public ocn_sflx_gh
real(rp), dimension(:,:), allocatable, public urb_sfc_z0h
module Lake admin
logical, public lake_do
module Land admin
logical, public land_do
module Ocean admin
logical, public ocean_do
module Urban admin
logical, public urban_do
module atmosphere / grid / cartesC index
integer, public ie
end point of inner domain: x, local
integer, public je
end point of inner domain: y, local
integer, public is
start point of inner domain: x, local
integer, public js
start point of inner domain: y, local
module atmosphere / hydrometeor
real(rp), public lhf
latent heat of fusion for use [J/kg]
real(rp), public cv_ice
CV for ice [J/kg/K].
real(rp), public cv_water
CV for water [J/kg/K].
module CONSTANT
Definition: scale_const.F90:11
real(rp), public const_eps
small number
Definition: scale_const.F90:35
real(rp), public const_undef
Definition: scale_const.F90:43
module coupler / surface-atmospehre
integer, parameter, public i_r_ir
integer, parameter, public i_r_direct
integer, parameter, public n_rad_rgn
integer, parameter, public i_r_diffuse
integer, parameter, public n_rad_dir
integer, parameter, public i_r_vis
module DEBUG
Definition: scale_debug.F90:11
module STDIO
Definition: scale_io.F90:10
module LANDUSE
real(rp), dimension(:,:), allocatable, public landuse_fact_lake
lake factor
real(rp), dimension(:,:), allocatable, public landuse_fact_urban
urban factor
real(rp), dimension(:,:), allocatable, public landuse_fact_ocean
ocean factor
real(rp), dimension(:,:), allocatable, public landuse_fact_land
land factor
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
module PRECISION
integer, parameter, public rp
module profiler
Definition: scale_prof.F90:11
module TRACER
real(rp), dimension(qa_max), public tracer_engi0
real(rp), dimension(qa_max), public tracer_cv
integer, public qa