SCALE-RM
scale_atmos_phy_rd_mm5sw.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
13 !-------------------------------------------------------------------------------
15  !-----------------------------------------------------------------------------
16  !
17  !++ used modules
18  !
19  use scale_precision
20  use scale_stdio
21  use scale_prof
23  use scale_tracer
24 
25  use scale_atmos_phy_rd_common, only: &
26  i_sw, &
27  i_lw, &
28  i_dn, &
29  i_up
30 
31  !-----------------------------------------------------------------------------
32  implicit none
33  private
34  !-----------------------------------------------------------------------------
35  !
36  !++ Public procedure
37  !
38  public :: swrad
39  public :: swinit
40 
41  !-----------------------------------------------------------------------------
42  !
43  !++ Public parameters & variables
44  !
45  !-----------------------------------------------------------------------------
46  !
47  !++ Private procedure
48  !
49  !-----------------------------------------------------------------------------
50  !
51  !++ Private parameters & variables
52  !
53  REAL(RP) :: CSSCA
54 
55  CONTAINS
56 
57 !------------------------------------------------------------------
58  SUBROUTINE swrad(dt,RTHRATEN,SDOWN3D,GSW,XLAT,XLONG,ALBEDO, &
59  rho_phy,T3D,P3D,pi3D,dz8w, &
60  solins,cosSZA, &
61  QV3D,QC3D,QR3D,QI3D,QS3D,QG3D, &
62  F_QV,F_QC,F_QR,F_QI,F_QS,F_QG, &
63  icloud,warm_rain )
64 
65  use scale_const, only: &
66  grav => const_grav, & ! 9.81
67  rdry => const_rdry, & ! 287.
68  cpdry => const_cpdry ! 1004.6
69  use scale_atmos_solarins, only: &
70  solarins => atmos_solarins_constant ! 1360.xxx (W/m2)
71 !------------------------------------------------------------------
72  IMPLICIT NONE
73 !------------------------------------------------------------------
74 
75  real(DP), INTENT(IN) :: dt
76  real(RP), INTENT(INOUT) :: rthraten(ia,ka,ja)
77  real(RP), INTENT(INOUT) :: sdown3d(ia,ka,ja)
78  real(RP), INTENT(INOUT) :: gsw(ia,ja)
79 
80  real(RP), INTENT(IN) :: xlat(ia,ja), xlong(ia,ja), albedo(ia, ja)
81  real(RP), INTENT(IN) :: rho_phy(ia,ka,ja)
82  real(RP), INTENT(IN) :: p3d(ia,ka,ja), &
83  t3d(ia,ka,ja), &
84  pi3d(ia,ka,ja), &
85  dz8w(ia,ka,ja)
86  real(RP), INTENT(IN) :: solins(ia,ja),cossza(ia,ja)
87  real(RP), OPTIONAL, INTENT(IN) :: qv3d (ia,ka,ja), &
88  qc3d (ia,ka,ja), &
89  qr3d (ia,ka,ja), &
90  qi3d (ia,ka,ja), &
91  qs3d (ia,ka,ja), &
92  qg3d (ia,ka,ja)
93 
94  LOGICAL, OPTIONAL, INTENT(IN ) :: f_qv,f_qc,f_qr,f_qi,f_qs,f_qg
95  INTEGER, INTENT(IN ) :: icloud
96  LOGICAL, INTENT(IN ) :: warm_rain
97 
98  !real, INTENT(IN ) :: RADFRQ,DEGRAD,XTIME,DECLIN
99 
100 
101  integer :: its,ite,jts,jte,kts,kte
102  real(RP) :: r, cp, g, solcon
103 
104 ! LOCAL VARS
105 
106  real(RP), DIMENSION( KS:KE ) :: tten1d, &
107  rho01d, &
108  p1d, &
109  dz, &
110  t1d, &
111  qv1d, &
112  qc1d, &
113  qr1d, &
114  qi1d, &
115  qs1d, &
116  qg1d
117 
118  real(RP) :: sdown1d(ks:ke+1)
119 !
120  real(RP) :: xlat0,xlong0,alb0,gsw0,cossza0,solins0
121  real(RP) :: aer_dry1(ks:ke),aer_water1(ks:ke)
122 
123 !
124  INTEGER :: i,j,k,nk
125  LOGICAL :: predicate , do_topo_shading
126 
127 
128 !------------------------------------------------------------------
129 
130  r = rdry
131  cp = cpdry
132  g = grav
133  solcon = solarins
134  its = is ; ite = ie
135  jts = js ; jte = je
136  kts = ks ; kte = ke
137 
138 
139  j_loop: DO j=jts,jte
140  i_loop: DO i=its,ite
141 
142 ! reverse vars
143  DO k=ks,ke
144  qv1d(k)=0.
145  qc1d(k)=0.
146  qr1d(k)=0.
147  qi1d(k)=0.
148  qs1d(k)=0.
149  qg1d(k)=0.
150  ENDDO
151 
152  DO k=ks,ke
153  !! NK=kme-1-K+kms
154  nk=ke-(k-ks)
155  tten1d(k)=0.
156 
157  t1d(k)=t3d(i,nk,j)
158  p1d(k)=p3d(i,nk,j)
159  rho01d(k)=rho_phy(i,nk,j)
160  dz(k)=dz8w(i,nk,j)
161  ENDDO
162 
163  !IF( PRESENT(pm2_5_dry) .AND. PRESENT(pm2_5_water) )THEN
164  ! DO K=kts,kte
165  ! NK=kme-1-K+kms
166  ! aer_dry1(k) = pm2_5_dry(i,nk,j)
167  ! aer_water1(k) = pm2_5_water(i,nk,j)
168  ! ENDDO
169  !ELSE
170  do k=ks,ke
171  aer_dry1(k) = 0.0
172  aer_water1(k) = 0.0
173  enddo
174  !ENDIF
175 
176  IF (PRESENT(f_qv) .AND. PRESENT(qv3d)) THEN
177  IF (f_qv) THEN
178  do k=ks,ke
179  !NK=kme-1-K+kms
180  nk=ke-(k-ks)
181  qv1d(k)=qv3d(i,nk,j)
182  qv1d(k)=max(0.0_rp,qv1d(k))
183  ENDDO
184  ENDIF
185  ENDIF
186 
187  IF (PRESENT(f_qc) .AND. PRESENT(qc3d)) THEN
188  IF (f_qc) THEN
189  do k=ks,ke
190  !NK=kme-1-K+kms
191  nk=ke-(k-ks)
192  qc1d(k)=qc3d(i,nk,j)
193  qc1d(k)=max(0.0_rp,qc1d(k))
194  enddo
195  ENDIF
196  ENDIF
197 
198  IF (PRESENT(f_qr) .AND. PRESENT(qr3d)) THEN
199  IF (f_qr) THEN
200  do k=kts,kte
201  !NK=kme-1-K+kms
202  nk=ke-(k-ks)
203  qr1d(k)=qr3d(i,nk,j)
204  qr1d(k)=max(0.0_rp,qr1d(k))
205  enddo
206  ENDIF
207  ENDIF
208 
209 !
210  IF ( PRESENT( f_qi ) ) THEN
211  predicate = f_qi
212  ELSE
213  predicate = .false.
214  ENDIF
215 
216  IF ( predicate .AND. PRESENT( qi3d ) ) THEN
217  do k=ks,ke
218  !NK=kme-1-K+kms
219  nk=ke-(k-ks)
220  qi1d(k)=qi3d(i,nk,j)
221  qi1d(k)=max(0.0_rp,qi1d(k))
222  enddo
223  ELSE
224  IF (.NOT. warm_rain) THEN
225  do k=ks,ke
226  IF(t1d(k) < 273.15) THEN
227  qi1d(k)=qc1d(k)
228  qc1d(k)=0.
229  qs1d(k)=qr1d(k)
230  qr1d(k)=0.
231  ENDIF
232  enddo
233  ENDIF
234  ENDIF
235 
236  IF (PRESENT(f_qs) .AND. PRESENT(qs3d)) THEN
237  IF (f_qs) THEN
238  do k=ks,ke
239  !NK=kme-1-K+kms
240  nk=ke-(k-ks)
241  qs1d(k)=qs3d(i,nk,j)
242  qs1d(k)=max(0.0_rp,qs1d(k))
243  enddo
244  ENDIF
245  ENDIF
246 
247  IF (PRESENT(f_qg) .AND. PRESENT(qg3d)) THEN
248  IF (f_qg) THEN
249  do k=ks,ke
250  !NK=kme-1-K+kms
251  nk=ke-(k-ks)
252  qg1d(k)=qg3d(i,nk,j)
253  qg1d(k)=max(0.0_rp,qg1d(k))
254  enddo
255  ENDIF
256  ENDIF
257 
258  xlat0=xlat(i,j)
259  xlong0=xlong(i,j)
260  solins0=solins(i,j)
261  cossza0=cossza(i,j)
262  alb0=albedo(i,j)
263 ! slope code removed - factor now done in surface driver
264  CALL swpara(tten1d,sdown1d,gsw0,alb0,cossza0, &
265  t1d,qv1d,qc1d,qr1d,qi1d,qs1d,qg1d,p1d, &
266  rho01d,dz, &
267  r,cp,g,solins0, &
268  xlat0,xlong0, &
269  icloud,aer_dry1,aer_water1, &
270  kts,kte )
271  gsw(i,j)=gsw0
272 
273 ! DO K=kts,kte
274 ! NK=kme-1-K+kms
275 ! NK=KA-1-K+1
276 ! RTHRATEN(I,K,J)=RTHRATEN(I,K,J)+TTEN1D(NK)/pi3D(I,K,J)
277 ! SDOWN3D(I,K,J)=SDOWN1D(NK)
278 ! if(k==kte)then
279 ! NK=KA-1-(K+1)+1
280 ! SDOWN3D(I,K+1,J)=SDOWN1D(NK)
281 ! endif
282 ! ENDDO
283 
284  do k=ks,ke
285  nk=ke-(k-ks)
286  rthraten(i,k,j)=rthraten(i,k,j)+tten1d(nk)/pi3d(i,k,j)
287  sdown3d(i,k,j)=sdown1d(nk)
288  !print *,KS,KE,K,"<=",NK
289  enddo
290  k=ks-1 ; nk=ke+1
291  sdown3d(i,ks-1,j)=sdown1d(ke+1)
292  !print *,KS,KE,K,"<=",NK
293 !
294  ENDDO i_loop
295  ENDDO j_loop
296 
297  return
298  END SUBROUTINE swrad
299 
300 !------------------------------------------------------------------
301  SUBROUTINE swpara(TTEN,SDOWN,GSW,ALBEDO,cosSZA, &
302  T,QV,QC,QR,QI,QS,QG,P, &
303  RHO0, DZ, &
304  R,CP,G,solins, &
305  XXLAT,XXLON, &
306  ICLOUD,aer_dry1,aer_water1, &
307  kts,kte )
309  use scale_time, only: &
310  nowdate => time_nowdate !< current time [yyyy mm dd hh mm ss]
311  use scale_const, only: &
312  d2r => const_d2r ! degree to radian
313 !------------------------------------------------------------------
314 ! TO CALCULATE SHORT-WAVE ABSORPTION AND SCATTERING IN CLEAR
315 ! AIR AND REFLECTION AND ABSORPTION IN CLOUD LAYERS (STEPHENS,
316 ! 1984)
317 ! CHANGES:
318 ! REDUCE EFFECTS OF ICE CLOUDS AND PRECIP ON LIQUID WATER PATH
319 ! ADD EFFECT OF GRAUPEL
320 !------------------------------------------------------------------
321  IMPLICIT NONE
322 
323  INTEGER, INTENT(IN ) :: kts,kte
324 
325  real(RP), DIMENSION( kts:kte ), INTENT(IN ) :: &
326  RHO0, &
327  T, &
328  P, &
329  DZ, &
330  QV, &
331  QC, &
332  QR, &
333  QI, &
334  QS, &
335  QG
336 
337  real(RP), DIMENSION( kts:kte ), INTENT(INOUT) :: TTEN
338  real(RP), DIMENSION( kts:kte+1 ), INTENT(OUT) :: SDOWN
339 !
340  real(RP), INTENT(IN ) :: R,CP,G
341  real(RP), INTENT(IN) :: ALBEDO,cosSZA,solins
342  integer , INTENT(IN) :: icloud
343  real(RP), INTENT(INOUT) :: GSW
344  real(RP), INTENT(IN ) :: XXLAT,XXLON
345 
346  !real, INTENT(IN ) :: GMT,RADFRQ,XTIME,XLAT,XLONG,DEGRAD
347 !
348 ! For slope-dependent radiation
349 ! INTEGER, OPTIONAL, INTENT(IN) :: slope_rad,shadow
350 ! real, OPTIONAL, INTENT(IN) :: slp_azi,slope
351 !
352 
353 ! LOCAL VARS
354 
355  real(RP) :: XLWP( kts:kte )
356  real(RP) :: XATP( kts:kte )
357  real(RP) :: XWVP( kts:kte )
358  real(RP) :: aer_dry1( kts:kte )
359  real(RP) :: aer_water1( kts:kte )
360  real(RP) :: RO( kts:kte )
361 !
362  real(RP) :: ALBTAB(4,5), ABSTAB(4,5)
363  real(RP) :: XMUVAL(4)
364 
365  real(RP) :: beta
366  real(RP) :: DayOfYear
367 !------------------------------------------------------------------
368 
369  DATA albtab/0.0,0.0,0.0,0.0, &
370  69.0,58.0,40.0,15.0, &
371  90.0,80.0,70.0,60.0, &
372  94.0,90.0,82.0,78.0, &
373  96.0,92.0,85.0,80.0/
374 
375  DATA abstab/0.0,0.0,0.0,0.0, &
376  0.0,2.5,4.0,5.0, &
377  0.0,2.6,7.0,10.0, &
378  0.0,3.3,10.0,14.0, &
379  0.0,3.7,10.0,15.0/
380 
381  DATA xmuval/0.0,0.2,0.5,1.0/
382 
383  real(RP) :: bext340, absc, alba, alw, csza,dabsa,dsca,dabs
384  real(RP) :: bexth2o, dscld, ff,oldalb,oldabs,oldabc
385  real(RP) :: soltop, totabs, ugcm, uv,xabs,xabsa,wv
386  real(RP) :: wgm, xalb, xi, xsca, xmu,xabsc,trans0,yj
387  real(RP) :: HRANG,XT24,TLOCTM,tloc,dsec,lon,lat
388  real(RP) :: DECLIN
389  real(RP) :: ww
390  INTEGER :: iil,ii,jjl,ju,k,iu
391 
392 ! For slope-dependent radiation
393 
394  real(RP) :: diffuse_frac, corr_fac, csza_slp
395 
396  gsw=0.0
397  sdown=0.0
398  bext340=5.e-6_rp
399  bexth2o=5.e-6_rp
400  soltop = solins
401 
402  !XT24=MOD(XTIME+RADFRQ*0.5,1440.)
403  !TLOCTM=GMT+XT24/60.+XLONG/15.
404  !HRANG=15.*(TLOCTM-12.)*DEGRAD
405  !XXLAT=XLAT*DEGRAD
406  !CSZA=SIN(XXLAT)*SIN(DECLIN)+COS(XXLAT)*COS(DECLIN)*COS(HRANG)
407 
408  ! local time
409  lat = xxlat / d2r
410  lon = xxlon / d2r
411 
412  tloc = mod((nowdate(4) + int(lon/15.0_rp)),24 )
413  dsec = real(NOWDATE(5)*60 + NOWDATE(6)) / 60.0_RP /60.0_RP
414  tloctm = real(NOWDATE(4)) + LON/15.0_RP + dsec
415  hrang = 15.*(tloctm-12.)*d2r
416 
417  csza=cossza
418 
419 ! RETURN IF NIGHT
420  IF(csza <= 1.0e-9_rp)GOTO 7
421 !
422  DO k=kts, kte
423 
424 ! P in the unit of 10mb
425  ro(k)=p(k)/(r*t(k))
426  xwvp(k)=ro(k)*qv(k)*dz(k)*1000.
427 ! KG/M**2
428  xatp(k)=ro(k)*dz(k)
429  ENDDO
430 !
431 ! G/M**2
432 ! REDUCE WEIGHT OF LIQUID AND ICE IN SHORT-WAVE SCHEME
433 ! ADD GRAUPEL EFFECT (ASSUMED SAME AS RAIN)
434 !
435  IF (icloud==0)THEN
436  DO k=kts, kte
437  xlwp(k)=0.
438  ENDDO
439  ELSE
440  DO k=kts, kte
441  xlwp(k)=ro(k)*1000.*dz(k)*(qc(k)+0.1*qi(k)+0.05* &
442  qr(k)+0.02*qs(k)+0.05*qg(k))
443  ENDDO
444  ENDIF
445 !
446  xmu=csza
447  ! SDOWN(1)=SOLTOP*XMU !adachi
448 
449  sdown(kts)=soltop * max(cossza,0.0_rp)
450 
451 ! SET WW (G/M**2) LIQUID WATER PATH INTEGRATED DOWN
452 ! SET UV (G/M**2) WATER VAPOR PATH INTEGRATED DOWN
453  ww=0.
454  uv=0.
455  oldalb=0.
456  oldabc=0.
457  totabs=0.
458 ! CONTRIBUTIONS DUE TO CLEAR AIR AND CLOUD
459  dsca=0.
460  dabs=0.
461  dscld=0.
462 !
463 ! CONTRIBUTION DUE TO AEROSOLS (FOR CHEMISTRY)
464  dabsa=0.
465 !
466  DO 200 k=kts,kte
467  ww=ww+xlwp(k)
468  uv=uv+xwvp(k)
469 ! WGM IS WW/COS(THETA) (G/M**2)
470 ! UGCM IS UV/COS(THETA) (G/CM**2)
471  wgm=ww/xmu
472  ugcm=uv*0.0001/xmu
473 !
474  oldabs=totabs
475 ! WATER VAPOR ABSORPTION AS IN LACIS AND HANSEN (1974)
476  totabs=2.9*ugcm/((1.+141.5*ugcm)**0.635+5.925*ugcm)
477 ! APPROXIMATE RAYLEIGH + AEROSOL SCATTERING
478 ! XSCA=1.E-5*XATP(K)/XMU
479 ! XSCA=(1.E-5*XATP(K)+aer_dry1(K)*bext340+aer_water1(K)*bexth2o)/XMU
480  beta=0.4*(1.0-xmu)+0.1
481 ! CSSCA - CLEAR-SKY SCATTERING SET FROM NAMELIST SWRAD_SCAT
482  xsca=(cssca*xatp(k)+beta*aer_dry1(k)*bext340*dz(k) &
483  +beta*aer_water1(k)*bexth2o*dz(k))/xmu
484 
485 ! LAYER VAPOR ABSORPTION DONE FIRST
486  ! XABS=(TOTABS-OLDABS)*(SDOWN(1)-DSCLD-DSCA-DABSA)/SDOWN(K)
487  xabs=(totabs-oldabs)*(sdown(kts)-dscld-dsca-dabsa)/sdown(k)
488 !rs AEROSOL ABSORB (would be elemental carbon). So far XABSA = 0.
489  xabsa=0.
490  IF(xabs<0.)xabs=0.
491 !
492  alw=log10(wgm+1.0_rp)
493  IF(alw>3.999)alw=3.999_rp
494 !
495  DO ii=1,3
496  IF(xmu>xmuval(ii))THEN
497  iil=ii
498  iu=ii+1
499  xi=(xmu-xmuval(ii))/(xmuval(ii+1)-xmuval(ii))+float(iil)
500  ENDIF
501  ENDDO
502 !
503  jjl = int(alw)+1
504  ju = jjl+1
505  yj = alw+1.
506 ! CLOUD ALBEDO
507  alba=(albtab(iu,ju)*(xi-iil)*(yj-jjl) &
508  +albtab(iil,ju)*(iu-xi)*(yj-jjl) &
509  +albtab(iu,jjl)*(xi-iil)*(ju-yj) &
510  +albtab(iil,jjl)*(iu-xi)*(ju-yj)) &
511  /((iu-iil)*(ju-jjl))
512 ! CLOUD ABSORPTION
513  absc=(abstab(iu,ju)*(xi-iil)*(yj-jjl) &
514  +abstab(iil,ju)*(iu-xi)*(yj-jjl) &
515  +abstab(iu,jjl)*(xi-iil)*(ju-yj) &
516  +abstab(iil,jjl)*(iu-xi)*(ju-yj)) &
517  /((iu-iil)*(ju-jjl))
518 ! LAYER ALBEDO AND ABSORPTION
519  !XALB=(ALBA-OLDALB)*(SDOWN(1)-DSCA-DABS)/SDOWN(K)
520  !XABSC=(ABSC-OLDABC)*(SDOWN(1)-DSCA-DABS)/SDOWN(K)
521  xalb=(alba-oldalb)*(sdown(kts)-dsca-dabs)/sdown(k)
522  xabsc=(absc-oldabc)*(sdown(kts)-dsca-dabs)/sdown(k)
523  IF(xalb<0.)xalb=0.
524  IF(xabsc<0.)xabsc=0.
525  dscld=dscld+(xalb+xabsc)*sdown(k)*0.01
526  dsca=dsca+xsca*sdown(k)
527  dabs=dabs+xabs*sdown(k)
528  dabsa=dabsa+xabsa*sdown(k)
529  oldalb=alba
530  oldabc=absc
531 ! LAYER TRANSMISSIVITY
532  trans0=100.0_rp-xalb-xabsc-xabs*100.0_rp-xsca*100.0_rp
533  IF(trans0<1.)THEN
534  ff=99.0_rp/(xalb+xabsc+xabs*100.0_rp+xsca*100.0_rp)
535  xalb=xalb*ff
536  xabsc=xabsc*ff
537  xabs=xabs*ff
538  xsca=xsca*ff
539  trans0=1.0_rp
540  ENDIF
541  sdown(k+1)=max(1.0e-9_rp,sdown(k)*trans0*0.01_rp)
542  tten(k)=sdown(k)*(xabsc+xabs*100._rp+xabsa*100._rp)*0.01_rp/(ro(k)*cp*dz(k))
543  200 CONTINUE
544 !
545  gsw=(1.-albedo)*sdown(kte+1)
546 
547 ! IF (PRESENT(slope_rad)) THEN
548 !! Slope-dependent solar radiation part
549 !
550 ! if (slope_rad==1) then
551 !
552 !! Parameterize diffuse fraction of global solar radiation as a function of the ratio between TOA radiation and surface global radiation
553 !
554 ! diffuse_frac = min(1.,1/(max(0.1,2.1-2.8*log(log(SDOWN(kts)/max(SDOWN(kte+1),1.e-3))))))
555 ! if ((slope==0).OR.(diffuse_frac==1).OR.(csza<1.e-2)) then ! no topographic effects when all radiation is diffuse or the sun is too close to the horizon
556 ! corr_fac = 1
557 ! goto 140
558 ! endif
559 !
560 !! cosine of zenith angle over sloping topography
561 !
562 ! csza_slp = ((SIN(XXLAT)*COS(HRANG))* &
563 ! (-cos(slp_azi)*sin(slope))-SIN(HRANG)*(sin(slp_azi)*sin(slope))+ &
564 ! (COS(XXLAT)*COS(HRANG))*cos(slope))* &
565 ! COS(DECLIN)+(COS(XXLAT)*(cos(slp_azi)*sin(slope))+ &
566 ! SIN(XXLAT)*cos(slope))*SIN(DECLIN)
567 !
568 ! csza_slp = 0
569 ! IF(csza_slp<=1.E-4) csza_slp = 0
570 !
571 ! Topographic shading
572 !
573 ! if (shadow==1) csza_slp = 0
574 !
575 ! Correction factor for sloping topography; the diffuse fraction of solar radiation is assumed to be unaffected by the slope
576 ! corr_fac = diffuse_frac + (1-diffuse_frac)*csza_slp/csza
577 !
578 ! 140 continue
579 !
580 ! GSW=(1.-ALBEDO)*SDOWN(kte+1)*corr_fac
581 !
582 ! endif
583 ! ENDIF
584 
585  !print *,"A1",TTEN(KE),SDOWN(KE),GSW,ALBEDO,cosSZA
586  !print *,"A2",QV(KE),QC(KE),QR(KE),QI(KE),QS(KE),QG(KE)
587  !print *,"A3",T(KE),P(KE),RHO0(KE),DZ(KE)
588  !print *,"A4",R,CP,G,solins
589  !print *,"A5",XXLAT,XXLON,XXLAT/D2R,XXLON/D2R
590  !print *,"A6",ICLOUD,aer_dry1(KE),aer_water1(KE)
591  !print *,"A7",kts,kte
592 
593  7 CONTINUE
594  return
595  END SUBROUTINE swpara
596 
597 !====================================================================
598  SUBROUTINE swinit
599 !--------------------------------------------------------------------
600  IMPLICIT NONE
601 !--------------------------------------------------------------------
602 ! LOGICAL , INTENT(IN) :: allowed_to_read
603 ! INTEGER , INTENT(IN) :: ids, ide, jds, jde, kds, kde, &
604 ! ims, ime, jms, jme, kms, kme, &
605 ! its, ite, jts, jte, kts, kte
606 ! real , INTENT(IN) :: swrad_scat
607 
608  LOGICAL :: allowed_to_read = .true.
609  real :: swrad_scat = 1
610 
611 
612 ! CSSCA - CLEAR-SKY SCATTERING SET FROM NAMELIST SWRAD_SCAT
613  cssca = swrad_scat * 1.e-5
614 
615  END SUBROUTINE swinit
616 
617 end module scale_atmos_phy_rd_mm5sw
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
real(rp), public atmos_solarins_constant
subroutine swpara(TTEN, SDOWN, GSW, ALBEDO, cosSZA, T, QV, QC, QR, QI, QS, QG, P, RHO0, DZ, R, CP, G, solins, XXLAT, XXLON, ICLOUD, aer_dry1, aer_water1, kts, kte)
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
integer, parameter, public i_lw
module grid index
integer, parameter, public i_sw
module TRACER
integer, public ia
of whole cells: x, local, with HALO
typedef float(real32_t)
integer, parameter, public i_dn
integer, public ka
of whole cells: z, local, with HALO
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
subroutine, public swrad(dt, RTHRATEN, SDOWN3D, GSW, XLAT, XLONG, ALBEDO, rho_phy, T3D, P3D, pi3D, dz8w, solins, cosSZA, QV3D, QC3D, QR3D, QI3D, QS3D, QG3D, F_QV, F_QC, F_QR, F_QI, F_QS, F_QG, icloud, warm_rain)
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
module PRECISION
module ATMOSPHERE / Physics Radiation
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:65
module ATMOSPHERE / Physics Radiation
integer, parameter, public i_up
integer, public ja
of whole cells: y, local, with HALO