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