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