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  !$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)
409  !$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_GH,URB_SFLX_QTRC,URB_SFLX_ENGI,URB_U10,URB_V10,URB_T2,URB_Q2)
410  !$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)
411  !$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)
412  allocate( urb_atm_temp(ia,ja) )
413  allocate( urb_atm_pres(ia,ja) )
414  allocate( urb_atm_w(ia,ja) )
415  allocate( urb_atm_u(ia,ja) )
416  allocate( urb_atm_v(ia,ja) )
417  allocate( urb_atm_dens(ia,ja) )
418  allocate( urb_atm_qv(ia,ja) )
419  allocate( urb_atm_pbl(ia,ja) )
420  allocate( urb_atm_sfc_dens(ia,ja) )
421  allocate( urb_atm_sfc_pres(ia,ja) )
423  allocate( urb_atm_cossza(ia,ja) )
424  allocate( urb_atm_sflx_water(ia,ja) )
425  allocate( urb_atm_sflx_engi(ia,ja) )
426  urb_atm_temp(:,:) = undef
427  urb_atm_pres(:,:) = undef
428  urb_atm_w(:,:) = undef
429  urb_atm_u(:,:) = undef
430  urb_atm_v(:,:) = undef
431  urb_atm_dens(:,:) = undef
432  urb_atm_qv(:,:) = undef
433  urb_atm_pbl(:,:) = undef
434  urb_atm_sfc_dens(:,:) = undef
435  urb_atm_sfc_pres(:,:) = undef
436  urb_atm_sflx_rad_dn(:,:,:,:) = undef
437  urb_atm_cossza(:,:) = undef
438  urb_atm_sflx_water(:,:) = undef
439  urb_atm_sflx_engi(:,:) = undef
440  !$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)
441 
442  ! counter intialize
443  cnt_putatm_ocn = 0.0_rp
444  cnt_putatm_lnd = 0.0_rp
445  cnt_putatm_urb = 0.0_rp
446  cnt_putocn = 0.0_rp
447  cnt_putlnd = 0.0_rp
448  cnt_puturb = 0.0_rp
449 
450  if ( atmos_hydrometeor_dry ) then
451  ocn_atm_qv = 0.0_rp
452  lnd_atm_qv = 0.0_rp
453  urb_atm_qv = 0.0_rp
454  endif
455 
456  return
457  end subroutine cpl_vars_setup
458 
459  !-----------------------------------------------------------------------------
461  subroutine cpl_vars_finalize
462  implicit none
463  !---------------------------------------------------------------------------
464 
465  log_newline
466  log_info("CPL_vars_finalize",*) 'Finalize'
467 
468  !$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)
469  deallocate( ocn_sfc_temp )
470  deallocate( ocn_sfc_albedo )
471  deallocate( ocn_sfc_z0m )
472  deallocate( ocn_sfc_z0h )
473  deallocate( ocn_sfc_z0e )
474  deallocate( ocn_sflx_mu )
475  deallocate( ocn_sflx_mv )
476  deallocate( ocn_sflx_mw )
477  deallocate( ocn_sflx_sh )
478  deallocate( ocn_sflx_lh )
479  deallocate( ocn_sflx_gh )
480  deallocate( ocn_sflx_qtrc )
481  deallocate( ocn_sflx_engi )
482  deallocate( ocn_u10 )
483  deallocate( ocn_v10 )
484  deallocate( ocn_t2 )
485  deallocate( ocn_q2 )
486 
487  !$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)
488  deallocate( lnd_sfc_temp )
489  deallocate( lnd_sfc_albedo )
490  deallocate( lnd_sfc_z0m )
491  deallocate( lnd_sfc_z0h )
492  deallocate( lnd_sfc_z0e )
493  deallocate( lnd_sflx_mu )
494  deallocate( lnd_sflx_mv )
495  deallocate( lnd_sflx_mw )
496  deallocate( lnd_sflx_sh )
497  deallocate( lnd_sflx_lh )
498  deallocate( lnd_sflx_gh )
499  deallocate( lnd_sflx_qtrc )
500  deallocate( lnd_sflx_engi )
501  deallocate( lnd_u10 )
502  deallocate( lnd_v10 )
503  deallocate( lnd_t2 )
504  deallocate( lnd_q2 )
505 
506  !$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)
507  deallocate( urb_sfc_temp )
508  deallocate( urb_sfc_albedo )
509  deallocate( urb_sfc_z0m )
510  deallocate( urb_sfc_z0h )
511  deallocate( urb_sfc_z0e )
512  deallocate( urb_sflx_mu )
513  deallocate( urb_sflx_mv )
514  deallocate( urb_sflx_mw )
515  deallocate( urb_sflx_sh )
516  deallocate( urb_sflx_lh )
517  deallocate( urb_sflx_shex )
518  deallocate( urb_sflx_lhex )
519  deallocate( urb_sflx_qvex )
520  deallocate( urb_sflx_gh )
521  deallocate( urb_sflx_qtrc )
522  deallocate( urb_sflx_engi )
523  deallocate( urb_u10 )
524  deallocate( urb_v10 )
525  deallocate( urb_t2 )
526  deallocate( urb_q2 )
527 
528  !$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)
529  deallocate( ocn_atm_temp )
530  deallocate( ocn_atm_pres )
531  deallocate( ocn_atm_w )
532  deallocate( ocn_atm_u )
533  deallocate( ocn_atm_v )
534  deallocate( ocn_atm_dens )
535  deallocate( ocn_atm_qv )
536  deallocate( ocn_atm_pbl )
537  deallocate( ocn_atm_sfc_dens )
538  deallocate( ocn_atm_sfc_pres )
539  deallocate( ocn_atm_sflx_rad_dn )
540  deallocate( ocn_atm_cossza )
541  deallocate( ocn_atm_sflx_water )
542  deallocate( ocn_atm_sflx_engi )
543 
544  !$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)
545  deallocate( lnd_atm_temp )
546  deallocate( lnd_atm_pres )
547  deallocate( lnd_atm_w )
548  deallocate( lnd_atm_u )
549  deallocate( lnd_atm_v )
550  deallocate( lnd_atm_dens )
551  deallocate( lnd_atm_qv )
552  deallocate( lnd_atm_pbl )
553  deallocate( lnd_atm_sfc_dens )
554  deallocate( lnd_atm_sfc_pres )
555  deallocate( lnd_atm_sflx_rad_dn )
556  deallocate( lnd_atm_cossza )
557  deallocate( lnd_atm_sflx_water )
558  deallocate( lnd_atm_sflx_engi )
559 
560  !$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)
561  deallocate( urb_atm_temp )
562  deallocate( urb_atm_pres )
563  deallocate( urb_atm_w )
564  deallocate( urb_atm_u )
565  deallocate( urb_atm_v )
566  deallocate( urb_atm_dens )
567  deallocate( urb_atm_qv )
568  deallocate( urb_atm_pbl )
569  deallocate( urb_atm_sfc_dens )
570  deallocate( urb_atm_sfc_pres )
571  deallocate( urb_atm_sflx_rad_dn )
572  deallocate( urb_atm_cossza )
573  deallocate( urb_atm_sflx_water )
574  deallocate( urb_atm_sflx_engi )
575 
576  return
577  end subroutine cpl_vars_finalize
578 
579  !-----------------------------------------------------------------------------
580  subroutine cpl_putatm( &
581  TEMP, &
582  PRES, &
583  W, &
584  U, &
585  V, &
586  DENS, &
587  QV, &
588  PBL, &
589  SFC_DENS, &
590  SFC_PRES, &
591  SFLX_rad_dn, &
592  cosSZA, &
593  SFLX_water, &
594  SFLX_ENGI, &
595  countup )
598  cv_water, &
599  cv_ice, &
600  lhf
601  implicit none
602 
603  real(rp), intent(in) :: temp (ia,ja)
604  real(rp), intent(in) :: pres (ia,ja)
605  real(rp), intent(in) :: w (ia,ja)
606  real(rp), intent(in) :: u (ia,ja)
607  real(rp), intent(in) :: v (ia,ja)
608  real(rp), intent(in) :: dens (ia,ja)
609  real(rp), intent(in) :: qv (ia,ja)
610  real(rp), intent(in) :: pbl (ia,ja)
611  real(rp), intent(in) :: sfc_dens (ia,ja)
612  real(rp), intent(in) :: sfc_pres (ia,ja)
613  real(rp), intent(in) :: sflx_rad_dn(ia,ja,n_rad_dir,n_rad_rgn)
614  real(rp), intent(in) :: cossza (ia,ja)
615  real(rp), intent(in) :: sflx_water (ia,ja)
616  real(rp), intent(in) :: sflx_engi (ia,ja)
617  logical, intent(in) :: countup
618 
619  integer :: i, j, idir, irgn
620  !---------------------------------------------------------------------------
621 
622  !$omp parallel do default(none) private(i,j,idir,irgn) OMP_SCHEDULE_ &
623  !$omp shared(JS,JE,IS,IE,OCN_ATM_TEMP,OCN_ATM_PRES,OCN_ATM_W,OCN_ATM_U) &
624  !$omp shared(OCN_ATM_V,OCN_ATM_DENS,CNT_putATM_OCN,TEMP,PRES,W,U,V,DENS) &
625  !$omp shared(ATMOS_HYDROMETEOR_dry,OCN_ATM_QV,OCN_ATM_PBL,OCN_ATM_SFC_DENS,OCN_ATM_SFC_PRES,OCN_ATM_SFLX_rad_dn) &
626  !$omp shared(OCN_ATM_cosSZA,OCN_ATM_SFLX_water,OCN_ATM_SFLX_ENGI,QV,PBL,SFC_DENS,SFC_PRES) &
627  !$omp shared(SFLX_rad_dn,cosSZA,SFLX_water,SFLX_ENGI) &
628  !$omp shared(LND_ATM_TEMP,LND_ATM_PRES,LND_ATM_W,LND_ATM_U,LND_ATM_V,LND_ATM_DENS) &
629  !$omp shared(LND_ATM_QV,LND_ATM_PBL,LND_ATM_SFC_DENS,LND_ATM_SFC_PRES,LND_ATM_SFLX_rad_dn,LND_ATM_cosSZA) &
630  !$omp shared(LND_ATM_SFLX_water,LND_ATM_SFLX_ENGI) &
631  !$omp shared(URB_ATM_TEMP,URB_ATM_PRES,URB_ATM_W,URB_ATM_U,URB_ATM_V,URB_ATM_DENS) &
632  !$omp shared(URB_ATM_QV,URB_ATM_PBL,URB_ATM_SFC_DENS,URB_ATM_SFC_PRES,URB_ATM_SFLX_rad_dn,URB_ATM_cosSZA) &
633  !$omp shared(URB_ATM_SFLX_water,URB_ATM_SFLX_ENGI,CNT_putATM_URB,CNT_putATM_LND)
634  !$acc kernels
635  do j = js, je
636  do i = is, ie
637  ! for ocean
638  ocn_atm_temp(i,j) = ocn_atm_temp(i,j) * cnt_putatm_ocn + temp(i,j)
639  ocn_atm_pres(i,j) = ocn_atm_pres(i,j) * cnt_putatm_ocn + pres(i,j)
640  ocn_atm_w(i,j) = ocn_atm_w(i,j) * cnt_putatm_ocn + w(i,j)
641  ocn_atm_u(i,j) = ocn_atm_u(i,j) * cnt_putatm_ocn + u(i,j)
642  ocn_atm_v(i,j) = ocn_atm_v(i,j) * cnt_putatm_ocn + v(i,j)
643  ocn_atm_dens(i,j) = ocn_atm_dens(i,j) * cnt_putatm_ocn + dens(i,j)
644  if ( .NOT. atmos_hydrometeor_dry ) &
645  ocn_atm_qv(i,j) = ocn_atm_qv(i,j) * cnt_putatm_ocn + qv(i,j)
646  ocn_atm_pbl(i,j) = ocn_atm_pbl(i,j) * cnt_putatm_ocn + pbl(i,j)
647  ocn_atm_sfc_dens(i,j) = ocn_atm_sfc_dens(i,j) * cnt_putatm_ocn + sfc_dens(i,j)
648  ocn_atm_sfc_pres(i,j) = ocn_atm_sfc_pres(i,j) * cnt_putatm_ocn + sfc_pres(i,j)
649  ocn_atm_cossza(i,j) = ocn_atm_cossza(i,j) * cnt_putatm_ocn + cossza(i,j)
650  ocn_atm_sflx_water(i,j) = ocn_atm_sflx_water(i,j) * cnt_putatm_ocn + sflx_water(i,j)
651  ocn_atm_sflx_engi(i,j) = ocn_atm_sflx_engi(i,j) * cnt_putatm_ocn + sflx_engi(i,j)
652  do irgn = i_r_ir, i_r_vis
653  do idir = i_r_direct, i_r_diffuse
654  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)
655  enddo
656  enddo
657 
658  ! for land
659  lnd_atm_temp(i,j) = lnd_atm_temp(i,j) * cnt_putatm_lnd + temp(i,j)
660  lnd_atm_pres(i,j) = lnd_atm_pres(i,j) * cnt_putatm_lnd + pres(i,j)
661  lnd_atm_w(i,j) = lnd_atm_w(i,j) * cnt_putatm_lnd + w(i,j)
662  lnd_atm_u(i,j) = lnd_atm_u(i,j) * cnt_putatm_lnd + u(i,j)
663  lnd_atm_v(i,j) = lnd_atm_v(i,j) * cnt_putatm_lnd + v(i,j)
664  lnd_atm_dens(i,j) = lnd_atm_dens(i,j) * cnt_putatm_lnd + dens(i,j)
665  if ( .NOT. atmos_hydrometeor_dry ) &
666  lnd_atm_qv(i,j) = lnd_atm_qv(i,j) * cnt_putatm_lnd + qv(i,j)
667  lnd_atm_pbl(i,j) = lnd_atm_pbl(i,j) * cnt_putatm_lnd + pbl(i,j)
668  lnd_atm_sfc_dens(i,j) = lnd_atm_sfc_dens(i,j) * cnt_putatm_lnd + sfc_dens(i,j)
669  lnd_atm_sfc_pres(i,j) = lnd_atm_sfc_pres(i,j) * cnt_putatm_lnd + sfc_pres(i,j)
670  lnd_atm_cossza(i,j) = lnd_atm_cossza(i,j) * cnt_putatm_lnd + cossza(i,j)
671  lnd_atm_sflx_water(i,j) = lnd_atm_sflx_water(i,j) * cnt_putatm_lnd + sflx_water(i,j)
672  lnd_atm_sflx_engi(i,j) = lnd_atm_sflx_engi(i,j) * cnt_putatm_lnd + sflx_engi(i,j)
673  do irgn = i_r_ir, i_r_vis
674  do idir = i_r_direct, i_r_diffuse
675  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)
676  enddo
677  enddo
678 
679  ! for urban
680  urb_atm_temp(i,j) = urb_atm_temp(i,j) * cnt_putatm_urb + temp(i,j)
681  urb_atm_pres(i,j) = urb_atm_pres(i,j) * cnt_putatm_urb + pres(i,j)
682  urb_atm_w(i,j) = urb_atm_w(i,j) * cnt_putatm_urb + w(i,j)
683  urb_atm_u(i,j) = urb_atm_u(i,j) * cnt_putatm_urb + u(i,j)
684  urb_atm_v(i,j) = urb_atm_v(i,j) * cnt_putatm_urb + v(i,j)
685  urb_atm_dens(i,j) = urb_atm_dens(i,j) * cnt_putatm_urb + dens(i,j)
686  if ( .NOT. atmos_hydrometeor_dry ) &
687  urb_atm_qv(i,j) = urb_atm_qv(i,j) * cnt_putatm_urb + qv(i,j)
688  urb_atm_pbl(i,j) = urb_atm_pbl(i,j) * cnt_putatm_urb + pbl(i,j)
689  urb_atm_sfc_dens(i,j) = urb_atm_sfc_dens(i,j) * cnt_putatm_urb + sfc_dens(i,j)
690  urb_atm_sfc_pres(i,j) = urb_atm_sfc_pres(i,j) * cnt_putatm_urb + sfc_pres(i,j)
691  urb_atm_cossza(i,j) = urb_atm_cossza(i,j) * cnt_putatm_urb + cossza(i,j)
692  urb_atm_sflx_water(i,j) = urb_atm_sflx_water(i,j) * cnt_putatm_urb + sflx_water(i,j)
693  urb_atm_sflx_engi(i,j) = urb_atm_sflx_engi(i,j) * cnt_putatm_urb + sflx_engi(i,j)
694  do irgn = i_r_ir, i_r_vis
695  do idir = i_r_direct, i_r_diffuse
696  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)
697  enddo
698  enddo
699 
700 
701  ! for ocean
702  ocn_atm_temp(i,j) = ocn_atm_temp(i,j) / ( cnt_putatm_ocn + 1.0_rp )
703  ocn_atm_pres(i,j) = ocn_atm_pres(i,j) / ( cnt_putatm_ocn + 1.0_rp )
704  ocn_atm_w(i,j) = ocn_atm_w(i,j) / ( cnt_putatm_ocn + 1.0_rp )
705  ocn_atm_u(i,j) = ocn_atm_u(i,j) / ( cnt_putatm_ocn + 1.0_rp )
706  ocn_atm_v(i,j) = ocn_atm_v(i,j) / ( cnt_putatm_ocn + 1.0_rp )
707  ocn_atm_dens(i,j) = ocn_atm_dens(i,j) / ( cnt_putatm_ocn + 1.0_rp )
708  ocn_atm_qv(i,j) = ocn_atm_qv(i,j) / ( cnt_putatm_ocn + 1.0_rp )
709  ocn_atm_pbl(i,j) = ocn_atm_pbl(i,j) / ( cnt_putatm_ocn + 1.0_rp )
710  ocn_atm_sfc_dens(i,j) = ocn_atm_sfc_dens(i,j) / ( cnt_putatm_ocn + 1.0_rp )
711  ocn_atm_sfc_pres(i,j) = ocn_atm_sfc_pres(i,j) / ( cnt_putatm_ocn + 1.0_rp )
712  ocn_atm_cossza(i,j) = ocn_atm_cossza(i,j) / ( cnt_putatm_ocn + 1.0_rp )
713  ocn_atm_sflx_water(i,j) = ocn_atm_sflx_water(i,j) / ( cnt_putatm_ocn + 1.0_rp )
714  ocn_atm_sflx_engi(i,j) = ocn_atm_sflx_engi(i,j) / ( cnt_putatm_ocn + 1.0_rp )
715  do irgn = i_r_ir, i_r_vis
716  do idir = i_r_direct, i_r_diffuse
717  ocn_atm_sflx_rad_dn(i,j,idir,irgn) = ocn_atm_sflx_rad_dn(i,j,idir,irgn) / ( cnt_putatm_ocn + 1.0_rp )
718  enddo
719  enddo
720 
721  ! for land
722  lnd_atm_temp(i,j) = lnd_atm_temp(i,j) / ( cnt_putatm_lnd + 1.0_rp )
723  lnd_atm_pres(i,j) = lnd_atm_pres(i,j) / ( cnt_putatm_lnd + 1.0_rp )
724  lnd_atm_w(i,j) = lnd_atm_w(i,j) / ( cnt_putatm_lnd + 1.0_rp )
725  lnd_atm_u(i,j) = lnd_atm_u(i,j) / ( cnt_putatm_lnd + 1.0_rp )
726  lnd_atm_v(i,j) = lnd_atm_v(i,j) / ( cnt_putatm_lnd + 1.0_rp )
727  lnd_atm_dens(i,j) = lnd_atm_dens(i,j) / ( cnt_putatm_lnd + 1.0_rp )
728  lnd_atm_qv(i,j) = lnd_atm_qv(i,j) / ( cnt_putatm_lnd + 1.0_rp )
729  lnd_atm_pbl(i,j) = lnd_atm_pbl(i,j) / ( cnt_putatm_lnd + 1.0_rp )
730  lnd_atm_sfc_dens(i,j) = lnd_atm_sfc_dens(i,j) / ( cnt_putatm_lnd + 1.0_rp )
731  lnd_atm_sfc_pres(i,j) = lnd_atm_sfc_pres(i,j) / ( cnt_putatm_lnd + 1.0_rp )
732  lnd_atm_cossza(i,j) = lnd_atm_cossza(i,j) / ( cnt_putatm_lnd + 1.0_rp )
733  lnd_atm_sflx_water(i,j) = lnd_atm_sflx_water(i,j) / ( cnt_putatm_lnd + 1.0_rp )
734  lnd_atm_sflx_engi(i,j) = lnd_atm_sflx_engi(i,j) / ( cnt_putatm_lnd + 1.0_rp )
735  do irgn = i_r_ir, i_r_vis
736  do idir = i_r_direct, i_r_diffuse
737  lnd_atm_sflx_rad_dn(i,j,idir,irgn) = lnd_atm_sflx_rad_dn(i,j,idir,irgn) / ( cnt_putatm_lnd + 1.0_rp )
738  enddo
739  enddo
740 
741  ! for urban
742  urb_atm_temp(i,j) = urb_atm_temp(i,j) / ( cnt_putatm_urb + 1.0_rp )
743  urb_atm_pres(i,j) = urb_atm_pres(i,j) / ( cnt_putatm_urb + 1.0_rp )
744  urb_atm_w(i,j) = urb_atm_w(i,j) / ( cnt_putatm_urb + 1.0_rp )
745  urb_atm_u(i,j) = urb_atm_u(i,j) / ( cnt_putatm_urb + 1.0_rp )
746  urb_atm_v(i,j) = urb_atm_v(i,j) / ( cnt_putatm_urb + 1.0_rp )
747  urb_atm_dens(i,j) = urb_atm_dens(i,j) / ( cnt_putatm_urb + 1.0_rp )
748  urb_atm_qv(i,j) = urb_atm_qv(i,j) / ( cnt_putatm_urb + 1.0_rp )
749  urb_atm_pbl(i,j) = urb_atm_pbl(i,j) / ( cnt_putatm_urb + 1.0_rp )
750  urb_atm_sfc_dens(i,j) = urb_atm_sfc_dens(i,j) / ( cnt_putatm_urb + 1.0_rp )
751  urb_atm_sfc_pres(i,j) = urb_atm_sfc_pres(i,j) / ( cnt_putatm_urb + 1.0_rp )
752  urb_atm_cossza(i,j) = urb_atm_cossza(i,j) / ( cnt_putatm_urb + 1.0_rp )
753  urb_atm_sflx_water(i,j) = urb_atm_sflx_water(i,j) / ( cnt_putatm_urb + 1.0_rp )
754  urb_atm_sflx_engi(i,j) = urb_atm_sflx_engi(i,j) / ( cnt_putatm_urb + 1.0_rp )
755  do irgn = i_r_ir, i_r_vis
756  do idir = i_r_direct, i_r_diffuse
757  urb_atm_sflx_rad_dn(i,j,idir,irgn) = urb_atm_sflx_rad_dn(i,j,idir,irgn) / ( cnt_putatm_urb + 1.0_rp )
758  enddo
759  enddo
760  enddo
761  enddo
762  !$acc end kernels
763 
764  if( countup ) then
765  cnt_putatm_ocn = cnt_putatm_ocn + 1.0_rp
766  cnt_putatm_lnd = cnt_putatm_lnd + 1.0_rp
767  cnt_putatm_urb = cnt_putatm_urb + 1.0_rp
768  endif
769 
770  return
771  end subroutine cpl_putatm
772 
773  !-----------------------------------------------------------------------------
774  subroutine cpl_putocn( &
775  SFC_TEMP, &
776  SFC_albedo, &
777  SFC_Z0M, &
778  SFC_Z0H, &
779  SFC_Z0E, &
780  SFLX_MW, &
781  SFLX_MU, &
782  SFLX_MV, &
783  SFLX_SH, &
784  SFLX_LH, &
785  SFLX_GH, &
786  SFLX_QTRC, &
787  U10, &
788  V10, &
789  T2, &
790  Q2, &
791  mask, &
792  countup )
793  implicit none
794 
795  real(rp), intent(in) :: sfc_temp (ia,ja)
796  real(rp), intent(in) :: sfc_albedo(ia,ja,n_rad_dir,n_rad_rgn)
797  real(rp), intent(in) :: sfc_z0m (ia,ja)
798  real(rp), intent(in) :: sfc_z0h (ia,ja)
799  real(rp), intent(in) :: sfc_z0e (ia,ja)
800  real(rp), intent(in) :: sflx_mw (ia,ja)
801  real(rp), intent(in) :: sflx_mu (ia,ja)
802  real(rp), intent(in) :: sflx_mv (ia,ja)
803  real(rp), intent(in) :: sflx_sh (ia,ja)
804  real(rp), intent(in) :: sflx_lh (ia,ja)
805  real(rp), intent(in) :: sflx_gh (ia,ja)
806  real(rp), intent(in) :: sflx_qtrc (ia,ja,qa)
807  real(rp), intent(in) :: u10 (ia,ja)
808  real(rp), intent(in) :: v10 (ia,ja)
809  real(rp), intent(in) :: t2 (ia,ja)
810  real(rp), intent(in) :: q2 (ia,ja)
811  logical, intent(in) :: mask (ia,ja)
812  logical, intent(in) :: countup
813 
814  integer :: i, j, iq, idir, irgn
815  !---------------------------------------------------------------------------
816 
817  !$omp parallel do default(none) private(i,j,idir,irgn) OMP_SCHEDULE_ &
818  !$omp shared(JS,JE,IS,IE,QA) &
819  !$omp shared(OCN_SFC_TEMP,OCN_SFC_albedo,OCN_SFC_Z0M,OCN_SFC_Z0H,OCN_SFC_Z0E) &
820  !$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) &
821  !$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) &
822  !$omp shared(mask) &
823  !$omp shared(CNT_putOCN) &
824  !$omp shared(TRACER_CV,TRACER_ENGI0)
825  !$acc kernels
826  do j = js, je
827  do i = is, ie
828  if ( mask(i,j) ) then
829  ocn_sfc_temp(i,j) = ocn_sfc_temp(i,j) * cnt_putocn + sfc_temp(i,j)
830  ocn_sfc_z0m(i,j) = ocn_sfc_z0m(i,j) * cnt_putocn + sfc_z0m(i,j)
831  ocn_sfc_z0h(i,j) = ocn_sfc_z0h(i,j) * cnt_putocn + sfc_z0h(i,j)
832  ocn_sfc_z0e(i,j) = ocn_sfc_z0e(i,j) * cnt_putocn + sfc_z0e(i,j)
833  ocn_sflx_mw(i,j) = ocn_sflx_mw(i,j) * cnt_putocn + sflx_mw(i,j)
834  ocn_sflx_mu(i,j) = ocn_sflx_mu(i,j) * cnt_putocn + sflx_mu(i,j)
835  ocn_sflx_mv(i,j) = ocn_sflx_mv(i,j) * cnt_putocn + sflx_mv(i,j)
836  ocn_sflx_sh(i,j) = ocn_sflx_sh(i,j) * cnt_putocn + sflx_sh(i,j)
837  ocn_sflx_lh(i,j) = ocn_sflx_lh(i,j) * cnt_putocn + sflx_lh(i,j)
838  ocn_sflx_gh(i,j) = ocn_sflx_gh(i,j) * cnt_putocn + sflx_gh(i,j)
840  do iq = 1, qa
841  ocn_sflx_qtrc(i,j,iq) = ocn_sflx_qtrc(i,j,iq) * cnt_putocn + sflx_qtrc(i,j,iq)
842  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) )
843  end do
844  ocn_u10(i,j) = ocn_u10(i,j) * cnt_putocn + u10(i,j)
845  ocn_v10(i,j) = ocn_v10(i,j) * cnt_putocn + v10(i,j)
846  ocn_t2(i,j) = ocn_t2(i,j) * cnt_putocn + t2(i,j)
847  ocn_q2(i,j) = ocn_q2(i,j) * cnt_putocn + q2(i,j)
848  do irgn = i_r_ir, i_r_vis
849  do idir = i_r_direct, i_r_diffuse
850  ocn_sfc_albedo(i,j,idir,irgn) = ocn_sfc_albedo(i,j,idir,irgn) * cnt_putocn + sfc_albedo(i,j,idir,irgn)
851  enddo
852  enddo
853 
854  ocn_sfc_temp(i,j) = ocn_sfc_temp(i,j) / ( cnt_putocn + 1.0_rp )
855  ocn_sfc_z0m(i,j) = ocn_sfc_z0m(i,j) / ( cnt_putocn + 1.0_rp )
856  ocn_sfc_z0h(i,j) = ocn_sfc_z0h(i,j) / ( cnt_putocn + 1.0_rp )
857  ocn_sfc_z0e(i,j) = ocn_sfc_z0e(i,j) / ( cnt_putocn + 1.0_rp )
858  ocn_sflx_mw(i,j) = ocn_sflx_mw(i,j) / ( cnt_putocn + 1.0_rp )
859  ocn_sflx_mu(i,j) = ocn_sflx_mu(i,j) / ( cnt_putocn + 1.0_rp )
860  ocn_sflx_mv(i,j) = ocn_sflx_mv(i,j) / ( cnt_putocn + 1.0_rp )
861  ocn_sflx_sh(i,j) = ocn_sflx_sh(i,j) / ( cnt_putocn + 1.0_rp )
862  ocn_sflx_lh(i,j) = ocn_sflx_lh(i,j) / ( cnt_putocn + 1.0_rp )
863  ocn_sflx_gh(i,j) = ocn_sflx_gh(i,j) / ( cnt_putocn + 1.0_rp )
864  ocn_sflx_qtrc(i,j,:) = ocn_sflx_qtrc(i,j,:) / ( cnt_putocn + 1.0_rp )
865  ocn_sflx_engi(i,j) = ocn_sflx_engi(i,j) / ( cnt_putocn + 1.0_rp )
866  ocn_u10(i,j) = ocn_u10(i,j) / ( cnt_putocn + 1.0_rp )
867  ocn_v10(i,j) = ocn_v10(i,j) / ( cnt_putocn + 1.0_rp )
868  ocn_t2(i,j) = ocn_t2(i,j) / ( cnt_putocn + 1.0_rp )
869  ocn_q2(i,j) = ocn_q2(i,j) / ( cnt_putocn + 1.0_rp )
870  do irgn = i_r_ir, i_r_vis
871  do idir = i_r_direct, i_r_diffuse
872  ocn_sfc_albedo(i,j,idir,irgn) = ocn_sfc_albedo(i,j,idir,irgn) / ( cnt_putocn + 1.0_rp )
873  enddo
874  enddo
875  end if
876  enddo
877  enddo
878  !$acc end kernels
879 
880  if( countup ) then
881  cnt_putocn = cnt_putocn + 1.0_rp
882  endif
883 
884  return
885  end subroutine cpl_putocn
886 
887  !-----------------------------------------------------------------------------
888  subroutine cpl_putlnd( &
889  SFC_TEMP, &
890  SFC_albedo, &
891  SFC_Z0M, &
892  SFC_Z0H, &
893  SFC_Z0E, &
894  SFLX_MW, &
895  SFLX_MU, &
896  SFLX_MV, &
897  SFLX_SH, &
898  SFLX_LH, &
899  SFLX_GH, &
900  SFLX_QTRC, &
901  U10, &
902  V10, &
903  T2, &
904  Q2, &
905  mask, &
906  countup )
907  implicit none
908 
909  real(rp), intent(in) :: sfc_temp (ia,ja)
910  real(rp), intent(in) :: sfc_albedo(ia,ja,n_rad_dir,n_rad_rgn)
911  real(rp), intent(in) :: sfc_z0m (ia,ja)
912  real(rp), intent(in) :: sfc_z0h (ia,ja)
913  real(rp), intent(in) :: sfc_z0e (ia,ja)
914  real(rp), intent(in) :: sflx_mw (ia,ja)
915  real(rp), intent(in) :: sflx_mu (ia,ja)
916  real(rp), intent(in) :: sflx_mv (ia,ja)
917  real(rp), intent(in) :: sflx_sh (ia,ja)
918  real(rp), intent(in) :: sflx_lh (ia,ja)
919  real(rp), intent(in) :: sflx_gh (ia,ja)
920  real(rp), intent(in) :: sflx_qtrc (ia,ja,qa)
921  real(rp), intent(in) :: u10 (ia,ja)
922  real(rp), intent(in) :: v10 (ia,ja)
923  real(rp), intent(in) :: t2 (ia,ja)
924  real(rp), intent(in) :: q2 (ia,ja)
925  logical, intent(in) :: mask (ia,ja)
926  logical, intent(in) :: countup
927 
928  integer :: i, j, iq, idir, irgn
929  !---------------------------------------------------------------------------
930 
931  !$omp parallel do default(none) &
932  !$omp shared(JS,JE,IS,IE,QA) &
933  !$omp shared(LND_SFC_TEMP,LND_SFC_albedo,LND_SFC_Z0M,LND_SFC_Z0H,LND_SFC_Z0E) &
934  !$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) &
935  !$omp shared(LND_U10,LND_V10,LND_T2,LND_Q2,CNT_putLND,SFC_TEMP,SFC_albedo,SFC_Z0M,SFC_Z0H) &
936  !$omp shared(SFC_Z0E,SFLX_MW,SFLX_MU,SFLX_MV,SFLX_SH,SFLX_LH,SFLX_GH,SFLX_QTRC,U10,V10,T2,Q2) &
937  !$omp shared(mask) &
938  !$omp shared(TRACER_CV,TRACER_ENGI0) &
939  !$omp private(i,j,idir,irgn) OMP_SCHEDULE_
940  !$acc kernels
941  do j = js, je
942  do i = is, ie
943  if ( mask(i,j) ) then
944  lnd_sfc_temp(i,j) = lnd_sfc_temp(i,j) * cnt_putlnd + sfc_temp(i,j)
945  lnd_sfc_z0m(i,j) = lnd_sfc_z0m(i,j) * cnt_putlnd + sfc_z0m(i,j)
946  lnd_sfc_z0h(i,j) = lnd_sfc_z0h(i,j) * cnt_putlnd + sfc_z0h(i,j)
947  lnd_sfc_z0e(i,j) = lnd_sfc_z0e(i,j) * cnt_putlnd + sfc_z0e(i,j)
948  lnd_sflx_mw(i,j) = lnd_sflx_mw(i,j) * cnt_putlnd + sflx_mw(i,j)
949  lnd_sflx_mu(i,j) = lnd_sflx_mu(i,j) * cnt_putlnd + sflx_mu(i,j)
950  lnd_sflx_mv(i,j) = lnd_sflx_mv(i,j) * cnt_putlnd + sflx_mv(i,j)
951  lnd_sflx_sh(i,j) = lnd_sflx_sh(i,j) * cnt_putlnd + sflx_sh(i,j)
952  lnd_sflx_lh(i,j) = lnd_sflx_lh(i,j) * cnt_putlnd + sflx_lh(i,j)
953  lnd_sflx_gh(i,j) = lnd_sflx_gh(i,j) * cnt_putlnd + sflx_gh(i,j)
955  do iq = 1, qa
956  lnd_sflx_qtrc(i,j,iq) = lnd_sflx_qtrc(i,j,iq) * cnt_putlnd + sflx_qtrc(i,j,iq)
957  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) )
958  end do
959  lnd_u10(i,j) = lnd_u10(i,j) * cnt_putlnd + u10(i,j)
960  lnd_v10(i,j) = lnd_v10(i,j) * cnt_putlnd + v10(i,j)
961  lnd_t2(i,j) = lnd_t2(i,j) * cnt_putlnd + t2(i,j)
962  lnd_q2(i,j) = lnd_q2(i,j) * cnt_putlnd + q2(i,j)
963  do irgn = i_r_ir, i_r_vis
964  do idir = i_r_direct, i_r_diffuse
965  lnd_sfc_albedo(i,j,idir,irgn) = lnd_sfc_albedo(i,j,idir,irgn) * cnt_putlnd + sfc_albedo(i,j,idir,irgn)
966  enddo
967  enddo
968 
969  lnd_sfc_temp(i,j) = lnd_sfc_temp(i,j) / ( cnt_putlnd + 1.0_rp )
970  lnd_sfc_z0m(i,j) = lnd_sfc_z0m(i,j) / ( cnt_putlnd + 1.0_rp )
971  lnd_sfc_z0h(i,j) = lnd_sfc_z0h(i,j) / ( cnt_putlnd + 1.0_rp )
972  lnd_sfc_z0e(i,j) = lnd_sfc_z0e(i,j) / ( cnt_putlnd + 1.0_rp )
973  lnd_sflx_mw(i,j) = lnd_sflx_mw(i,j) / ( cnt_putlnd + 1.0_rp )
974  lnd_sflx_mu(i,j) = lnd_sflx_mu(i,j) / ( cnt_putlnd + 1.0_rp )
975  lnd_sflx_mv(i,j) = lnd_sflx_mv(i,j) / ( cnt_putlnd + 1.0_rp )
976  lnd_sflx_sh(i,j) = lnd_sflx_sh(i,j) / ( cnt_putlnd + 1.0_rp )
977  lnd_sflx_lh(i,j) = lnd_sflx_lh(i,j) / ( cnt_putlnd + 1.0_rp )
978  lnd_sflx_gh(i,j) = lnd_sflx_gh(i,j) / ( cnt_putlnd + 1.0_rp )
979  lnd_sflx_qtrc(i,j,:) = lnd_sflx_qtrc(i,j,:) / ( cnt_putlnd + 1.0_rp )
980  lnd_sflx_engi(i,j) = lnd_sflx_engi(i,j) / ( cnt_putlnd + 1.0_rp )
981  lnd_u10(i,j) = lnd_u10(i,j) / ( cnt_putlnd + 1.0_rp )
982  lnd_v10(i,j) = lnd_v10(i,j) / ( cnt_putlnd + 1.0_rp )
983  lnd_t2(i,j) = lnd_t2(i,j) / ( cnt_putlnd + 1.0_rp )
984  lnd_q2(i,j) = lnd_q2(i,j) / ( cnt_putlnd + 1.0_rp )
985  do irgn = i_r_ir, i_r_vis
986  do idir = i_r_direct, i_r_diffuse
987  lnd_sfc_albedo(i,j,idir,irgn) = lnd_sfc_albedo(i,j,idir,irgn) / ( cnt_putlnd + 1.0_rp )
988  enddo
989  enddo
990  end if
991  enddo
992  enddo
993  !$acc end kernels
994 
995  if( countup ) then
996  cnt_putlnd = cnt_putlnd + 1.0_rp
997  endif
998 
999  return
1000  end subroutine cpl_putlnd
1001 
1002  !-----------------------------------------------------------------------------
1003  subroutine cpl_puturb( &
1004  SFC_TEMP, &
1005  SFC_albedo, &
1006  SFC_Z0M, &
1007  SFC_Z0H, &
1008  SFC_Z0E, &
1009  SFLX_MW, &
1010  SFLX_MU, &
1011  SFLX_MV, &
1012  SFLX_SH, &
1013  SFLX_LH, &
1014  SFLX_SHEX, &
1015  SFLX_LHEX, &
1016  SFLX_QVEX, &
1017  SFLX_GH, &
1018  SFLX_QTRC, &
1019  U10, &
1020  V10, &
1021  T2, &
1022  Q2, &
1023  mask, &
1024  countup )
1025  implicit none
1026 
1027  real(rp), intent(in) :: sfc_temp (ia,ja)
1028  real(rp), intent(in) :: sfc_albedo(ia,ja,n_rad_dir,n_rad_rgn)
1029  real(rp), intent(in) :: sfc_z0m (ia,ja)
1030  real(rp), intent(in) :: sfc_z0h (ia,ja)
1031  real(rp), intent(in) :: sfc_z0e (ia,ja)
1032  real(rp), intent(in) :: sflx_mw (ia,ja)
1033  real(rp), intent(in) :: sflx_mu (ia,ja)
1034  real(rp), intent(in) :: sflx_mv (ia,ja)
1035  real(rp), intent(in) :: sflx_sh (ia,ja)
1036  real(rp), intent(in) :: sflx_lh (ia,ja)
1037  real(rp), intent(in) :: sflx_shex (ia,ja)
1038  real(rp), intent(in) :: sflx_lhex (ia,ja)
1039  real(rp), intent(in) :: sflx_qvex (ia,ja)
1040  real(rp), intent(in) :: sflx_gh (ia,ja)
1041  real(rp), intent(in) :: sflx_qtrc (ia,ja,qa)
1042  real(rp), intent(in) :: u10 (ia,ja)
1043  real(rp), intent(in) :: v10 (ia,ja)
1044  real(rp), intent(in) :: t2 (ia,ja)
1045  real(rp), intent(in) :: q2 (ia,ja)
1046  logical, intent(in) :: mask (ia,ja)
1047  logical, intent(in) :: countup
1048 
1049  integer :: i, j, iq, idir, irgn
1050  !---------------------------------------------------------------------------
1051 
1052  !$omp parallel do default(none) OMP_SCHEDULE_ &
1053  !$omp shared(JS,JE,IS,IE,QA, &
1054  !$omp URB_SFC_TEMP,URB_SFC_albedo,URB_SFC_Z0M,URB_SFC_Z0H,URB_SFC_Z0E, &
1055  !$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, &
1056  !$omp URB_SFLX_QTRC,URB_SFLX_ENGI,URB_U10,URB_V10,URB_T2,URB_Q2,CNT_putURB, &
1057  !$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, &
1058  !$omp mask, &
1059  !$omp TRACER_CV,TRACER_ENGI0)
1060  !$acc kernels
1061  do j = js, je
1062  do i = is, ie
1063  if ( mask(i,j) ) then
1064  urb_sfc_temp(i,j) = urb_sfc_temp(i,j) * cnt_puturb + sfc_temp(i,j)
1065  urb_sfc_z0m(i,j) = urb_sfc_z0m(i,j) * cnt_puturb + sfc_z0m(i,j)
1066  urb_sfc_z0h(i,j) = urb_sfc_z0h(i,j) * cnt_puturb + sfc_z0h(i,j)
1067  urb_sfc_z0e(i,j) = urb_sfc_z0e(i,j) * cnt_puturb + sfc_z0e(i,j)
1068  urb_sflx_mw(i,j) = urb_sflx_mw(i,j) * cnt_puturb + sflx_mw(i,j)
1069  urb_sflx_mu(i,j) = urb_sflx_mu(i,j) * cnt_puturb + sflx_mu(i,j)
1070  urb_sflx_mv(i,j) = urb_sflx_mv(i,j) * cnt_puturb + sflx_mv(i,j)
1071  urb_sflx_sh(i,j) = urb_sflx_sh(i,j) * cnt_puturb + sflx_sh(i,j)
1072  urb_sflx_lh(i,j) = urb_sflx_lh(i,j) * cnt_puturb + sflx_lh(i,j)
1073  urb_sflx_shex(i,j) = urb_sflx_shex(i,j) * cnt_puturb + sflx_shex(i,j)
1074  urb_sflx_lhex(i,j) = urb_sflx_lhex(i,j) * cnt_puturb + sflx_lhex(i,j)
1075  urb_sflx_qvex(i,j) = urb_sflx_qvex(i,j) * cnt_puturb + sflx_qvex(i,j)
1076  urb_sflx_gh(i,j) = urb_sflx_gh(i,j) * cnt_puturb + sflx_gh(i,j)
1077  urb_sflx_engi(i,j) = urb_sflx_engi(i,j) * cnt_puturb
1078  do iq = 1, qa
1079  urb_sflx_qtrc(i,j,iq) = urb_sflx_qtrc(i,j,iq) * cnt_puturb + sflx_qtrc(i,j,iq)
1080  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) )
1081  end do
1082  urb_u10(i,j) = urb_u10(i,j) * cnt_puturb + u10(i,j)
1083  urb_v10(i,j) = urb_v10(i,j) * cnt_puturb + v10(i,j)
1084  urb_t2(i,j) = urb_t2(i,j) * cnt_puturb + t2(i,j)
1085  urb_q2(i,j) = urb_q2(i,j) * cnt_puturb + q2(i,j)
1086  do irgn = i_r_ir, i_r_vis
1087  do idir = i_r_direct, i_r_diffuse
1088  urb_sfc_albedo(i,j,idir,irgn) = urb_sfc_albedo(i,j,idir,irgn) * cnt_puturb + sfc_albedo(i,j,idir,irgn)
1089  enddo
1090  enddo
1091 
1092  urb_sfc_temp(i,j) = urb_sfc_temp(i,j) / ( cnt_puturb + 1.0_rp )
1093  urb_sfc_z0m(i,j) = urb_sfc_z0m(i,j) / ( cnt_puturb + 1.0_rp )
1094  urb_sfc_z0h(i,j) = urb_sfc_z0h(i,j) / ( cnt_puturb + 1.0_rp )
1095  urb_sfc_z0e(i,j) = urb_sfc_z0e(i,j) / ( cnt_puturb + 1.0_rp )
1096  urb_sflx_mw(i,j) = urb_sflx_mw(i,j) / ( cnt_puturb + 1.0_rp )
1097  urb_sflx_mu(i,j) = urb_sflx_mu(i,j) / ( cnt_puturb + 1.0_rp )
1098  urb_sflx_mv(i,j) = urb_sflx_mv(i,j) / ( cnt_puturb + 1.0_rp )
1099  urb_sflx_sh(i,j) = urb_sflx_sh(i,j) / ( cnt_puturb + 1.0_rp )
1100  urb_sflx_lh(i,j) = urb_sflx_lh(i,j) / ( cnt_puturb + 1.0_rp )
1101  urb_sflx_shex(i,j) = urb_sflx_shex(i,j) / ( cnt_puturb + 1.0_rp )
1102  urb_sflx_lhex(i,j) = urb_sflx_lhex(i,j) / ( cnt_puturb + 1.0_rp )
1103  urb_sflx_qvex(i,j) = urb_sflx_qvex(i,j) / ( cnt_puturb + 1.0_rp )
1104  urb_sflx_gh(i,j) = urb_sflx_gh(i,j) / ( cnt_puturb + 1.0_rp )
1105  urb_sflx_qtrc(i,j,:) = urb_sflx_qtrc(i,j,:) / ( cnt_puturb + 1.0_rp )
1106  urb_sflx_engi(i,j) = urb_sflx_engi(i,j) / ( cnt_puturb + 1.0_rp )
1107  urb_u10(i,j) = urb_u10(i,j) / ( cnt_puturb + 1.0_rp )
1108  urb_v10(i,j) = urb_v10(i,j) / ( cnt_puturb + 1.0_rp )
1109  urb_t2(i,j) = urb_t2(i,j) / ( cnt_puturb + 1.0_rp )
1110  urb_q2(i,j) = urb_q2(i,j) / ( cnt_puturb + 1.0_rp )
1111  do irgn = i_r_ir, i_r_vis
1112  do idir = i_r_direct, i_r_diffuse
1113  urb_sfc_albedo(i,j,idir,irgn) = urb_sfc_albedo(i,j,idir,irgn) / ( cnt_puturb + 1.0_rp )
1114  enddo
1115  enddo
1116  end if
1117  enddo
1118  enddo
1119  !$acc end kernels
1120 
1121  if( countup ) then
1122  cnt_puturb = cnt_puturb + 1.0_rp
1123  endif
1124 
1125  return
1126  end subroutine cpl_puturb
1127 
1128  !-----------------------------------------------------------------------------
1129  subroutine cpl_getsfc_atm( &
1130  SFC_TEMP, &
1131  SFC_albedo, &
1132  SFC_Z0M, &
1133  SFC_Z0H, &
1134  SFC_Z0E, &
1135  SFLX_MW, &
1136  SFLX_MU, &
1137  SFLX_MV, &
1138  SFLX_SH, &
1139  SFLX_LH, &
1140  SFLX_SHEX, &
1141  SFLX_LHEX, &
1142  SFLX_QVEX, &
1143  SFLX_GH, &
1144  SFLX_QTRC, &
1145  SFLX_ENGI, &
1146  U10, &
1147  V10, &
1148  T2, &
1149  Q2 )
1150  use scale_landuse, only: &
1151  fact_ocean => landuse_fact_ocean, &
1152  fact_land => landuse_fact_land, &
1153  fact_urban => landuse_fact_urban
1154  use scale_const, only: &
1155  eps => const_eps
1156  use scale_atmos_hydrometeor, only: &
1157  i_qv
1158  implicit none
1159 
1160  real(rp), intent(out) :: sfc_temp (ia,ja)
1161  real(rp), intent(out) :: sfc_albedo(ia,ja,n_rad_dir,n_rad_rgn)
1162  real(rp), intent(out) :: sfc_z0m (ia,ja)
1163  real(rp), intent(out) :: sfc_z0h (ia,ja)
1164  real(rp), intent(out) :: sfc_z0e (ia,ja)
1165  real(rp), intent(out) :: sflx_mw (ia,ja)
1166  real(rp), intent(out) :: sflx_mu (ia,ja)
1167  real(rp), intent(out) :: sflx_mv (ia,ja)
1168  real(rp), intent(out) :: sflx_sh (ia,ja)
1169  real(rp), intent(out) :: sflx_lh (ia,ja)
1170  real(rp), intent(out) :: sflx_shex (ia,ja)
1171  real(rp), intent(out) :: sflx_lhex (ia,ja)
1172  real(rp), intent(out) :: sflx_qvex (ia,ja)
1173  real(rp), intent(out) :: sflx_gh (ia,ja)
1174  real(rp), intent(out) :: sflx_qtrc (ia,ja,qa)
1175  real(rp), intent(out) :: sflx_engi (ia,ja)
1176  real(rp), intent(out) :: u10 (ia,ja)
1177  real(rp), intent(out) :: v10 (ia,ja)
1178  real(rp), intent(out) :: t2 (ia,ja)
1179  real(rp), intent(out) :: q2 (ia,ja)
1180 
1181  integer :: i, j, idir, irgn, iq
1182  !---------------------------------------------------------------------------
1183 
1184  !$omp parallel do default(none) &
1185  !$omp shared(JS,JE,IS,IE,QA,I_QV,EPS,SFLX_QTRC,SFLX_ENGI,SFC_TEMP,SFC_albedo,SFC_Z0M,SFC_Z0H,SFC_Z0E) &
1186  !$omp shared(SFLX_MW,SFLX_MU,SFLX_MV,SFLX_SH,SFLX_LH,SFLX_SHEX,SFLX_LHEX,SFLX_QVEX,SFLX_GH,U10,V10,T2,Q2) &
1187  !$omp shared(fact_ocean,fact_land,fact_urban,OCN_SFC_TEMP,LND_SFC_TEMP,URB_SFC_TEMP,OCN_SFC_albedo) &
1188  !$omp shared(LND_SFC_albedo,URB_SFC_albedo,OCN_SFC_Z0M,LND_SFC_Z0M,URB_SFC_Z0M) &
1189  !$omp shared(OCN_SFC_Z0H,LND_SFC_Z0H,URB_SFC_Z0H,OCN_SFC_Z0E,LND_SFC_Z0E,URB_SFC_Z0E) &
1190  !$omp shared(OCN_SFLX_MW,LND_SFLX_MW,URB_SFLX_MW,OCN_SFLX_MU,LND_SFLX_MU,URB_SFLX_MU) &
1191  !$omp shared(OCN_SFLX_MV,LND_SFLX_MV,URB_SFLX_MV) &
1192  !$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) &
1193  !$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) &
1194  !$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) &
1195  !$omp private(i,j,iq) OMP_SCHEDULE_
1196  !$acc kernels
1197  do j = js, je
1198  do i = is, ie
1199  sfc_temp(i,j) = fact_ocean(i,j) * ocn_sfc_temp(i,j) &
1200  + fact_land(i,j) * lnd_sfc_temp(i,j) &
1201  + fact_urban(i,j) * urb_sfc_temp(i,j)
1202 
1203  do irgn = i_r_ir, i_r_vis
1204  do idir = i_r_direct, i_r_diffuse
1205  sfc_albedo(i,j,idir,irgn) = fact_ocean(i,j) * ocn_sfc_albedo(i,j,idir,irgn) &
1206  + fact_land(i,j) * lnd_sfc_albedo(i,j,idir,irgn) &
1207  + fact_urban(i,j) * urb_sfc_albedo(i,j,idir,irgn)
1208  enddo
1209  enddo
1210 
1211  sfc_z0m(i,j) = fact_ocean(i,j) * ocn_sfc_z0m(i,j) &
1212  + fact_land(i,j) * lnd_sfc_z0m(i,j) &
1213  + fact_urban(i,j) * urb_sfc_z0m(i,j)
1214 
1215  sfc_z0h(i,j) = fact_ocean(i,j) * ocn_sfc_z0h(i,j) &
1216  + fact_land(i,j) * lnd_sfc_z0h(i,j) &
1217  + fact_urban(i,j) * urb_sfc_z0h(i,j)
1218 
1219  sfc_z0e(i,j) = fact_ocean(i,j) * ocn_sfc_z0e(i,j) &
1220  + fact_land(i,j) * lnd_sfc_z0e(i,j) &
1221  + fact_urban(i,j) * urb_sfc_z0e(i,j)
1222 
1223  sflx_mw(i,j) = fact_ocean(i,j) * ocn_sflx_mw(i,j) &
1224  + fact_land(i,j) * lnd_sflx_mw(i,j) &
1225  + fact_urban(i,j) * urb_sflx_mw(i,j)
1226 
1227  sflx_mu(i,j) = fact_ocean(i,j) * ocn_sflx_mu(i,j) &
1228  + fact_land(i,j) * lnd_sflx_mu(i,j) &
1229  + fact_urban(i,j) * urb_sflx_mu(i,j)
1230 
1231  sflx_mv(i,j) = fact_ocean(i,j) * ocn_sflx_mv(i,j) &
1232  + fact_land(i,j) * lnd_sflx_mv(i,j) &
1233  + fact_urban(i,j) * urb_sflx_mv(i,j)
1234 
1235  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
1236 
1237  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
1238 
1239  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
1240 
1241  sflx_sh(i,j) = fact_ocean(i,j) * ocn_sflx_sh(i,j) &
1242  + fact_land(i,j) * lnd_sflx_sh(i,j) &
1243  + fact_urban(i,j) * urb_sflx_sh(i,j) &
1244  + sflx_shex(i,j) ! grid average
1245 
1246  sflx_lh(i,j) = fact_ocean(i,j) * ocn_sflx_lh(i,j) &
1247  + fact_land(i,j) * lnd_sflx_lh(i,j) &
1248  + fact_urban(i,j) * urb_sflx_lh(i,j) &
1249  + sflx_lhex(i,j) ! grid average
1250 
1251  sflx_gh(i,j) = fact_ocean(i,j) * ocn_sflx_gh(i,j) &
1252  + fact_land(i,j) * lnd_sflx_gh(i,j) &
1253  + fact_urban(i,j) * urb_sflx_gh(i,j)
1254 
1255  do iq = 1, qa
1256  sflx_qtrc(i,j,iq) = fact_ocean(i,j) * ocn_sflx_qtrc(i,j,iq) &
1257  + fact_land(i,j) * lnd_sflx_qtrc(i,j,iq) &
1258  + fact_urban(i,j) * urb_sflx_qtrc(i,j,iq)
1259  enddo
1260  sflx_qtrc(i,j,i_qv) = sflx_qtrc(i,j,i_qv) + sflx_qvex(i,j)
1261 
1262  sflx_engi(i,j) = fact_ocean(i,j) * ocn_sflx_engi(i,j) &
1263  + fact_land(i,j) * lnd_sflx_engi(i,j) &
1264  + fact_urban(i,j) * urb_sflx_engi(i,j)
1265 
1266  u10(i,j) = fact_ocean(i,j) * ocn_u10(i,j) &
1267  + fact_land(i,j) * lnd_u10(i,j) &
1268  + fact_urban(i,j) * urb_u10(i,j)
1269 
1270  v10(i,j) = fact_ocean(i,j) * ocn_v10(i,j) &
1271  + fact_land(i,j) * lnd_v10(i,j) &
1272  + fact_urban(i,j) * urb_v10(i,j)
1273 
1274  t2(i,j) = fact_ocean(i,j) * ocn_t2(i,j) &
1275  + fact_land(i,j) * lnd_t2(i,j) &
1276  + fact_urban(i,j) * urb_t2(i,j)
1277 
1278  q2(i,j) = fact_ocean(i,j) * ocn_q2(i,j) &
1279  + fact_land(i,j) * lnd_q2(i,j) &
1280  + fact_urban(i,j) * urb_q2(i,j)
1281  enddo
1282  enddo
1283  !$acc end kernels
1284 
1285  cnt_putocn = 0.0_rp
1286  cnt_putlnd = 0.0_rp
1287  cnt_puturb = 0.0_rp
1288 
1289  return
1290  end subroutine cpl_getsfc_atm
1291 
1292  !-----------------------------------------------------------------------------
1293  subroutine cpl_getatm_ocn( &
1294  TEMP, &
1295  PRES, &
1296  W, &
1297  U, &
1298  V, &
1299  DENS, &
1300  QV, &
1301  PBL, &
1302  SFC_DENS, &
1303  SFC_PRES, &
1304  SFLX_rad_dn, &
1305  cosSZA, &
1306  SFLX_water, &
1307  SFLX_ENGI )
1308  implicit none
1309 
1310  real(rp), intent(out) :: temp (ia,ja)
1311  real(rp), intent(out) :: pres (ia,ja)
1312  real(rp), intent(out) :: w (ia,ja)
1313  real(rp), intent(out) :: u (ia,ja)
1314  real(rp), intent(out) :: v (ia,ja)
1315  real(rp), intent(out) :: dens (ia,ja)
1316  real(rp), intent(out) :: qv (ia,ja)
1317  real(rp), intent(out) :: pbl (ia,ja)
1318  real(rp), intent(out) :: sfc_dens (ia,ja)
1319  real(rp), intent(out) :: sfc_pres (ia,ja)
1320  real(rp), intent(out) :: sflx_rad_dn(ia,ja,n_rad_dir,n_rad_rgn)
1321  real(rp), intent(out) :: cossza (ia,ja)
1322  real(rp), intent(out) :: sflx_water (ia,ja)
1323  real(rp), intent(out) :: sflx_engi (ia,ja)
1324 
1325  integer :: i, j, idir, irgn
1326  !---------------------------------------------------------------------------
1327 
1328 !OCL XFILL
1329  !$omp parallel do default(none) private(i,j,idir,irgn) OMP_SCHEDULE_ &
1330  !$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) &
1331  !$omp shared(OCN_ATM_TEMP,OCN_ATM_PRES,OCN_ATM_W,OCN_ATM_U,OCN_ATM_V,OCN_ATM_DENS,OCN_ATM_QV) &
1332  !$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)
1333  !$acc kernels
1334  do j = js, je
1335  do i = is, ie
1336  temp(i,j) = ocn_atm_temp(i,j)
1337  pres(i,j) = ocn_atm_pres(i,j)
1338  w(i,j) = ocn_atm_w(i,j)
1339  u(i,j) = ocn_atm_u(i,j)
1340  v(i,j) = ocn_atm_v(i,j)
1341  dens(i,j) = ocn_atm_dens(i,j)
1342  qv(i,j) = ocn_atm_qv(i,j)
1343  pbl(i,j) = ocn_atm_pbl(i,j)
1344  sfc_dens(i,j) = ocn_atm_sfc_dens(i,j)
1345  sfc_pres(i,j) = ocn_atm_sfc_pres(i,j)
1346  cossza(i,j) = ocn_atm_cossza(i,j)
1347  sflx_water(i,j) = ocn_atm_sflx_water(i,j)
1348  sflx_engi(i,j) = ocn_atm_sflx_engi(i,j)
1349  do irgn = i_r_ir, i_r_vis
1350  do idir = i_r_direct, i_r_diffuse
1351  sflx_rad_dn(i,j,idir,irgn) = ocn_atm_sflx_rad_dn(i,j,idir,irgn)
1352  enddo
1353  enddo
1354  enddo
1355  enddo
1356  !$acc end kernels
1357 
1358  cnt_putatm_ocn = 0.0_rp
1359 
1360  return
1361  end subroutine cpl_getatm_ocn
1362 
1363  !-----------------------------------------------------------------------------
1364  subroutine cpl_getatm_lnd( &
1365  TEMP, &
1366  PRES, &
1367  W, &
1368  U, &
1369  V, &
1370  DENS, &
1371  QV, &
1372  PBL, &
1373  SFC_DENS, &
1374  SFC_PRES, &
1375  SFLX_rad_dn, &
1376  cosSZA, &
1377  SFLX_water, &
1378  SFLX_ENGI )
1379  implicit none
1380 
1381  real(rp), intent(out) :: temp (ia,ja)
1382  real(rp), intent(out) :: pres (ia,ja)
1383  real(rp), intent(out) :: w (ia,ja)
1384  real(rp), intent(out) :: u (ia,ja)
1385  real(rp), intent(out) :: v (ia,ja)
1386  real(rp), intent(out) :: dens (ia,ja)
1387  real(rp), intent(out) :: qv (ia,ja)
1388  real(rp), intent(out) :: pbl (ia,ja)
1389  real(rp), intent(out) :: sfc_dens (ia,ja)
1390  real(rp), intent(out) :: sfc_pres (ia,ja)
1391  real(rp), intent(out) :: sflx_rad_dn(ia,ja,n_rad_dir,n_rad_rgn)
1392  real(rp), intent(out) :: cossza (ia,ja)
1393  real(rp), intent(out) :: sflx_water (ia,ja)
1394  real(rp), intent(out) :: sflx_engi (ia,ja)
1395 
1396  integer :: i, j, idir, irgn
1397  !---------------------------------------------------------------------------
1398 
1399 !OCL XFILL
1400  !$omp parallel do default(none) private(i,j,idir,irgn) OMP_SCHEDULE_ &
1401  !$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) &
1402  !$omp shared(LND_ATM_TEMP,LND_ATM_PRES,LND_ATM_W,LND_ATM_U,LND_ATM_V,LND_ATM_DENS,LND_ATM_QV) &
1403  !$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)
1404  !$acc kernels
1405  do j = js, je
1406  do i = is, ie
1407  temp(i,j) = lnd_atm_temp(i,j)
1408  pres(i,j) = lnd_atm_pres(i,j)
1409  w(i,j) = lnd_atm_w(i,j)
1410  u(i,j) = lnd_atm_u(i,j)
1411  v(i,j) = lnd_atm_v(i,j)
1412  dens(i,j) = lnd_atm_dens(i,j)
1413  qv(i,j) = lnd_atm_qv(i,j)
1414  pbl(i,j) = lnd_atm_pbl(i,j)
1415  sfc_dens(i,j) = lnd_atm_sfc_dens(i,j)
1416  sfc_pres(i,j) = lnd_atm_sfc_pres(i,j)
1417  cossza(i,j) = lnd_atm_cossza(i,j)
1418  sflx_water(i,j) = lnd_atm_sflx_water(i,j)
1419  sflx_engi(i,j) = lnd_atm_sflx_engi(i,j)
1420  do irgn = i_r_ir, i_r_vis
1421  do idir = i_r_direct, i_r_diffuse
1422  sflx_rad_dn(i,j,idir,irgn) = lnd_atm_sflx_rad_dn(i,j,idir,irgn)
1423  enddo
1424  enddo
1425  enddo
1426  enddo
1427  !$acc end kernels
1428 
1429  cnt_putatm_lnd = 0.0_rp
1430 
1431  return
1432  end subroutine cpl_getatm_lnd
1433 
1434  !-----------------------------------------------------------------------------
1435  subroutine cpl_getatm_urb( &
1436  TEMP, &
1437  PRES, &
1438  W, &
1439  U, &
1440  V, &
1441  DENS, &
1442  QV, &
1443  PBL, &
1444  SFC_DENS, &
1445  SFC_PRES, &
1446  SFLX_rad_dn, &
1447  cosSZA, &
1448  SFLX_water, &
1449  SFLX_ENGI )
1450  implicit none
1451 
1452  real(rp), intent(out) :: temp (ia,ja)
1453  real(rp), intent(out) :: pres (ia,ja)
1454  real(rp), intent(out) :: w (ia,ja)
1455  real(rp), intent(out) :: u (ia,ja)
1456  real(rp), intent(out) :: v (ia,ja)
1457  real(rp), intent(out) :: dens (ia,ja)
1458  real(rp), intent(out) :: qv (ia,ja)
1459  real(rp), intent(out) :: pbl (ia,ja)
1460  real(rp), intent(out) :: sfc_dens (ia,ja)
1461  real(rp), intent(out) :: sfc_pres (ia,ja)
1462  real(rp), intent(out) :: sflx_rad_dn(ia,ja,n_rad_dir,n_rad_rgn)
1463  real(rp), intent(out) :: cossza (ia,ja)
1464  real(rp), intent(out) :: sflx_water (ia,ja)
1465  real(rp), intent(out) :: sflx_engi (ia,ja)
1466 
1467  integer :: i, j, idir, irgn
1468  !---------------------------------------------------------------------------
1469 
1470 !OCL XFILL
1471  !$omp parallel do default(none) private(i,j,idir,irgn) OMP_SCHEDULE_ &
1472  !$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) &
1473  !$omp shared(URB_ATM_TEMP,URB_ATM_PRES,URB_ATM_W,URB_ATM_U,URB_ATM_V,URB_ATM_DENS,URB_ATM_QV) &
1474  !$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)
1475  !$acc kernels
1476  do j = js, je
1477  do i = is, ie
1478  temp(i,j) = urb_atm_temp(i,j)
1479  pres(i,j) = urb_atm_pres(i,j)
1480  w(i,j) = urb_atm_w(i,j)
1481  u(i,j) = urb_atm_u(i,j)
1482  v(i,j) = urb_atm_v(i,j)
1483  dens(i,j) = urb_atm_dens(i,j)
1484  qv(i,j) = urb_atm_qv(i,j)
1485  pbl(i,j) = urb_atm_pbl(i,j)
1486  sfc_dens(i,j) = urb_atm_sfc_dens(i,j)
1487  sfc_pres(i,j) = urb_atm_sfc_pres(i,j)
1488  cossza(i,j) = urb_atm_cossza(i,j)
1489  sflx_water(i,j) = urb_atm_sflx_water(i,j)
1490  sflx_engi(i,j) = urb_atm_sflx_engi(i,j)
1491  do irgn = i_r_ir, i_r_vis
1492  do idir = i_r_direct, i_r_diffuse
1493  sflx_rad_dn(i,j,idir,irgn) = urb_atm_sflx_rad_dn(i,j,idir,irgn)
1494  enddo
1495  enddo
1496  enddo
1497  enddo
1498  !$acc end kernels
1499 
1500  cnt_putatm_urb = 0.0_rp
1501 
1502  return
1503  end subroutine cpl_getatm_urb
1504 
1505 end module mod_cpl_vars
mod_cpl_vars::ocn_atm_cossza
real(rp), dimension(:,:), allocatable, public ocn_atm_cossza
Definition: mod_cpl_vars.F90:120
mod_cpl_vars::lnd_atm_qv
real(rp), dimension(:,:), allocatable, public lnd_atm_qv
Definition: mod_cpl_vars.F90:131
mod_cpl_vars::ocn_atm_pres
real(rp), dimension(:,:), allocatable, public ocn_atm_pres
Definition: mod_cpl_vars.F90:110
scale_landuse::landuse_fact_lake
real(rp), dimension(:,:), allocatable, public landuse_fact_lake
lake factor
Definition: scale_landuse.F90:48
scale_cpl_sfc_index::n_rad_dir
integer, parameter, public n_rad_dir
Definition: scale_cpl_sfc_index.F90:36
mod_cpl_vars::lnd_sflx_engi
real(rp), dimension(:,:), allocatable, public lnd_sflx_engi
Definition: mod_cpl_vars.F90:80
mod_cpl_vars::urb_sflx_shex
real(rp), dimension(:,:), allocatable, public urb_sflx_shex
Definition: mod_cpl_vars.F90:97
mod_cpl_vars::urb_atm_w
real(rp), dimension(:,:), allocatable, public urb_atm_w
Definition: mod_cpl_vars.F90:143
mod_cpl_vars::ocn_sflx_mu
real(rp), dimension(:,:), allocatable, public ocn_sflx_mu
Definition: mod_cpl_vars.F90:55
scale_landuse::landuse_fact_ocean
real(rp), dimension(:,:), allocatable, public landuse_fact_ocean
ocean factor
Definition: scale_landuse.F90:45
mod_urban_admin::urban_do
logical, public urban_do
Definition: mod_urban_admin.F90:32
mod_cpl_vars::lnd_atm_sfc_dens
real(rp), dimension(:,:), allocatable, public lnd_atm_sfc_dens
Definition: mod_cpl_vars.F90:133
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_cpl_sfc_index::i_r_direct
integer, parameter, public i_r_direct
Definition: scale_cpl_sfc_index.F90:37
mod_cpl_vars::ocn_sflx_lh
real(rp), dimension(:,:), allocatable, public ocn_sflx_lh
Definition: mod_cpl_vars.F90:58
scale_tracer::qa
integer, public qa
Definition: scale_tracer.F90:35
mod_cpl_vars::ocn_sflx_mv
real(rp), dimension(:,:), allocatable, public ocn_sflx_mv
Definition: mod_cpl_vars.F90:56
mod_cpl_vars::lnd_sfc_temp
real(rp), dimension(:,:), allocatable, public lnd_sfc_temp
Definition: mod_cpl_vars.F90:68
mod_cpl_vars::cnt_putocn
real(rp), public cnt_putocn
Definition: mod_cpl_vars.F90:160
mod_cpl_vars::urb_sflx_mv
real(rp), dimension(:,:), allocatable, public urb_sflx_mv
Definition: mod_cpl_vars.F90:94
mod_cpl_vars::cnt_putatm_urb
real(rp), public cnt_putatm_urb
Definition: mod_cpl_vars.F90:159
mod_cpl_vars::urb_sflx_lh
real(rp), dimension(:,:), allocatable, public urb_sflx_lh
Definition: mod_cpl_vars.F90:96
scale_cpl_sfc_index::i_r_diffuse
integer, parameter, public i_r_diffuse
Definition: scale_cpl_sfc_index.F90:38
mod_cpl_vars
module COUPLER Variables
Definition: mod_cpl_vars.F90:12
mod_cpl_vars::urb_atm_v
real(rp), dimension(:,:), allocatable, public urb_atm_v
Definition: mod_cpl_vars.F90:145
scale_precision
module PRECISION
Definition: scale_precision.F90:14
mod_cpl_vars::urb_atm_sfc_dens
real(rp), dimension(:,:), allocatable, public urb_atm_sfc_dens
Definition: mod_cpl_vars.F90:149
scale_cpl_sfc_index::i_r_ir
integer, parameter, public i_r_ir
Definition: scale_cpl_sfc_index.F90:29
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
mod_cpl_vars::ocn_atm_pbl
real(rp), dimension(:,:), allocatable, public ocn_atm_pbl
Definition: mod_cpl_vars.F90:116
mod_cpl_vars::ocn_atm_u
real(rp), dimension(:,:), allocatable, public ocn_atm_u
Definition: mod_cpl_vars.F90:112
mod_cpl_vars::lnd_sflx_mv
real(rp), dimension(:,:), allocatable, public lnd_sflx_mv
Definition: mod_cpl_vars.F90:75
mod_cpl_vars::lnd_atm_pres
real(rp), dimension(:,:), allocatable, public lnd_atm_pres
Definition: mod_cpl_vars.F90:126
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
mod_cpl_vars::lnd_atm_sflx_rad_dn
real(rp), dimension(:,:,:,:), allocatable, public lnd_atm_sflx_rad_dn
Definition: mod_cpl_vars.F90:135
mod_cpl_vars::lnd_q2
real(rp), dimension(:,:), allocatable, public lnd_q2
Definition: mod_cpl_vars.F90:84
mod_cpl_vars::ocn_atm_w
real(rp), dimension(:,:), allocatable, public ocn_atm_w
Definition: mod_cpl_vars.F90:111
mod_cpl_vars::urb_sfc_z0h
real(rp), dimension(:,:), allocatable, public urb_sfc_z0h
Definition: mod_cpl_vars.F90:90
mod_cpl_vars::lnd_u10
real(rp), dimension(:,:), allocatable, public lnd_u10
Definition: mod_cpl_vars.F90:81
mod_cpl_vars::urb_atm_sflx_water
real(rp), dimension(:,:), allocatable, public urb_atm_sflx_water
Definition: mod_cpl_vars.F90:153
mod_cpl_vars::cnt_putlnd
real(rp), public cnt_putlnd
Definition: mod_cpl_vars.F90:161
scale_atmos_hydrometeor::atmos_hydrometeor_dry
logical, public atmos_hydrometeor_dry
Definition: scale_atmos_hydrometeor.F90:114
mod_cpl_vars::ocn_q2
real(rp), dimension(:,:), allocatable, public ocn_q2
Definition: mod_cpl_vars.F90:65
scale_tracer::tracer_engi0
real(rp), dimension(qa_max), public tracer_engi0
Definition: scale_tracer.F90:45
mod_cpl_vars::cpl_getsfc_atm
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)
Definition: mod_cpl_vars.F90:1150
mod_cpl_vars::lnd_sfc_z0e
real(rp), dimension(:,:), allocatable, public lnd_sfc_z0e
Definition: mod_cpl_vars.F90:72
mod_cpl_vars::urb_sfc_z0m
real(rp), dimension(:,:), allocatable, public urb_sfc_z0m
Definition: mod_cpl_vars.F90:89
mod_cpl_vars::ocn_atm_qv
real(rp), dimension(:,:), allocatable, public ocn_atm_qv
Definition: mod_cpl_vars.F90:115
scale_prc
module PROCESS
Definition: scale_prc.F90:11
mod_cpl_vars::lnd_t2
real(rp), dimension(:,:), allocatable, public lnd_t2
Definition: mod_cpl_vars.F90:83
mod_cpl_vars::lnd_atm_w
real(rp), dimension(:,:), allocatable, public lnd_atm_w
Definition: mod_cpl_vars.F90:127
mod_cpl_vars::ocn_atm_sfc_pres
real(rp), dimension(:,:), allocatable, public ocn_atm_sfc_pres
Definition: mod_cpl_vars.F90:118
mod_cpl_vars::urb_sflx_qvex
real(rp), dimension(:,:), allocatable, public urb_sflx_qvex
Definition: mod_cpl_vars.F90:99
mod_cpl_vars::cpl_putatm
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)
Definition: mod_cpl_vars.F90:596
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
mod_cpl_vars::lnd_sflx_gh
real(rp), dimension(:,:), allocatable, public lnd_sflx_gh
Definition: mod_cpl_vars.F90:78
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
scale_io
module STDIO
Definition: scale_io.F90:10
mod_cpl_vars::ocn_sflx_sh
real(rp), dimension(:,:), allocatable, public ocn_sflx_sh
Definition: mod_cpl_vars.F90:57
mod_cpl_vars::lnd_sfc_z0h
real(rp), dimension(:,:), allocatable, public lnd_sfc_z0h
Definition: mod_cpl_vars.F90:71
mod_ocean_admin
module Ocean admin
Definition: mod_ocean_admin.F90:11
mod_cpl_vars::urb_sfc_albedo
real(rp), dimension(:,:,:,:), allocatable, public urb_sfc_albedo
Definition: mod_cpl_vars.F90:88
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
mod_cpl_vars::ocn_sfc_z0m
real(rp), dimension(:,:), allocatable, public ocn_sfc_z0m
Definition: mod_cpl_vars.F90:51
mod_land_admin
module Land admin
Definition: mod_land_admin.F90:11
mod_lake_admin::lake_do
logical, public lake_do
Definition: mod_lake_admin.F90:32
mod_cpl_vars::cpl_putocn
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)
Definition: mod_cpl_vars.F90:793
mod_cpl_vars::lnd_sflx_sh
real(rp), dimension(:,:), allocatable, public lnd_sflx_sh
Definition: mod_cpl_vars.F90:76
mod_cpl_vars::lnd_sflx_lh
real(rp), dimension(:,:), allocatable, public lnd_sflx_lh
Definition: mod_cpl_vars.F90:77
mod_cpl_vars::urb_atm_u
real(rp), dimension(:,:), allocatable, public urb_atm_u
Definition: mod_cpl_vars.F90:144
mod_cpl_vars::lnd_sflx_mw
real(rp), dimension(:,:), allocatable, public lnd_sflx_mw
Definition: mod_cpl_vars.F90:73
scale_tracer::tracer_cv
real(rp), dimension(qa_max), public tracer_cv
Definition: scale_tracer.F90:42
mod_cpl_vars::cpl_getatm_urb
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)
Definition: mod_cpl_vars.F90:1450
mod_cpl_vars::lnd_sflx_qtrc
real(rp), dimension(:,:,:), allocatable, public lnd_sflx_qtrc
Definition: mod_cpl_vars.F90:79
mod_cpl_vars::urb_atm_pbl
real(rp), dimension(:,:), allocatable, public urb_atm_pbl
Definition: mod_cpl_vars.F90:148
mod_cpl_vars::cnt_putatm_ocn
real(rp), public cnt_putatm_ocn
Definition: mod_cpl_vars.F90:157
mod_cpl_vars::lnd_atm_cossza
real(rp), dimension(:,:), allocatable, public lnd_atm_cossza
Definition: mod_cpl_vars.F90:136
scale_prof
module profiler
Definition: scale_prof.F90:11
mod_cpl_vars::ocn_atm_sflx_water
real(rp), dimension(:,:), allocatable, public ocn_atm_sflx_water
Definition: mod_cpl_vars.F90:121
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
mod_cpl_vars::ocn_atm_v
real(rp), dimension(:,:), allocatable, public ocn_atm_v
Definition: mod_cpl_vars.F90:113
mod_cpl_vars::ocn_sflx_mw
real(rp), dimension(:,:), allocatable, public ocn_sflx_mw
Definition: mod_cpl_vars.F90:54
mod_cpl_vars::cpl_vars_setup
subroutine, public cpl_vars_setup
Setup.
Definition: mod_cpl_vars.F90:177
mod_cpl_vars::lnd_sflx_mu
real(rp), dimension(:,:), allocatable, public lnd_sflx_mu
Definition: mod_cpl_vars.F90:74
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
mod_cpl_vars::lnd_atm_sflx_water
real(rp), dimension(:,:), allocatable, public lnd_atm_sflx_water
Definition: mod_cpl_vars.F90:137
mod_cpl_vars::lnd_atm_v
real(rp), dimension(:,:), allocatable, public lnd_atm_v
Definition: mod_cpl_vars.F90:129
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_atmos_hydrometeor::i_qv
integer, public i_qv
Definition: scale_atmos_hydrometeor.F90:93
mod_cpl_vars::ocn_u10
real(rp), dimension(:,:), allocatable, public ocn_u10
Definition: mod_cpl_vars.F90:62
mod_cpl_vars::urb_u10
real(rp), dimension(:,:), allocatable, public urb_u10
Definition: mod_cpl_vars.F90:103
scale_landuse::landuse_fact_land
real(rp), dimension(:,:), allocatable, public landuse_fact_land
land factor
Definition: scale_landuse.F90:46
mod_cpl_vars::ocn_sflx_engi
real(rp), dimension(:,:), allocatable, public ocn_sflx_engi
Definition: mod_cpl_vars.F90:61
mod_cpl_vars::urb_sflx_lhex
real(rp), dimension(:,:), allocatable, public urb_sflx_lhex
Definition: mod_cpl_vars.F90:98
mod_cpl_vars::urb_atm_cossza
real(rp), dimension(:,:), allocatable, public urb_atm_cossza
Definition: mod_cpl_vars.F90:152
scale_debug
module DEBUG
Definition: scale_debug.F90:11
mod_cpl_vars::lnd_sfc_z0m
real(rp), dimension(:,:), allocatable, public lnd_sfc_z0m
Definition: mod_cpl_vars.F90:70
mod_cpl_vars::urb_sfc_z0e
real(rp), dimension(:,:), allocatable, public urb_sfc_z0e
Definition: mod_cpl_vars.F90:91
mod_ocean_admin::ocean_do
logical, public ocean_do
Definition: mod_ocean_admin.F90:32
mod_cpl_vars::urb_atm_pres
real(rp), dimension(:,:), allocatable, public urb_atm_pres
Definition: mod_cpl_vars.F90:142
mod_cpl_vars::urb_sfc_temp
real(rp), dimension(:,:), allocatable, public urb_sfc_temp
Definition: mod_cpl_vars.F90:87
mod_cpl_vars::lnd_atm_sfc_pres
real(rp), dimension(:,:), allocatable, public lnd_atm_sfc_pres
Definition: mod_cpl_vars.F90:134
mod_cpl_vars::urb_sflx_sh
real(rp), dimension(:,:), allocatable, public urb_sflx_sh
Definition: mod_cpl_vars.F90:95
mod_cpl_vars::ocn_atm_dens
real(rp), dimension(:,:), allocatable, public ocn_atm_dens
Definition: mod_cpl_vars.F90:114
mod_cpl_vars::lnd_atm_dens
real(rp), dimension(:,:), allocatable, public lnd_atm_dens
Definition: mod_cpl_vars.F90:130
scale_cpl_sfc_index
module coupler / surface-atmospehre
Definition: scale_cpl_sfc_index.F90:11
mod_cpl_vars::ocn_sfc_z0e
real(rp), dimension(:,:), allocatable, public ocn_sfc_z0e
Definition: mod_cpl_vars.F90:53
mod_cpl_vars::cpl_puturb
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)
Definition: mod_cpl_vars.F90:1025
mod_land_admin::land_do
logical, public land_do
Definition: mod_land_admin.F90:41
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
mod_cpl_vars::ocn_t2
real(rp), dimension(:,:), allocatable, public ocn_t2
Definition: mod_cpl_vars.F90:64
scale_cpl_sfc_index::i_r_vis
integer, parameter, public i_r_vis
Definition: scale_cpl_sfc_index.F90:31
mod_cpl_vars::lnd_atm_pbl
real(rp), dimension(:,:), allocatable, public lnd_atm_pbl
Definition: mod_cpl_vars.F90:132
mod_cpl_vars::urb_sflx_gh
real(rp), dimension(:,:), allocatable, public urb_sflx_gh
Definition: mod_cpl_vars.F90:100
mod_cpl_vars::cpl_getatm_ocn
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)
Definition: mod_cpl_vars.F90:1308
mod_cpl_vars::lnd_atm_sflx_engi
real(rp), dimension(:,:), allocatable, public lnd_atm_sflx_engi
Definition: mod_cpl_vars.F90:138
scale_atmos_hydrometeor::lhf
real(rp), public lhf
latent heat of fusion for use [J/kg]
Definition: scale_atmos_hydrometeor.F90:146
mod_cpl_vars::cpl_vars_finalize
subroutine, public cpl_vars_finalize
Finalize.
Definition: mod_cpl_vars.F90:462
mod_cpl_vars::urb_t2
real(rp), dimension(:,:), allocatable, public urb_t2
Definition: mod_cpl_vars.F90:105
mod_cpl_vars::urb_atm_dens
real(rp), dimension(:,:), allocatable, public urb_atm_dens
Definition: mod_cpl_vars.F90:146
mod_cpl_vars::urb_sflx_mu
real(rp), dimension(:,:), allocatable, public urb_sflx_mu
Definition: mod_cpl_vars.F90:93
mod_cpl_vars::lnd_atm_temp
real(rp), dimension(:,:), allocatable, public lnd_atm_temp
Definition: mod_cpl_vars.F90:125
mod_cpl_vars::lnd_v10
real(rp), dimension(:,:), allocatable, public lnd_v10
Definition: mod_cpl_vars.F90:82
mod_cpl_vars::lnd_atm_u
real(rp), dimension(:,:), allocatable, public lnd_atm_u
Definition: mod_cpl_vars.F90:128
scale_landuse::landuse_fact_urban
real(rp), dimension(:,:), allocatable, public landuse_fact_urban
urban factor
Definition: scale_landuse.F90:47
mod_cpl_vars::ocn_sflx_gh
real(rp), dimension(:,:), allocatable, public ocn_sflx_gh
Definition: mod_cpl_vars.F90:59
mod_cpl_vars::ocn_atm_temp
real(rp), dimension(:,:), allocatable, public ocn_atm_temp
Definition: mod_cpl_vars.F90:109
mod_cpl_vars::ocn_sflx_qtrc
real(rp), dimension(:,:,:), allocatable, public ocn_sflx_qtrc
Definition: mod_cpl_vars.F90:60
mod_cpl_vars::urb_sflx_qtrc
real(rp), dimension(:,:,:), allocatable, public urb_sflx_qtrc
Definition: mod_cpl_vars.F90:101
mod_cpl_vars::ocn_sfc_temp
real(rp), dimension(:,:), allocatable, public ocn_sfc_temp
Definition: mod_cpl_vars.F90:49
scale_landuse
module LANDUSE
Definition: scale_landuse.F90:19
mod_cpl_vars::urb_atm_temp
real(rp), dimension(:,:), allocatable, public urb_atm_temp
Definition: mod_cpl_vars.F90:141
mod_urban_admin
module Urban admin
Definition: mod_urban_admin.F90:11
mod_cpl_vars::urb_q2
real(rp), dimension(:,:), allocatable, public urb_q2
Definition: mod_cpl_vars.F90:106
scale_cpl_sfc_index::n_rad_rgn
integer, parameter, public n_rad_rgn
Definition: scale_cpl_sfc_index.F90:28
mod_cpl_vars::urb_atm_qv
real(rp), dimension(:,:), allocatable, public urb_atm_qv
Definition: mod_cpl_vars.F90:147
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
mod_cpl_vars::ocn_atm_sfc_dens
real(rp), dimension(:,:), allocatable, public ocn_atm_sfc_dens
Definition: mod_cpl_vars.F90:117
mod_cpl_vars::urb_atm_sfc_pres
real(rp), dimension(:,:), allocatable, public urb_atm_sfc_pres
Definition: mod_cpl_vars.F90:150
mod_lake_admin
module Lake admin
Definition: mod_lake_admin.F90:11
mod_cpl_vars::cpl_putlnd
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)
Definition: mod_cpl_vars.F90:907
mod_cpl_vars::ocn_atm_sflx_rad_dn
real(rp), dimension(:,:,:,:), allocatable, public ocn_atm_sflx_rad_dn
Definition: mod_cpl_vars.F90:119
mod_cpl_vars::ocn_sfc_albedo
real(rp), dimension(:,:,:,:), allocatable, public ocn_sfc_albedo
Definition: mod_cpl_vars.F90:50
mod_cpl_vars::ocn_atm_sflx_engi
real(rp), dimension(:,:), allocatable, public ocn_atm_sflx_engi
Definition: mod_cpl_vars.F90:122
mod_cpl_vars::cnt_putatm_lnd
real(rp), public cnt_putatm_lnd
Definition: mod_cpl_vars.F90:158
mod_cpl_vars::urb_sflx_mw
real(rp), dimension(:,:), allocatable, public urb_sflx_mw
Definition: mod_cpl_vars.F90:92
mod_cpl_vars::ocn_sfc_z0h
real(rp), dimension(:,:), allocatable, public ocn_sfc_z0h
Definition: mod_cpl_vars.F90:52
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56
mod_cpl_vars::cpl_getatm_lnd
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)
Definition: mod_cpl_vars.F90:1379
mod_cpl_vars::urb_sflx_engi
real(rp), dimension(:,:), allocatable, public urb_sflx_engi
Definition: mod_cpl_vars.F90:102
scale_atmos_hydrometeor::cv_water
real(rp), public cv_water
CV for water [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:151
mod_cpl_vars::cnt_puturb
real(rp), public cnt_puturb
Definition: mod_cpl_vars.F90:162
mod_cpl_vars::lnd_sfc_albedo
real(rp), dimension(:,:,:,:), allocatable, public lnd_sfc_albedo
Definition: mod_cpl_vars.F90:69
scale_atmos_hydrometeor::cv_ice
real(rp), public cv_ice
CV for ice [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:153
mod_cpl_vars::urb_v10
real(rp), dimension(:,:), allocatable, public urb_v10
Definition: mod_cpl_vars.F90:104
mod_cpl_vars::ocn_v10
real(rp), dimension(:,:), allocatable, public ocn_v10
Definition: mod_cpl_vars.F90:63
mod_cpl_vars::urb_atm_sflx_engi
real(rp), dimension(:,:), allocatable, public urb_atm_sflx_engi
Definition: mod_cpl_vars.F90:154
mod_cpl_vars::urb_atm_sflx_rad_dn
real(rp), dimension(:,:,:,:), allocatable, public urb_atm_sflx_rad_dn
Definition: mod_cpl_vars.F90:151