SCALE-RM
mod_cpl_vars.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_stdio
18  use scale_prof
19  use scale_debug
21  use scale_tracer
22 
23  use scale_const, only: &
24  i_sw => const_i_sw, &
25  i_lw => const_i_lw
26  !-----------------------------------------------------------------------------
27  implicit none
28  private
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public procedure
32  !
33  public :: cpl_vars_setup
34 
35  public :: cpl_putatm
36  public :: cpl_putocn
37  public :: cpl_putlnd
38  public :: cpl_puturb
39  public :: cpl_getsfc_atm
40  public :: cpl_getatm_ocn
41  public :: cpl_getatm_lnd
42  public :: cpl_getatm_urb
43 
44  !-----------------------------------------------------------------------------
45  !
46  !++ Public parameters & variables
47  !
48 
49  ! Input from ocean model
50  real(RP), public, allocatable :: ocn_sfc_temp (:,:) ! ocean surface skin temperature [K]
51  real(RP), public, allocatable :: ocn_sfc_albedo(:,:,:) ! ocean surface albedo [0-1]
52  real(RP), public, allocatable :: ocn_sfc_z0m (:,:) ! ocean surface roughness length for momemtum [m]
53  real(RP), public, allocatable :: ocn_sfc_z0h (:,:) ! ocean surface roughness length for heat [m]
54  real(RP), public, allocatable :: ocn_sfc_z0e (:,:) ! ocean surface roughness length for vapor [m]
55  real(RP), public, allocatable :: ocn_sflx_mw (:,:) ! ocean surface w-momentum flux [kg/m2/s]
56  real(RP), public, allocatable :: ocn_sflx_mu (:,:) ! ocean surface u-momentum flux [kg/m2/s]
57  real(RP), public, allocatable :: ocn_sflx_mv (:,:) ! ocean surface v-momentum flux [kg/m2/s]
58  real(RP), public, allocatable :: ocn_sflx_sh (:,:) ! ocean surface sensible heat flux [J/m2/s]
59  real(RP), public, allocatable :: ocn_sflx_lh (:,:) ! ocean surface latent heat flux [J/m2/s]
60  real(RP), public, allocatable :: ocn_sflx_wh (:,:) ! ocean surface water heat flux [J/m2/s]
61  real(RP), public, allocatable :: ocn_sflx_evap (:,:) ! ocean surface water vapor flux [kg/m2/s]
62  real(RP), public, allocatable :: ocn_u10 (:,:) ! ocean velocity u at 10m [m/s]
63  real(RP), public, allocatable :: ocn_v10 (:,:) ! ocean velocity v at 10m [m/s]
64  real(RP), public, allocatable :: ocn_t2 (:,:) ! ocean temperature at 2m [K]
65  real(RP), public, allocatable :: ocn_q2 (:,:) ! ocean water vapor at 2m [kg/kg]
66 
67  ! Input from land model
68  real(RP), public, allocatable :: lnd_sfc_temp (:,:) ! land surface skin temperature [K]
69  real(RP), public, allocatable :: lnd_sfc_albedo(:,:,:) ! land surface albedo [0-1]
70  real(RP), public, allocatable :: lnd_sfc_z0m (:,:) ! land surface roughness length for momemtum [m]
71  real(RP), public, allocatable :: lnd_sfc_z0h (:,:) ! land surface roughness length for heat [m]
72  real(RP), public, allocatable :: lnd_sfc_z0e (:,:) ! land surface roughness length for vapor [m]
73  real(RP), public, allocatable :: lnd_sflx_mw (:,:) ! land surface w-momentum flux [kg/m2/s]
74  real(RP), public, allocatable :: lnd_sflx_mu (:,:) ! land surface u-momentum flux [kg/m2/s]
75  real(RP), public, allocatable :: lnd_sflx_mv (:,:) ! land surface v-momentum flux [kg/m2/s]
76  real(RP), public, allocatable :: lnd_sflx_sh (:,:) ! land surface sensible heat flux [J/m2/s]
77  real(RP), public, allocatable :: lnd_sflx_lh (:,:) ! land surface latent heat flux [J/m2/s]
78  real(RP), public, allocatable :: lnd_sflx_gh (:,:) ! land surface ground heat flux [J/m2/s]
79  real(RP), public, allocatable :: lnd_sflx_evap (:,:) ! land surface water vapor flux [kg/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 [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/m2/s]
92  real(RP), public, allocatable :: urb_sflx_mu (:,:) ! urban surface u-momentum flux [kg/m2/s]
93  real(RP), public, allocatable :: urb_sflx_mv (:,:) ! urban surface v-momentum flux [kg/m2/s]
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_gh (:,:) ! urban surface ground heat flux [J/m2/s]
97  real(RP), public, allocatable :: urb_sflx_evap (:,:) ! urban surface water vapor flux [kg/m2/s]
98  real(RP), public, allocatable :: urb_u10 (:,:) ! urban velocity u at 10m [m/s]
99  real(RP), public, allocatable :: urb_v10 (:,:) ! urban velocity v at 10m [m/s]
100  real(RP), public, allocatable :: urb_t2 (:,:) ! urban temperature at 2m [K]
101  real(RP), public, allocatable :: urb_q2 (:,:) ! urban water vapor at 2m [kg/kg]
102 
103  ! Output to ocean model
104  real(RP), public, allocatable :: ocn_atm_temp (:,:) ! temperature at the lowermost atmosphere layer [K]
105  real(RP), public, allocatable :: ocn_atm_pres (:,:) ! pressure at the lowermost atmosphere layer [Pa]
106  real(RP), public, allocatable :: ocn_atm_w (:,:) ! velocity w at the lowermost atmosphere layer [m/s]
107  real(RP), public, allocatable :: ocn_atm_u (:,:) ! velocity u at the lowermost atmosphere layer [m/s]
108  real(RP), public, allocatable :: ocn_atm_v (:,:) ! velocity v at the lowermost atmosphere layer [m/s]
109  real(RP), public, allocatable :: ocn_atm_dens (:,:) ! density at the lowermost atmosphere layer [kg/m3]
110  real(RP), public, allocatable :: ocn_atm_qv (:,:) ! water vapor at the lowermost atmosphere layer [kg/kg]
111  real(RP), public, allocatable :: ocn_atm_pbl (:,:) ! the top of atmospheric mixing layer [m]
112  real(RP), public, allocatable :: ocn_atm_sfc_pres (:,:) ! surface pressure [Pa]
113  real(RP), public, allocatable :: ocn_atm_sflx_rad_dn(:,:,:,:) ! downward radiation flux (SW/LW,direct/diffuse) [J/m2/s]
114  real(RP), public, allocatable :: ocn_atm_cossza (:,:) ! cos(solar zenith angle) [0-1]
115  real(RP), public, allocatable :: ocn_atm_sflx_rain (:,:) ! liquid water flux [kg/m2/s]
116  real(RP), public, allocatable :: ocn_atm_sflx_snow (:,:) ! ice water flux [kg/m2/s]
117 
118  ! Output to land model
119  real(RP), public, allocatable :: lnd_atm_temp (:,:) ! temperature at the lowermost atmosphere layer [K]
120  real(RP), public, allocatable :: lnd_atm_pres (:,:) ! pressure at the lowermost atmosphere layer [Pa]
121  real(RP), public, allocatable :: lnd_atm_w (:,:) ! velocity w at the lowermost atmosphere layer [m/s]
122  real(RP), public, allocatable :: lnd_atm_u (:,:) ! velocity u at the lowermost atmosphere layer [m/s]
123  real(RP), public, allocatable :: lnd_atm_v (:,:) ! velocity v at the lowermost atmosphere layer [m/s]
124  real(RP), public, allocatable :: lnd_atm_dens (:,:) ! density at the lowermost atmosphere layer [kg/m3]
125  real(RP), public, allocatable :: lnd_atm_qv (:,:) ! water vapor at the lowermost atmosphere layer [kg/kg]
126  real(RP), public, allocatable :: lnd_atm_pbl (:,:) ! the top of atmospheric mixing layer [m]
127  real(RP), public, allocatable :: lnd_atm_sfc_pres (:,:) ! surface pressure [Pa]
128  real(RP), public, allocatable :: lnd_atm_sflx_rad_dn(:,:,:,:) ! downward radiation flux (SW/LW,direct/diffuse) [J/m2/s]
129  real(RP), public, allocatable :: lnd_atm_cossza (:,:) ! cos(solar zenith angle) [0-1]
130  real(RP), public, allocatable :: lnd_atm_sflx_rain (:,:) ! liquid water flux [kg/m2/s]
131  real(RP), public, allocatable :: lnd_atm_sflx_snow (:,:) ! ice water flux [kg/m2/s]
132 
133  ! Output to urban model
134  real(RP), public, allocatable :: urb_atm_temp (:,:) ! temperature at the lowermost atmosphere layer [K]
135  real(RP), public, allocatable :: urb_atm_pres (:,:) ! pressure at the lowermost atmosphere layer [Pa]
136  real(RP), public, allocatable :: urb_atm_w (:,:) ! velocity w at the lowermost atmosphere layer [m/s]
137  real(RP), public, allocatable :: urb_atm_u (:,:) ! velocity u at the lowermost atmosphere layer [m/s]
138  real(RP), public, allocatable :: urb_atm_v (:,:) ! velocity v at the lowermost atmosphere layer [m/s]
139  real(RP), public, allocatable :: urb_atm_dens (:,:) ! density at the lowermost atmosphere layer [kg/m3]
140  real(RP), public, allocatable :: urb_atm_qv (:,:) ! water vapor at the lowermost atmosphere layer [kg/kg]
141  real(RP), public, allocatable :: urb_atm_pbl (:,:) ! the top of atmospheric mixing layer [m]
142  real(RP), public, allocatable :: urb_atm_sfc_pres (:,:) ! surface pressure [Pa]
143  real(RP), public, allocatable :: urb_atm_sflx_rad_dn(:,:,:,:) ! downward radiation flux (SW/LW,direct/diffuse) [J/m2/s]
144  real(RP), public, allocatable :: urb_atm_cossza (:,:) ! cos(solar zenith angle) [0-1]
145  real(RP), public, allocatable :: urb_atm_sflx_rain (:,:) ! liquid water flux [kg/m2/s]
146  real(RP), public, allocatable :: urb_atm_sflx_snow (:,:) ! ice water flux [kg/m2/s]
147 
148  ! counter
149  real(RP), public :: cnt_putatm_ocn ! put counter for atmos to ocean
150  real(RP), public :: cnt_putatm_lnd ! put counter for atmos to land
151  real(RP), public :: cnt_putatm_urb ! put counter for atmos to urban
152  real(RP), public :: cnt_putocn ! put counter for ocean
153  real(RP), public :: cnt_putlnd ! put counter for land
154  real(RP), public :: cnt_puturb ! put counter for urban
155 
156  !-----------------------------------------------------------------------------
157  !
158  !++ Private procedure
159  !
160  !-----------------------------------------------------------------------------
161  !
162  !++ Private parameters & variables
163  !
164  !-----------------------------------------------------------------------------
165 contains
166  !-----------------------------------------------------------------------------
168  subroutine cpl_vars_setup
169  use scale_const, only: &
170  undef => const_undef
171  implicit none
172  !---------------------------------------------------------------------------
173 
174  if( io_l ) write(io_fid_log,*)
175  if( io_l ) write(io_fid_log,*) '++++++ Module[VARS] / Categ[CPL] / Origin[SCALE-RM]'
176 
177  allocate( ocn_sfc_temp(ia,ja) )
178  allocate( ocn_sfc_albedo(ia,ja,2) )
179  allocate( ocn_sfc_z0m(ia,ja) )
180  allocate( ocn_sfc_z0h(ia,ja) )
181  allocate( ocn_sfc_z0e(ia,ja) )
182  allocate( ocn_sflx_mu(ia,ja) )
183  allocate( ocn_sflx_mv(ia,ja) )
184  allocate( ocn_sflx_mw(ia,ja) )
185  allocate( ocn_sflx_sh(ia,ja) )
186  allocate( ocn_sflx_lh(ia,ja) )
187  allocate( ocn_sflx_wh(ia,ja) )
188  allocate( ocn_sflx_evap(ia,ja) )
189  allocate( ocn_u10(ia,ja) )
190  allocate( ocn_v10(ia,ja) )
191  allocate( ocn_t2(ia,ja) )
192  allocate( ocn_q2(ia,ja) )
193  ocn_sfc_temp(:,:) = undef
194  ocn_sfc_albedo(:,:,:) = undef
195  ocn_sfc_z0m(:,:) = undef
196  ocn_sfc_z0h(:,:) = undef
197  ocn_sfc_z0e(:,:) = undef
198  ocn_sflx_mu(:,:) = undef
199  ocn_sflx_mv(:,:) = undef
200  ocn_sflx_mw(:,:) = undef
201  ocn_sflx_sh(:,:) = undef
202  ocn_sflx_lh(:,:) = undef
203  ocn_sflx_wh(:,:) = undef
204  ocn_sflx_evap(:,:) = undef
205  ocn_u10(:,:) = undef
206  ocn_v10(:,:) = undef
207  ocn_t2(:,:) = undef
208  ocn_q2(:,:) = undef
209 
210  allocate( lnd_sfc_temp(ia,ja) )
211  allocate( lnd_sfc_albedo(ia,ja,2) )
212  allocate( lnd_sfc_z0m(ia,ja) )
213  allocate( lnd_sfc_z0h(ia,ja) )
214  allocate( lnd_sfc_z0e(ia,ja) )
215  allocate( lnd_sflx_mu(ia,ja) )
216  allocate( lnd_sflx_mv(ia,ja) )
217  allocate( lnd_sflx_mw(ia,ja) )
218  allocate( lnd_sflx_sh(ia,ja) )
219  allocate( lnd_sflx_lh(ia,ja) )
220  allocate( lnd_sflx_gh(ia,ja) )
221  allocate( lnd_sflx_evap(ia,ja) )
222  allocate( lnd_u10(ia,ja) )
223  allocate( lnd_v10(ia,ja) )
224  allocate( lnd_t2(ia,ja) )
225  allocate( lnd_q2(ia,ja) )
226  lnd_sfc_temp(:,:) = undef
227  lnd_sfc_albedo(:,:,:) = undef
228  lnd_sfc_z0m(:,:) = undef
229  lnd_sfc_z0h(:,:) = undef
230  lnd_sfc_z0e(:,:) = undef
231  lnd_sflx_mu(:,:) = undef
232  lnd_sflx_mv(:,:) = undef
233  lnd_sflx_mw(:,:) = undef
234  lnd_sflx_sh(:,:) = undef
235  lnd_sflx_lh(:,:) = undef
236  lnd_sflx_gh(:,:) = undef
237  lnd_sflx_evap(:,:) = undef
238  lnd_u10(:,:) = undef
239  lnd_v10(:,:) = undef
240  lnd_t2(:,:) = undef
241  lnd_q2(:,:) = undef
242 
243  allocate( urb_sfc_temp(ia,ja) )
244  allocate( urb_sfc_albedo(ia,ja,2) )
245  allocate( urb_sfc_z0m(ia,ja) )
246  allocate( urb_sfc_z0h(ia,ja) )
247  allocate( urb_sfc_z0e(ia,ja) )
248  allocate( urb_sflx_mu(ia,ja) )
249  allocate( urb_sflx_mv(ia,ja) )
250  allocate( urb_sflx_mw(ia,ja) )
251  allocate( urb_sflx_sh(ia,ja) )
252  allocate( urb_sflx_lh(ia,ja) )
253  allocate( urb_sflx_gh(ia,ja) )
254  allocate( urb_sflx_evap(ia,ja) )
255  allocate( urb_u10(ia,ja) )
256  allocate( urb_v10(ia,ja) )
257  allocate( urb_t2(ia,ja) )
258  allocate( urb_q2(ia,ja) )
259  urb_sfc_temp(:,:) = undef
260  urb_sfc_albedo(:,:,:) = undef
261  urb_sfc_z0m(:,:) = undef
262  urb_sfc_z0h(:,:) = undef
263  urb_sfc_z0e(:,:) = undef
264  urb_sflx_mu(:,:) = undef
265  urb_sflx_mv(:,:) = undef
266  urb_sflx_mw(:,:) = undef
267  urb_sflx_sh(:,:) = undef
268  urb_sflx_lh(:,:) = undef
269  urb_sflx_gh(:,:) = undef
270  urb_sflx_evap(:,:) = undef
271  urb_u10(:,:) = undef
272  urb_v10(:,:) = undef
273  urb_t2(:,:) = undef
274  urb_q2(:,:) = undef
275 
276  allocate( ocn_atm_temp(ia,ja) )
277  allocate( ocn_atm_pres(ia,ja) )
278  allocate( ocn_atm_w(ia,ja) )
279  allocate( ocn_atm_u(ia,ja) )
280  allocate( ocn_atm_v(ia,ja) )
281  allocate( ocn_atm_dens(ia,ja) )
282  allocate( ocn_atm_qv(ia,ja) )
283  allocate( ocn_atm_pbl(ia,ja) )
284  allocate( ocn_atm_sfc_pres(ia,ja) )
285  allocate( ocn_atm_sflx_rad_dn(ia,ja,2,2) )
286  allocate( ocn_atm_cossza(ia,ja) )
287  allocate( ocn_atm_sflx_rain(ia,ja) )
288  allocate( ocn_atm_sflx_snow(ia,ja) )
289  ocn_atm_temp(:,:) = undef
290  ocn_atm_pres(:,:) = undef
291  ocn_atm_w(:,:) = undef
292  ocn_atm_u(:,:) = undef
293  ocn_atm_v(:,:) = undef
294  ocn_atm_dens(:,:) = undef
295  ocn_atm_qv(:,:) = undef
296  ocn_atm_pbl(:,:) = undef
297  ocn_atm_sfc_pres(:,:) = undef
298  ocn_atm_sflx_rad_dn(:,:,2,2) = undef
299  ocn_atm_cossza(:,:) = undef
300  ocn_atm_sflx_rain(:,:) = undef
301  ocn_atm_sflx_snow(:,:) = undef
302 
303  allocate( lnd_atm_temp(ia,ja) )
304  allocate( lnd_atm_pres(ia,ja) )
305  allocate( lnd_atm_w(ia,ja) )
306  allocate( lnd_atm_u(ia,ja) )
307  allocate( lnd_atm_v(ia,ja) )
308  allocate( lnd_atm_dens(ia,ja) )
309  allocate( lnd_atm_qv(ia,ja) )
310  allocate( lnd_atm_pbl(ia,ja) )
311  allocate( lnd_atm_sfc_pres(ia,ja) )
312  allocate( lnd_atm_sflx_rad_dn(ia,ja,2,2) )
313  allocate( lnd_atm_cossza(ia,ja) )
314  allocate( lnd_atm_sflx_rain(ia,ja) )
315  allocate( lnd_atm_sflx_snow(ia,ja) )
316  lnd_atm_temp(:,:) = undef
317  lnd_atm_pres(:,:) = undef
318  lnd_atm_w(:,:) = undef
319  lnd_atm_u(:,:) = undef
320  lnd_atm_v(:,:) = undef
321  lnd_atm_dens(:,:) = undef
322  lnd_atm_qv(:,:) = undef
323  lnd_atm_pbl(:,:) = undef
324  lnd_atm_sfc_pres(:,:) = undef
325  lnd_atm_sflx_rad_dn(:,:,2,2) = undef
326  lnd_atm_cossza(:,:) = undef
327  lnd_atm_sflx_rain(:,:) = undef
328  lnd_atm_sflx_snow(:,:) = undef
329 
330  allocate( urb_atm_temp(ia,ja) )
331  allocate( urb_atm_pres(ia,ja) )
332  allocate( urb_atm_w(ia,ja) )
333  allocate( urb_atm_u(ia,ja) )
334  allocate( urb_atm_v(ia,ja) )
335  allocate( urb_atm_dens(ia,ja) )
336  allocate( urb_atm_qv(ia,ja) )
337  allocate( urb_atm_pbl(ia,ja) )
338  allocate( urb_atm_sfc_pres(ia,ja) )
339  allocate( urb_atm_sflx_rad_dn(ia,ja,2,2) )
340  allocate( urb_atm_cossza(ia,ja) )
341  allocate( urb_atm_sflx_rain(ia,ja) )
342  allocate( urb_atm_sflx_snow(ia,ja) )
343  urb_atm_temp(:,:) = undef
344  urb_atm_pres(:,:) = undef
345  urb_atm_w(:,:) = undef
346  urb_atm_u(:,:) = undef
347  urb_atm_v(:,:) = undef
348  urb_atm_dens(:,:) = undef
349  urb_atm_qv(:,:) = undef
350  urb_atm_pbl(:,:) = undef
351  urb_atm_sfc_pres(:,:) = undef
352  urb_atm_sflx_rad_dn(:,:,2,2) = undef
353  urb_atm_cossza(:,:) = undef
354  urb_atm_sflx_rain(:,:) = undef
355  urb_atm_sflx_snow(:,:) = undef
356 
357  ! counter intialize
358  cnt_putatm_ocn = 0.0_rp
359  cnt_putatm_lnd = 0.0_rp
360  cnt_putatm_urb = 0.0_rp
361  cnt_putocn = 0.0_rp
362  cnt_putlnd = 0.0_rp
363  cnt_puturb = 0.0_rp
364 
365  return
366  end subroutine cpl_vars_setup
367 
368  !-----------------------------------------------------------------------------
369  subroutine cpl_putatm( &
370  TEMP, &
371  PRES, &
372  W, &
373  U, &
374  V, &
375  DENS, &
376  QTRC, &
377  PBL, &
378  SFC_PRES, &
379  SFLX_rad_dn, &
380  cosSZA, &
381  SFLX_rain, &
382  SFLX_snow, &
383  countup )
384  implicit none
385 
386  ! arguments
387  real(RP), intent(in) :: TEMP (ia,ja)
388  real(RP), intent(in) :: PRES (ia,ja)
389  real(RP), intent(in) :: W (ia,ja)
390  real(RP), intent(in) :: U (ia,ja)
391  real(RP), intent(in) :: V (ia,ja)
392  real(RP), intent(in) :: DENS (ia,ja)
393  real(RP), intent(in) :: QTRC (ia,ja,qa)
394  real(RP), intent(in) :: PBL (ia,ja)
395  real(RP), intent(in) :: SFC_PRES (ia,ja)
396  real(RP), intent(in) :: SFLX_rad_dn(ia,ja,2,2)
397  real(RP), intent(in) :: cosSZA (ia,ja)
398  real(RP), intent(in) :: SFLX_rain (ia,ja)
399  real(RP), intent(in) :: SFLX_snow (ia,ja)
400 
401  logical, intent(in) :: countup
402 
403  ! works
404  integer :: i, j
405  !---------------------------------------------------------------------------
406 
407  do j = js, je
408  do i = is, ie
409  ! for ocean
410  ocn_atm_temp(i,j) = ocn_atm_temp(i,j) * cnt_putatm_ocn + temp(i,j)
411  ocn_atm_pres(i,j) = ocn_atm_pres(i,j) * cnt_putatm_ocn + pres(i,j)
412  ocn_atm_w(i,j) = ocn_atm_w(i,j) * cnt_putatm_ocn + w(i,j)
413  ocn_atm_u(i,j) = ocn_atm_u(i,j) * cnt_putatm_ocn + u(i,j)
414  ocn_atm_v(i,j) = ocn_atm_v(i,j) * cnt_putatm_ocn + v(i,j)
415  ocn_atm_dens(i,j) = ocn_atm_dens(i,j) * cnt_putatm_ocn + dens(i,j)
416  ocn_atm_qv(i,j) = ocn_atm_qv(i,j) * cnt_putatm_ocn + qtrc(i,j,i_qv)
417  ocn_atm_pbl(i,j) = ocn_atm_pbl(i,j) * cnt_putatm_ocn + pbl(i,j)
418  ocn_atm_sfc_pres(i,j) = ocn_atm_sfc_pres(i,j) * cnt_putatm_ocn + sfc_pres(i,j)
419  ocn_atm_sflx_rad_dn(i,j,:,:) = ocn_atm_sflx_rad_dn(i,j,:,:) * cnt_putatm_ocn + sflx_rad_dn(i,j,:,:)
420  ocn_atm_cossza(i,j) = ocn_atm_cossza(i,j) * cnt_putatm_ocn + cossza(i,j)
421  ocn_atm_sflx_rain(i,j) = ocn_atm_sflx_rain(i,j) * cnt_putatm_ocn + sflx_rain(i,j)
422  ocn_atm_sflx_snow(i,j) = ocn_atm_sflx_snow(i,j) * cnt_putatm_ocn + sflx_snow(i,j)
423  ! for land
424  lnd_atm_temp(i,j) = lnd_atm_temp(i,j) * cnt_putatm_lnd + temp(i,j)
425  lnd_atm_pres(i,j) = lnd_atm_pres(i,j) * cnt_putatm_lnd + pres(i,j)
426  lnd_atm_w(i,j) = lnd_atm_w(i,j) * cnt_putatm_lnd + w(i,j)
427  lnd_atm_u(i,j) = lnd_atm_u(i,j) * cnt_putatm_lnd + u(i,j)
428  lnd_atm_v(i,j) = lnd_atm_v(i,j) * cnt_putatm_lnd + v(i,j)
429  lnd_atm_dens(i,j) = lnd_atm_dens(i,j) * cnt_putatm_lnd + dens(i,j)
430  lnd_atm_qv(i,j) = lnd_atm_qv(i,j) * cnt_putatm_lnd + qtrc(i,j,i_qv)
431  lnd_atm_pbl(i,j) = lnd_atm_pbl(i,j) * cnt_putatm_lnd + pbl(i,j)
432  lnd_atm_sfc_pres(i,j) = lnd_atm_sfc_pres(i,j) * cnt_putatm_lnd + sfc_pres(i,j)
433  lnd_atm_sflx_rad_dn(i,j,:,:) = lnd_atm_sflx_rad_dn(i,j,:,:) * cnt_putatm_lnd + sflx_rad_dn(i,j,:,:)
434  lnd_atm_cossza(i,j) = lnd_atm_cossza(i,j) * cnt_putatm_lnd + cossza(i,j)
435  lnd_atm_sflx_rain(i,j) = lnd_atm_sflx_rain(i,j) * cnt_putatm_lnd + sflx_rain(i,j)
436  lnd_atm_sflx_snow(i,j) = lnd_atm_sflx_snow(i,j) * cnt_putatm_lnd + sflx_snow(i,j)
437  ! for urban
438  urb_atm_temp(i,j) = urb_atm_temp(i,j) * cnt_putatm_urb + temp(i,j)
439  urb_atm_pres(i,j) = urb_atm_pres(i,j) * cnt_putatm_urb + pres(i,j)
440  urb_atm_w(i,j) = urb_atm_w(i,j) * cnt_putatm_urb + w(i,j)
441  urb_atm_u(i,j) = urb_atm_u(i,j) * cnt_putatm_urb + u(i,j)
442  urb_atm_v(i,j) = urb_atm_v(i,j) * cnt_putatm_urb + v(i,j)
443  urb_atm_dens(i,j) = urb_atm_dens(i,j) * cnt_putatm_urb + dens(i,j)
444  urb_atm_qv(i,j) = urb_atm_qv(i,j) * cnt_putatm_urb + qtrc(i,j,i_qv)
445  urb_atm_pbl(i,j) = urb_atm_pbl(i,j) * cnt_putatm_urb + pbl(i,j)
446  urb_atm_sfc_pres(i,j) = urb_atm_sfc_pres(i,j) * cnt_putatm_urb + sfc_pres(i,j)
447  urb_atm_sflx_rad_dn(i,j,:,:) = urb_atm_sflx_rad_dn(i,j,:,:) * cnt_putatm_urb + sflx_rad_dn(i,j,:,:)
448  urb_atm_cossza(i,j) = urb_atm_cossza(i,j) * cnt_putatm_urb + cossza(i,j)
449  urb_atm_sflx_rain(i,j) = urb_atm_sflx_rain(i,j) * cnt_putatm_urb + sflx_rain(i,j)
450  urb_atm_sflx_snow(i,j) = urb_atm_sflx_snow(i,j) * cnt_putatm_urb + sflx_snow(i,j)
451  enddo
452  enddo
453 
454  do j = js, je
455  do i = is, ie
456  ! for ocean
457  ocn_atm_temp(i,j) = ocn_atm_temp(i,j) / ( cnt_putatm_ocn + 1.0_rp )
458  ocn_atm_pres(i,j) = ocn_atm_pres(i,j) / ( cnt_putatm_ocn + 1.0_rp )
459  ocn_atm_w(i,j) = ocn_atm_w(i,j) / ( cnt_putatm_ocn + 1.0_rp )
460  ocn_atm_u(i,j) = ocn_atm_u(i,j) / ( cnt_putatm_ocn + 1.0_rp )
461  ocn_atm_v(i,j) = ocn_atm_v(i,j) / ( cnt_putatm_ocn + 1.0_rp )
462  ocn_atm_dens(i,j) = ocn_atm_dens(i,j) / ( cnt_putatm_ocn + 1.0_rp )
463  ocn_atm_qv(i,j) = ocn_atm_qv(i,j) / ( cnt_putatm_ocn + 1.0_rp )
464  ocn_atm_pbl(i,j) = ocn_atm_pbl(i,j) / ( cnt_putatm_ocn + 1.0_rp )
465  ocn_atm_sfc_pres(i,j) = ocn_atm_sfc_pres(i,j) / ( cnt_putatm_ocn + 1.0_rp )
466  ocn_atm_sflx_rad_dn(i,j,:,:) = ocn_atm_sflx_rad_dn(i,j,:,:) / ( cnt_putatm_ocn + 1.0_rp )
467  ocn_atm_cossza(i,j) = ocn_atm_cossza(i,j) / ( cnt_putatm_ocn + 1.0_rp )
468  ocn_atm_sflx_rain(i,j) = ocn_atm_sflx_rain(i,j) / ( cnt_putatm_ocn + 1.0_rp )
469  ocn_atm_sflx_snow(i,j) = ocn_atm_sflx_snow(i,j) / ( cnt_putatm_ocn + 1.0_rp )
470  ! for ocean
471  lnd_atm_temp(i,j) = lnd_atm_temp(i,j) / ( cnt_putatm_lnd + 1.0_rp )
472  lnd_atm_pres(i,j) = lnd_atm_pres(i,j) / ( cnt_putatm_lnd + 1.0_rp )
473  lnd_atm_w(i,j) = lnd_atm_w(i,j) / ( cnt_putatm_lnd + 1.0_rp )
474  lnd_atm_u(i,j) = lnd_atm_u(i,j) / ( cnt_putatm_lnd + 1.0_rp )
475  lnd_atm_v(i,j) = lnd_atm_v(i,j) / ( cnt_putatm_lnd + 1.0_rp )
476  lnd_atm_dens(i,j) = lnd_atm_dens(i,j) / ( cnt_putatm_lnd + 1.0_rp )
477  lnd_atm_qv(i,j) = lnd_atm_qv(i,j) / ( cnt_putatm_lnd + 1.0_rp )
478  lnd_atm_pbl(i,j) = lnd_atm_pbl(i,j) / ( cnt_putatm_lnd + 1.0_rp )
479  lnd_atm_sfc_pres(i,j) = lnd_atm_sfc_pres(i,j) / ( cnt_putatm_lnd + 1.0_rp )
480  lnd_atm_sflx_rad_dn(i,j,:,:) = lnd_atm_sflx_rad_dn(i,j,:,:) / ( cnt_putatm_lnd + 1.0_rp )
481  lnd_atm_cossza(i,j) = lnd_atm_cossza(i,j) / ( cnt_putatm_lnd + 1.0_rp )
482  lnd_atm_sflx_rain(i,j) = lnd_atm_sflx_rain(i,j) / ( cnt_putatm_lnd + 1.0_rp )
483  lnd_atm_sflx_snow(i,j) = lnd_atm_sflx_snow(i,j) / ( cnt_putatm_lnd + 1.0_rp )
484  ! for ocean
485  urb_atm_temp(i,j) = urb_atm_temp(i,j) / ( cnt_putatm_urb + 1.0_rp )
486  urb_atm_pres(i,j) = urb_atm_pres(i,j) / ( cnt_putatm_urb + 1.0_rp )
487  urb_atm_w(i,j) = urb_atm_w(i,j) / ( cnt_putatm_urb + 1.0_rp )
488  urb_atm_u(i,j) = urb_atm_u(i,j) / ( cnt_putatm_urb + 1.0_rp )
489  urb_atm_v(i,j) = urb_atm_v(i,j) / ( cnt_putatm_urb + 1.0_rp )
490  urb_atm_dens(i,j) = urb_atm_dens(i,j) / ( cnt_putatm_urb + 1.0_rp )
491  urb_atm_qv(i,j) = urb_atm_qv(i,j) / ( cnt_putatm_urb + 1.0_rp )
492  urb_atm_pbl(i,j) = urb_atm_pbl(i,j) / ( cnt_putatm_urb + 1.0_rp )
493  urb_atm_sfc_pres(i,j) = urb_atm_sfc_pres(i,j) / ( cnt_putatm_urb + 1.0_rp )
494  urb_atm_sflx_rad_dn(i,j,:,:) = urb_atm_sflx_rad_dn(i,j,:,:) / ( cnt_putatm_urb + 1.0_rp )
495  urb_atm_cossza(i,j) = urb_atm_cossza(i,j) / ( cnt_putatm_urb + 1.0_rp )
496  urb_atm_sflx_rain(i,j) = urb_atm_sflx_rain(i,j) / ( cnt_putatm_urb + 1.0_rp )
497  urb_atm_sflx_snow(i,j) = urb_atm_sflx_snow(i,j) / ( cnt_putatm_urb + 1.0_rp )
498  enddo
499  enddo
500 
501  if( countup ) then
502  cnt_putatm_ocn = cnt_putatm_ocn + 1.0_rp
503  cnt_putatm_lnd = cnt_putatm_lnd + 1.0_rp
504  cnt_putatm_urb = cnt_putatm_urb + 1.0_rp
505  end if
506 
507  return
508  end subroutine cpl_putatm
509 
510  !-----------------------------------------------------------------------------
511  subroutine cpl_putocn( &
512  SFC_TEMP, &
513  SFC_albedo, &
514  SFC_Z0M, &
515  SFC_Z0H, &
516  SFC_Z0E, &
517  SFLX_MW, &
518  SFLX_MU, &
519  SFLX_MV, &
520  SFLX_SH, &
521  SFLX_LH, &
522  SFLX_WH, &
523  SFLX_evap, &
524  U10, &
525  V10, &
526  T2, &
527  Q2, &
528  countup )
529  implicit none
530 
531  ! arguments
532  real(RP), intent(in) :: SFC_TEMP (ia,ja)
533  real(RP), intent(in) :: SFC_albedo(ia,ja,2)
534  real(RP), intent(in) :: SFC_Z0M (ia,ja)
535  real(RP), intent(in) :: SFC_Z0H (ia,ja)
536  real(RP), intent(in) :: SFC_Z0E (ia,ja)
537  real(RP), intent(in) :: SFLX_MW (ia,ja)
538  real(RP), intent(in) :: SFLX_MU (ia,ja)
539  real(RP), intent(in) :: SFLX_MV (ia,ja)
540  real(RP), intent(in) :: SFLX_SH (ia,ja)
541  real(RP), intent(in) :: SFLX_LH (ia,ja)
542  real(RP), intent(in) :: SFLX_WH (ia,ja)
543  real(RP), intent(in) :: SFLX_evap (ia,ja)
544  real(RP), intent(in) :: U10 (ia,ja)
545  real(RP), intent(in) :: V10 (ia,ja)
546  real(RP), intent(in) :: T2 (ia,ja)
547  real(RP), intent(in) :: Q2 (ia,ja)
548 
549  logical, intent(in) :: countup
550 
551  ! works
552  integer :: i, j
553  !---------------------------------------------------------------------------
554 
555  do j = js, je
556  do i = is, ie
557  ocn_sfc_temp(i,j) = ocn_sfc_temp(i,j) * cnt_putocn + sfc_temp(i,j)
558  ocn_sfc_albedo(i,j,i_lw) = ocn_sfc_albedo(i,j,i_lw) * cnt_putocn + sfc_albedo(i,j,i_lw)
559  ocn_sfc_albedo(i,j,i_sw) = ocn_sfc_albedo(i,j,i_sw) * cnt_putocn + sfc_albedo(i,j,i_sw)
560  ocn_sfc_z0m(i,j) = ocn_sfc_z0m(i,j) * cnt_putocn + sfc_z0m(i,j)
561  ocn_sfc_z0h(i,j) = ocn_sfc_z0h(i,j) * cnt_putocn + sfc_z0h(i,j)
562  ocn_sfc_z0e(i,j) = ocn_sfc_z0e(i,j) * cnt_putocn + sfc_z0e(i,j)
563  ocn_sflx_mw(i,j) = ocn_sflx_mw(i,j) * cnt_putocn + sflx_mw(i,j)
564  ocn_sflx_mu(i,j) = ocn_sflx_mu(i,j) * cnt_putocn + sflx_mu(i,j)
565  ocn_sflx_mv(i,j) = ocn_sflx_mv(i,j) * cnt_putocn + sflx_mv(i,j)
566  ocn_sflx_sh(i,j) = ocn_sflx_sh(i,j) * cnt_putocn + sflx_sh(i,j)
567  ocn_sflx_lh(i,j) = ocn_sflx_lh(i,j) * cnt_putocn + sflx_lh(i,j)
568  ocn_sflx_wh(i,j) = ocn_sflx_wh(i,j) * cnt_putocn + sflx_wh(i,j)
569  ocn_sflx_evap(i,j) = ocn_sflx_evap(i,j) * cnt_putocn + sflx_evap(i,j)
570  ocn_u10(i,j) = ocn_u10(i,j) * cnt_putocn + u10(i,j)
571  ocn_v10(i,j) = ocn_v10(i,j) * cnt_putocn + v10(i,j)
572  ocn_t2(i,j) = ocn_t2(i,j) * cnt_putocn + t2(i,j)
573  ocn_q2(i,j) = ocn_q2(i,j) * cnt_putocn + q2(i,j)
574  enddo
575  enddo
576 
577  do j = js, je
578  do i = is, ie
579  ocn_sfc_temp(i,j) = ocn_sfc_temp(i,j) / ( cnt_putocn + 1.0_rp )
580  ocn_sfc_albedo(i,j,i_lw) = ocn_sfc_albedo(i,j,i_lw) / ( cnt_putocn + 1.0_rp )
581  ocn_sfc_albedo(i,j,i_sw) = ocn_sfc_albedo(i,j,i_sw) / ( cnt_putocn + 1.0_rp )
582  ocn_sfc_z0m(i,j) = ocn_sfc_z0m(i,j) / ( cnt_putocn + 1.0_rp )
583  ocn_sfc_z0h(i,j) = ocn_sfc_z0h(i,j) / ( cnt_putocn + 1.0_rp )
584  ocn_sfc_z0e(i,j) = ocn_sfc_z0e(i,j) / ( cnt_putocn + 1.0_rp )
585  ocn_sflx_mw(i,j) = ocn_sflx_mw(i,j) / ( cnt_putocn + 1.0_rp )
586  ocn_sflx_mu(i,j) = ocn_sflx_mu(i,j) / ( cnt_putocn + 1.0_rp )
587  ocn_sflx_mv(i,j) = ocn_sflx_mv(i,j) / ( cnt_putocn + 1.0_rp )
588  ocn_sflx_sh(i,j) = ocn_sflx_sh(i,j) / ( cnt_putocn + 1.0_rp )
589  ocn_sflx_lh(i,j) = ocn_sflx_lh(i,j) / ( cnt_putocn + 1.0_rp )
590  ocn_sflx_wh(i,j) = ocn_sflx_wh(i,j) / ( cnt_putocn + 1.0_rp )
591  ocn_sflx_evap(i,j) = ocn_sflx_evap(i,j) / ( cnt_putocn + 1.0_rp )
592  ocn_u10(i,j) = ocn_u10(i,j) / ( cnt_putocn + 1.0_rp )
593  ocn_v10(i,j) = ocn_v10(i,j) / ( cnt_putocn + 1.0_rp )
594  ocn_t2(i,j) = ocn_t2(i,j) / ( cnt_putocn + 1.0_rp )
595  ocn_q2(i,j) = ocn_q2(i,j) / ( cnt_putocn + 1.0_rp )
596  enddo
597  enddo
598 
599  if( countup ) then
600  cnt_putocn = cnt_putocn + 1.0_rp
601  end if
602 
603  return
604  end subroutine cpl_putocn
605 
606  !-----------------------------------------------------------------------------
607  subroutine cpl_putlnd( &
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_GH, &
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_GH (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  lnd_sfc_temp(i,j) = lnd_sfc_temp(i,j) * cnt_putlnd + sfc_temp(i,j)
654  lnd_sfc_albedo(i,j,i_lw) = lnd_sfc_albedo(i,j,i_lw) * cnt_putlnd + sfc_albedo(i,j,i_lw)
655  lnd_sfc_albedo(i,j,i_sw) = lnd_sfc_albedo(i,j,i_sw) * cnt_putlnd + sfc_albedo(i,j,i_sw)
656  lnd_sfc_z0m(i,j) = lnd_sfc_z0m(i,j) * cnt_putlnd + sfc_z0m(i,j)
657  lnd_sfc_z0h(i,j) = lnd_sfc_z0h(i,j) * cnt_putlnd + sfc_z0h(i,j)
658  lnd_sfc_z0e(i,j) = lnd_sfc_z0e(i,j) * cnt_putlnd + sfc_z0e(i,j)
659  lnd_sflx_mw(i,j) = lnd_sflx_mw(i,j) * cnt_putlnd + sflx_mw(i,j)
660  lnd_sflx_mu(i,j) = lnd_sflx_mu(i,j) * cnt_putlnd + sflx_mu(i,j)
661  lnd_sflx_mv(i,j) = lnd_sflx_mv(i,j) * cnt_putlnd + sflx_mv(i,j)
662  lnd_sflx_sh(i,j) = lnd_sflx_sh(i,j) * cnt_putlnd + sflx_sh(i,j)
663  lnd_sflx_lh(i,j) = lnd_sflx_lh(i,j) * cnt_putlnd + sflx_lh(i,j)
664  lnd_sflx_gh(i,j) = lnd_sflx_gh(i,j) * cnt_putlnd + sflx_gh(i,j)
665  lnd_sflx_evap(i,j) = lnd_sflx_evap(i,j) * cnt_putlnd + sflx_evap(i,j)
666  lnd_u10(i,j) = lnd_u10(i,j) * cnt_putlnd + u10(i,j)
667  lnd_v10(i,j) = lnd_v10(i,j) * cnt_putlnd + v10(i,j)
668  lnd_t2(i,j) = lnd_t2(i,j) * cnt_putlnd + t2(i,j)
669  lnd_q2(i,j) = lnd_q2(i,j) * cnt_putlnd + q2(i,j)
670  enddo
671  enddo
672 
673  do j = js, je
674  do i = is, ie
675  lnd_sfc_temp(i,j) = lnd_sfc_temp(i,j) / ( cnt_putlnd + 1.0_rp )
676  lnd_sfc_albedo(i,j,i_lw) = lnd_sfc_albedo(i,j,i_lw) / ( cnt_putlnd + 1.0_rp )
677  lnd_sfc_albedo(i,j,i_sw) = lnd_sfc_albedo(i,j,i_sw) / ( cnt_putlnd + 1.0_rp )
678  lnd_sfc_z0m(i,j) = lnd_sfc_z0m(i,j) / ( cnt_putlnd + 1.0_rp )
679  lnd_sfc_z0h(i,j) = lnd_sfc_z0h(i,j) / ( cnt_putlnd + 1.0_rp )
680  lnd_sfc_z0e(i,j) = lnd_sfc_z0e(i,j) / ( cnt_putlnd + 1.0_rp )
681  lnd_sflx_mw(i,j) = lnd_sflx_mw(i,j) / ( cnt_putlnd + 1.0_rp )
682  lnd_sflx_mu(i,j) = lnd_sflx_mu(i,j) / ( cnt_putlnd + 1.0_rp )
683  lnd_sflx_mv(i,j) = lnd_sflx_mv(i,j) / ( cnt_putlnd + 1.0_rp )
684  lnd_sflx_sh(i,j) = lnd_sflx_sh(i,j) / ( cnt_putlnd + 1.0_rp )
685  lnd_sflx_lh(i,j) = lnd_sflx_lh(i,j) / ( cnt_putlnd + 1.0_rp )
686  lnd_sflx_gh(i,j) = lnd_sflx_gh(i,j) / ( cnt_putlnd + 1.0_rp )
687  lnd_sflx_evap(i,j) = lnd_sflx_evap(i,j) / ( cnt_putlnd + 1.0_rp )
688  lnd_u10(i,j) = lnd_u10(i,j) / ( cnt_putlnd + 1.0_rp )
689  lnd_v10(i,j) = lnd_v10(i,j) / ( cnt_putlnd + 1.0_rp )
690  lnd_t2(i,j) = lnd_t2(i,j) / ( cnt_putlnd + 1.0_rp )
691  lnd_q2(i,j) = lnd_q2(i,j) / ( cnt_putlnd + 1.0_rp )
692  enddo
693  enddo
694 
695  if( countup ) then
696  cnt_putlnd = cnt_putlnd + 1.0_rp
697  end if
698 
699  return
700  end subroutine cpl_putlnd
701 
702  !-----------------------------------------------------------------------------
703  subroutine cpl_puturb( &
704  SFC_TEMP, &
705  SFC_albedo, &
706  SFC_Z0M, &
707  SFC_Z0H, &
708  SFC_Z0E, &
709  SFLX_MW, &
710  SFLX_MU, &
711  SFLX_MV, &
712  SFLX_SH, &
713  SFLX_LH, &
714  SFLX_GH, &
715  SFLX_evap, &
716  U10, &
717  V10, &
718  T2, &
719  Q2, &
720  countup )
721  implicit none
722 
723  ! arguments
724  real(RP), intent(in) :: SFC_TEMP (ia,ja)
725  real(RP), intent(in) :: SFC_albedo(ia,ja,2)
726  real(RP), intent(in) :: SFC_Z0M (ia,ja)
727  real(RP), intent(in) :: SFC_Z0H (ia,ja)
728  real(RP), intent(in) :: SFC_Z0E (ia,ja)
729  real(RP), intent(in) :: SFLX_MW (ia,ja)
730  real(RP), intent(in) :: SFLX_MU (ia,ja)
731  real(RP), intent(in) :: SFLX_MV (ia,ja)
732  real(RP), intent(in) :: SFLX_SH (ia,ja)
733  real(RP), intent(in) :: SFLX_LH (ia,ja)
734  real(RP), intent(in) :: SFLX_GH (ia,ja)
735  real(RP), intent(in) :: SFLX_evap (ia,ja)
736  real(RP), intent(in) :: U10 (ia,ja)
737  real(RP), intent(in) :: V10 (ia,ja)
738  real(RP), intent(in) :: T2 (ia,ja)
739  real(RP), intent(in) :: Q2 (ia,ja)
740 
741  logical, intent(in) :: countup
742 
743  ! works
744  integer :: i, j
745  !---------------------------------------------------------------------------
746 
747  do j = js, je
748  do i = is, ie
749  urb_sfc_temp(i,j) = urb_sfc_temp(i,j) * cnt_puturb + sfc_temp(i,j)
750  urb_sfc_albedo(i,j,i_lw) = urb_sfc_albedo(i,j,i_lw) * cnt_puturb + sfc_albedo(i,j,i_lw)
751  urb_sfc_albedo(i,j,i_sw) = urb_sfc_albedo(i,j,i_sw) * cnt_puturb + sfc_albedo(i,j,i_sw)
752  urb_sfc_z0m(i,j) = urb_sfc_z0m(i,j) * cnt_puturb + sfc_z0m(i,j)
753  urb_sfc_z0h(i,j) = urb_sfc_z0h(i,j) * cnt_puturb + sfc_z0h(i,j)
754  urb_sfc_z0e(i,j) = urb_sfc_z0e(i,j) * cnt_puturb + sfc_z0e(i,j)
755  urb_sflx_mw(i,j) = urb_sflx_mw(i,j) * cnt_puturb + sflx_mw(i,j)
756  urb_sflx_mu(i,j) = urb_sflx_mu(i,j) * cnt_puturb + sflx_mu(i,j)
757  urb_sflx_mv(i,j) = urb_sflx_mv(i,j) * cnt_puturb + sflx_mv(i,j)
758  urb_sflx_sh(i,j) = urb_sflx_sh(i,j) * cnt_puturb + sflx_sh(i,j)
759  urb_sflx_lh(i,j) = urb_sflx_lh(i,j) * cnt_puturb + sflx_lh(i,j)
760  urb_sflx_gh(i,j) = urb_sflx_gh(i,j) * cnt_puturb + sflx_gh(i,j)
761  urb_sflx_evap(i,j) = urb_sflx_evap(i,j) * cnt_puturb + sflx_evap(i,j)
762  urb_u10(i,j) = urb_u10(i,j) * cnt_puturb + u10(i,j)
763  urb_v10(i,j) = urb_v10(i,j) * cnt_puturb + v10(i,j)
764  urb_t2(i,j) = urb_t2(i,j) * cnt_puturb + t2(i,j)
765  urb_q2(i,j) = urb_q2(i,j) * cnt_puturb + q2(i,j)
766  enddo
767  enddo
768 
769  do j = js, je
770  do i = is, ie
771  urb_sfc_temp(i,j) = urb_sfc_temp(i,j) / ( cnt_puturb + 1.0_rp )
772  urb_sfc_albedo(i,j,i_lw) = urb_sfc_albedo(i,j,i_lw) / ( cnt_puturb + 1.0_rp )
773  urb_sfc_albedo(i,j,i_sw) = urb_sfc_albedo(i,j,i_sw) / ( cnt_puturb + 1.0_rp )
774  urb_sfc_z0m(i,j) = urb_sfc_z0m(i,j) / ( cnt_puturb + 1.0_rp )
775  urb_sfc_z0h(i,j) = urb_sfc_z0h(i,j) / ( cnt_puturb + 1.0_rp )
776  urb_sfc_z0e(i,j) = urb_sfc_z0e(i,j) / ( cnt_puturb + 1.0_rp )
777  urb_sflx_mw(i,j) = urb_sflx_mw(i,j) / ( cnt_puturb + 1.0_rp )
778  urb_sflx_mu(i,j) = urb_sflx_mu(i,j) / ( cnt_puturb + 1.0_rp )
779  urb_sflx_mv(i,j) = urb_sflx_mv(i,j) / ( cnt_puturb + 1.0_rp )
780  urb_sflx_sh(i,j) = urb_sflx_sh(i,j) / ( cnt_puturb + 1.0_rp )
781  urb_sflx_lh(i,j) = urb_sflx_lh(i,j) / ( cnt_puturb + 1.0_rp )
782  urb_sflx_gh(i,j) = urb_sflx_gh(i,j) / ( cnt_puturb + 1.0_rp )
783  urb_sflx_evap(i,j) = urb_sflx_evap(i,j) / ( cnt_puturb + 1.0_rp )
784  urb_u10(i,j) = urb_u10(i,j) / ( cnt_puturb + 1.0_rp )
785  urb_v10(i,j) = urb_v10(i,j) / ( cnt_puturb + 1.0_rp )
786  urb_t2(i,j) = urb_t2(i,j) / ( cnt_puturb + 1.0_rp )
787  urb_q2(i,j) = urb_q2(i,j) / ( cnt_puturb + 1.0_rp )
788  enddo
789  enddo
790 
791  if( countup ) then
792  cnt_puturb = cnt_puturb + 1.0_rp
793  end if
794 
795  return
796  end subroutine cpl_puturb
797 
798  !-----------------------------------------------------------------------------
799  subroutine cpl_getsfc_atm( &
800  SFC_TEMP, &
801  SFC_albedo, &
802  SFC_Z0M, &
803  SFC_Z0H, &
804  SFC_Z0E, &
805  SFLX_MW, &
806  SFLX_MU, &
807  SFLX_MV, &
808  SFLX_SH, &
809  SFLX_LH, &
810  SFLX_GH, &
811  SFLX_QTRC, &
812  U10, &
813  V10, &
814  T2, &
815  Q2 )
816  use scale_landuse, only: &
817  fact_ocean => landuse_fact_ocean, &
818  fact_land => landuse_fact_land, &
819  fact_urban => landuse_fact_urban
820  implicit none
821 
822  real(RP), intent(out) :: SFC_TEMP (ia,ja)
823  real(RP), intent(out) :: SFC_albedo(ia,ja,2)
824  real(RP), intent(out) :: SFC_Z0M (ia,ja)
825  real(RP), intent(out) :: SFC_Z0H (ia,ja)
826  real(RP), intent(out) :: SFC_Z0E (ia,ja)
827  real(RP), intent(out) :: SFLX_MW (ia,ja)
828  real(RP), intent(out) :: SFLX_MU (ia,ja)
829  real(RP), intent(out) :: SFLX_MV (ia,ja)
830  real(RP), intent(out) :: SFLX_SH (ia,ja)
831  real(RP), intent(out) :: SFLX_LH (ia,ja)
832  real(RP), intent(out) :: SFLX_GH (ia,ja)
833  real(RP), intent(out) :: SFLX_QTRC (ia,ja,qa)
834  real(RP), intent(out) :: U10 (ia,ja)
835  real(RP), intent(out) :: V10 (ia,ja)
836  real(RP), intent(out) :: T2 (ia,ja)
837  real(RP), intent(out) :: Q2 (ia,ja)
838 
839  integer :: i, j, iq
840  !---------------------------------------------------------------------------
841 
842  do j = js, je
843  do i = is, ie
844  do iq = 1, qa
845  sflx_qtrc(i,j,iq) = 0.0_rp ! tentative
846  enddo
847 
848  sfc_temp(i,j) = fact_ocean(i,j) * ocn_sfc_temp(i,j) &
849  + fact_land(i,j) * lnd_sfc_temp(i,j) &
850  + fact_urban(i,j) * urb_sfc_temp(i,j)
851 
852  sfc_albedo(i,j,i_lw) = fact_ocean(i,j) * ocn_sfc_albedo(i,j,i_lw) &
853  + fact_land(i,j) * lnd_sfc_albedo(i,j,i_lw) &
854  + fact_urban(i,j) * urb_sfc_albedo(i,j,i_lw)
855 
856  sfc_albedo(i,j,i_sw) = fact_ocean(i,j) * ocn_sfc_albedo(i,j,i_sw) &
857  + fact_land(i,j) * lnd_sfc_albedo(i,j,i_sw) &
858  + fact_urban(i,j) * urb_sfc_albedo(i,j,i_sw)
859 
860  sfc_z0m(i,j) = fact_ocean(i,j) * ocn_sfc_z0m(i,j) &
861  + fact_land(i,j) * lnd_sfc_z0m(i,j) &
862  + fact_urban(i,j) * urb_sfc_z0m(i,j)
863 
864  sfc_z0h(i,j) = fact_ocean(i,j) * ocn_sfc_z0h(i,j) &
865  + fact_land(i,j) * lnd_sfc_z0h(i,j) &
866  + fact_urban(i,j) * urb_sfc_z0h(i,j)
867 
868  sfc_z0e(i,j) = fact_ocean(i,j) * ocn_sfc_z0e(i,j) &
869  + fact_land(i,j) * lnd_sfc_z0e(i,j) &
870  + fact_urban(i,j) * urb_sfc_z0e(i,j)
871 
872  sflx_mw(i,j) = fact_ocean(i,j) * ocn_sflx_mw(i,j) &
873  + fact_land(i,j) * lnd_sflx_mw(i,j) &
874  + fact_urban(i,j) * urb_sflx_mw(i,j)
875 
876  sflx_mu(i,j) = fact_ocean(i,j) * ocn_sflx_mu(i,j) &
877  + fact_land(i,j) * lnd_sflx_mu(i,j) &
878  + fact_urban(i,j) * urb_sflx_mu(i,j)
879 
880  sflx_mv(i,j) = fact_ocean(i,j) * ocn_sflx_mv(i,j) &
881  + fact_land(i,j) * lnd_sflx_mv(i,j) &
882  + fact_urban(i,j) * urb_sflx_mv(i,j)
883 
884  sflx_sh(i,j) = fact_ocean(i,j) * ocn_sflx_sh(i,j) &
885  + fact_land(i,j) * lnd_sflx_sh(i,j) &
886  + fact_urban(i,j) * urb_sflx_sh(i,j)
887 
888  sflx_lh(i,j) = fact_ocean(i,j) * ocn_sflx_lh(i,j) &
889  + fact_land(i,j) * lnd_sflx_lh(i,j) &
890  + fact_urban(i,j) * urb_sflx_lh(i,j)
891 
892  sflx_gh(i,j) = fact_ocean(i,j) * ocn_sflx_wh(i,j) &
893  + fact_land(i,j) * lnd_sflx_gh(i,j) &
894  + fact_urban(i,j) * urb_sflx_gh(i,j)
895 
896  sflx_qtrc(i,j,i_qv) = fact_ocean(i,j) * ocn_sflx_evap(i,j) &
897  + fact_land(i,j) * lnd_sflx_evap(i,j) &
898  + fact_urban(i,j) * urb_sflx_evap(i,j)
899 
900  u10(i,j) = fact_ocean(i,j) * ocn_u10(i,j) &
901  + fact_land(i,j) * lnd_u10(i,j) &
902  + fact_urban(i,j) * urb_u10(i,j)
903 
904  v10(i,j) = fact_ocean(i,j) * ocn_v10(i,j) &
905  + fact_land(i,j) * lnd_v10(i,j) &
906  + fact_urban(i,j) * urb_v10(i,j)
907 
908  t2(i,j) = fact_ocean(i,j) * ocn_t2(i,j) &
909  + fact_land(i,j) * lnd_t2(i,j) &
910  + fact_urban(i,j) * urb_t2(i,j)
911 
912  q2(i,j) = fact_ocean(i,j) * ocn_q2(i,j) &
913  + fact_land(i,j) * lnd_q2(i,j) &
914  + fact_urban(i,j) * urb_q2(i,j)
915  enddo
916  enddo
917 
918  cnt_putocn = 0.0_rp
919  cnt_putlnd = 0.0_rp
920  cnt_puturb = 0.0_rp
921 
922  return
923  end subroutine cpl_getsfc_atm
924 
925  !-----------------------------------------------------------------------------
926  subroutine cpl_getatm_ocn( &
927  TEMP, &
928  PRES, &
929  W, &
930  U, &
931  V, &
932  DENS, &
933  QV, &
934  PBL, &
935  SFC_PRES, &
936  SFLX_rad_dn, &
937  cosSZA, &
938  SFLX_rain, &
939  SFLX_snow )
940  implicit none
941 
942  real(RP), intent(out) :: TEMP (ia,ja)
943  real(RP), intent(out) :: PRES (ia,ja)
944  real(RP), intent(out) :: W (ia,ja)
945  real(RP), intent(out) :: U (ia,ja)
946  real(RP), intent(out) :: V (ia,ja)
947  real(RP), intent(out) :: DENS (ia,ja)
948  real(RP), intent(out) :: QV (ia,ja)
949  real(RP), intent(out) :: PBL (ia,ja)
950  real(RP), intent(out) :: SFC_PRES (ia,ja)
951  real(RP), intent(out) :: SFLX_rad_dn(ia,ja,2,2)
952  real(RP), intent(out) :: cosSZA (ia,ja)
953  real(RP), intent(out) :: SFLX_rain (ia,ja)
954  real(RP), intent(out) :: SFLX_snow (ia,ja)
955 
956  integer :: i, j
957  !---------------------------------------------------------------------------
958 
959 !OCL XFILL
960  do j = js, je
961  do i = is, ie
962  temp(i,j) = ocn_atm_temp(i,j)
963  pres(i,j) = ocn_atm_pres(i,j)
964  w(i,j) = ocn_atm_w(i,j)
965  u(i,j) = ocn_atm_u(i,j)
966  v(i,j) = ocn_atm_v(i,j)
967  dens(i,j) = ocn_atm_dens(i,j)
968  qv(i,j) = ocn_atm_qv(i,j)
969  pbl(i,j) = ocn_atm_pbl(i,j)
970  sfc_pres(i,j) = ocn_atm_sfc_pres(i,j)
971  sflx_rad_dn(i,j,:,:) = ocn_atm_sflx_rad_dn(i,j,:,:)
972  cossza(i,j) = ocn_atm_cossza(i,j)
973  sflx_rain(i,j) = ocn_atm_sflx_rain(i,j)
974  sflx_snow(i,j) = ocn_atm_sflx_snow(i,j)
975  enddo
976  enddo
977 
978  cnt_putatm_ocn = 0.0_rp
979 
980  return
981  end subroutine cpl_getatm_ocn
982 
983  !-----------------------------------------------------------------------------
984  subroutine cpl_getatm_lnd( &
985  TEMP, &
986  PRES, &
987  W, &
988  U, &
989  V, &
990  DENS, &
991  QV, &
992  PBL, &
993  SFC_PRES, &
994  SFLX_rad_dn, &
995  cosSZA, &
996  SFLX_rain, &
997  SFLX_snow )
998  implicit none
999 
1000  real(RP), intent(out) :: TEMP (ia,ja)
1001  real(RP), intent(out) :: PRES (ia,ja)
1002  real(RP), intent(out) :: W (ia,ja)
1003  real(RP), intent(out) :: U (ia,ja)
1004  real(RP), intent(out) :: V (ia,ja)
1005  real(RP), intent(out) :: DENS (ia,ja)
1006  real(RP), intent(out) :: QV (ia,ja)
1007  real(RP), intent(out) :: PBL (ia,ja)
1008  real(RP), intent(out) :: SFC_PRES (ia,ja)
1009  real(RP), intent(out) :: SFLX_rad_dn(ia,ja,2,2)
1010  real(RP), intent(out) :: cosSZA (ia,ja)
1011  real(RP), intent(out) :: SFLX_rain (ia,ja)
1012  real(RP), intent(out) :: SFLX_snow (ia,ja)
1013 
1014  integer :: i, j
1015  !---------------------------------------------------------------------------
1016 
1017 !OCL XFILL
1018  do j = js, je
1019  do i = is, ie
1020  temp(i,j) = lnd_atm_temp(i,j)
1021  pres(i,j) = lnd_atm_pres(i,j)
1022  w(i,j) = lnd_atm_w(i,j)
1023  u(i,j) = lnd_atm_u(i,j)
1024  v(i,j) = lnd_atm_v(i,j)
1025  dens(i,j) = lnd_atm_dens(i,j)
1026  qv(i,j) = lnd_atm_qv(i,j)
1027  pbl(i,j) = lnd_atm_pbl(i,j)
1028  sfc_pres(i,j) = lnd_atm_sfc_pres(i,j)
1029  sflx_rad_dn(i,j,:,:) = lnd_atm_sflx_rad_dn(i,j,:,:)
1030  cossza(i,j) = lnd_atm_cossza(i,j)
1031  sflx_rain(i,j) = lnd_atm_sflx_rain(i,j)
1032  sflx_snow(i,j) = lnd_atm_sflx_snow(i,j)
1033  enddo
1034  enddo
1035 
1036  cnt_putatm_lnd = 0.0_rp
1037 
1038  return
1039  end subroutine cpl_getatm_lnd
1040 
1041  !-----------------------------------------------------------------------------
1042  subroutine cpl_getatm_urb( &
1043  TEMP, &
1044  PRES, &
1045  W, &
1046  U, &
1047  V, &
1048  DENS, &
1049  QV, &
1050  PBL, &
1051  SFC_PRES, &
1052  SFLX_rad_dn, &
1053  cosSZA, &
1054  SFLX_rain, &
1055  SFLX_snow )
1056  implicit none
1057 
1058  real(RP), intent(out) :: TEMP (ia,ja)
1059  real(RP), intent(out) :: PRES (ia,ja)
1060  real(RP), intent(out) :: W (ia,ja)
1061  real(RP), intent(out) :: U (ia,ja)
1062  real(RP), intent(out) :: V (ia,ja)
1063  real(RP), intent(out) :: DENS (ia,ja)
1064  real(RP), intent(out) :: QV (ia,ja)
1065  real(RP), intent(out) :: PBL (ia,ja)
1066  real(RP), intent(out) :: SFC_PRES (ia,ja)
1067  real(RP), intent(out) :: SFLX_rad_dn(ia,ja,2,2)
1068  real(RP), intent(out) :: cosSZA (ia,ja)
1069  real(RP), intent(out) :: SFLX_rain (ia,ja)
1070  real(RP), intent(out) :: SFLX_snow (ia,ja)
1071 
1072  integer :: i, j
1073  !---------------------------------------------------------------------------
1074 
1075 !OCL XFILL
1076  do j = js, je
1077  do i = is, ie
1078  temp(i,j) = urb_atm_temp(i,j)
1079  pres(i,j) = urb_atm_pres(i,j)
1080  w(i,j) = urb_atm_w(i,j)
1081  u(i,j) = urb_atm_u(i,j)
1082  v(i,j) = urb_atm_v(i,j)
1083  dens(i,j) = urb_atm_dens(i,j)
1084  qv(i,j) = urb_atm_qv(i,j)
1085  pbl(i,j) = urb_atm_pbl(i,j)
1086  sfc_pres(i,j) = urb_atm_sfc_pres(i,j)
1087  sflx_rad_dn(i,j,:,:) = urb_atm_sflx_rad_dn(i,j,:,:)
1088  cossza(i,j) = urb_atm_cossza(i,j)
1089  sflx_rain(i,j) = urb_atm_sflx_rain(i,j)
1090  sflx_snow(i,j) = urb_atm_sflx_snow(i,j)
1091  enddo
1092  enddo
1093 
1094  cnt_putatm_urb = 0.0_rp
1095 
1096  return
1097  end subroutine cpl_getatm_urb
1098 
1099 end module mod_cpl_vars
integer, public is
start point of inner domain: x, local
module DEBUG
Definition: scale_debug.F90:13
real(rp), dimension(:,:), allocatable, public ocn_sflx_evap
real(rp), dimension(:,:), allocatable, public urb_atm_w
real(rp), dimension(:,:), allocatable, public urb_sflx_mv
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:98
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:59
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 x whole cells (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
integer, public i_qv
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
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:99
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
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)
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 y whole cells (local, with HALO)
real(rp), public cnt_putatm_urb