SCALE-RM
Functions/Subroutines
scale_atmos_phy_rd_mm5sw Module Reference

module ATMOSPHERE / Physics Radiation More...

Functions/Subroutines

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)
 
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)
 
subroutine, public swinit
 

Detailed Description

module ATMOSPHERE / Physics Radiation

Description
Atmospheric radiation transfer process SWRAD: MM5 SW (Dudhia) scheme from WRF
Author
Team SCALE
History

Function/Subroutine Documentation

◆ swrad()

subroutine, public scale_atmos_phy_rd_mm5sw::swrad ( real(dp), intent(in)  dt,
real(rp), dimension(ia,ka,ja), intent(inout)  RTHRATEN,
real(rp), dimension(ia,ka,ja), intent(inout)  SDOWN3D,
real(rp), dimension(ia,ja), intent(inout)  GSW,
real(rp), dimension(ia,ja), intent(in)  XLAT,
real(rp), dimension(ia,ja), intent(in)  XLONG,
real(rp), dimension(ia, ja), intent(in)  ALBEDO,
real(rp), dimension(ia,ka,ja), intent(in)  rho_phy,
real(rp), dimension(ia,ka,ja), intent(in)  T3D,
real(rp), dimension(ia,ka,ja), intent(in)  P3D,
real(rp), dimension(ia,ka,ja), intent(in)  pi3D,
real(rp), dimension(ia,ka,ja), intent(in)  dz8w,
real(rp), dimension(ia,ja), intent(in)  solins,
real(rp), dimension(ia,ja), intent(in)  cosSZA,
real(rp), dimension (ia,ka,ja), intent(in), optional  QV3D,
real(rp), dimension (ia,ka,ja), intent(in), optional  QC3D,
real(rp), dimension (ia,ka,ja), intent(in), optional  QR3D,
real(rp), dimension (ia,ka,ja), intent(in), optional  QI3D,
real(rp), dimension (ia,ka,ja), intent(in), optional  QS3D,
  QG3D,
real(rp), intent(in), optional  F_QV,
real(rp), intent(in), optional  F_QC,
real(rp), intent(in), optional  F_QR,
real(rp), intent(in), optional  F_QI,
real(rp), intent(in), optional  F_QS,
real(rp), intent(in), optional  F_QG,
integer, intent(in)  icloud,
logical, intent(in)  warm_rain 
)

Definition at line 58 of file scale_atmos_phy_rd_mm5sw.F90.

58 
59  use scale_const, only: &
60  grav => const_grav, & ! 9.81
61  rdry => const_rdry, & ! 287.
62  cpdry => const_cpdry ! 1004.6
63  use scale_atmos_solarins, only: &
64  solarins => atmos_solarins_constant ! 1360.xxx (W/m2)
65 !------------------------------------------------------------------
66  IMPLICIT NONE
67 !------------------------------------------------------------------
68 
69  real(DP), INTENT(IN) :: dt
70  real(RP), INTENT(INOUT) :: RTHRATEN(IA,KA,JA)
71  real(RP), INTENT(INOUT) :: SDOWN3D(IA,KA,JA)
72  real(RP), INTENT(INOUT) :: GSW(IA,JA)
73 
74  real(RP), INTENT(IN) :: XLAT(IA,JA), XLONG(IA,JA), ALBEDO(IA, JA)
75  real(RP), INTENT(IN) :: rho_phy(IA,KA,JA)
76  real(RP), INTENT(IN) :: P3D(IA,KA,JA), &
77  T3D(IA,KA,JA), &
78  pi3D(IA,KA,JA), &
79  dz8w(IA,KA,JA)
80  real(RP), INTENT(IN) :: solins(IA,JA),cosSZA(IA,JA)
81  real(RP), OPTIONAL, INTENT(IN) :: QV3D (IA,KA,JA), &
82  QC3D (IA,KA,JA), &
83  QR3D (IA,KA,JA), &
84  QI3D (IA,KA,JA), &
85  QS3D (IA,KA,JA), &
86  QG3D (IA,KA,JA)
87 
88  LOGICAL, OPTIONAL, INTENT(IN ) :: F_QV,F_QC,F_QR,F_QI,F_QS,F_QG
89  INTEGER, INTENT(IN ) :: icloud
90  LOGICAL, INTENT(IN ) :: warm_rain
91 
92  !real, INTENT(IN ) :: RADFRQ,DEGRAD,XTIME,DECLIN
93 
94 
95  integer :: its,ite,jts,jte,kts,kte
96  real(RP) :: R, CP, G, SOLCON
97 
98 ! LOCAL VARS
99 
100  real(RP), DIMENSION( KS:KE ) :: TTEN1D, &
101  RHO01D, &
102  P1D, &
103  DZ, &
104  T1D, &
105  QV1D, &
106  QC1D, &
107  QR1D, &
108  QI1D, &
109  QS1D, &
110  QG1D
111 
112  real(RP) :: SDOWN1D(KS:KE+1)
113 !
114  real(RP) :: XLAT0,XLONG0,ALB0,GSW0,cosSZA0,solins0
115  real(RP) :: aer_dry1(KS:KE),aer_water1(KS:KE)
116 
117 !
118  INTEGER :: i,j,K,NK
119  LOGICAL :: predicate , do_topo_shading
120 
121 
122 !------------------------------------------------------------------
123 
124  r = rdry
125  cp = cpdry
126  g = grav
127  solcon = solarins
128  its = is ; ite = ie
129  jts = js ; jte = je
130  kts = ks ; kte = ke
131 
132 
133  j_loop: DO j=jts,jte
134  i_loop: DO i=its,ite
135 
136 ! reverse vars
137  DO k=ks,ke
138  qv1d(k)=0.
139  qc1d(k)=0.
140  qr1d(k)=0.
141  qi1d(k)=0.
142  qs1d(k)=0.
143  qg1d(k)=0.
144  ENDDO
145 
146  DO k=ks,ke
147  !! NK=kme-1-K+kms
148  nk=ke-(k-ks)
149  tten1d(k)=0.
150 
151  t1d(k)=t3d(i,nk,j)
152  p1d(k)=p3d(i,nk,j)
153  rho01d(k)=rho_phy(i,nk,j)
154  dz(k)=dz8w(i,nk,j)
155  ENDDO
156 
157  !IF( PRESENT(pm2_5_dry) .AND. PRESENT(pm2_5_water) )THEN
158  ! DO K=kts,kte
159  ! NK=kme-1-K+kms
160  ! aer_dry1(k) = pm2_5_dry(i,nk,j)
161  ! aer_water1(k) = pm2_5_water(i,nk,j)
162  ! ENDDO
163  !ELSE
164  do k=ks,ke
165  aer_dry1(k) = 0.0
166  aer_water1(k) = 0.0
167  enddo
168  !ENDIF
169 
170  IF (PRESENT(f_qv) .AND. PRESENT(qv3d)) THEN
171  IF (f_qv) THEN
172  do k=ks,ke
173  !NK=kme-1-K+kms
174  nk=ke-(k-ks)
175  qv1d(k)=qv3d(i,nk,j)
176  qv1d(k)=max(0.0_rp,qv1d(k))
177  ENDDO
178  ENDIF
179  ENDIF
180 
181  IF (PRESENT(f_qc) .AND. PRESENT(qc3d)) THEN
182  IF (f_qc) THEN
183  do k=ks,ke
184  !NK=kme-1-K+kms
185  nk=ke-(k-ks)
186  qc1d(k)=qc3d(i,nk,j)
187  qc1d(k)=max(0.0_rp,qc1d(k))
188  enddo
189  ENDIF
190  ENDIF
191 
192  IF (PRESENT(f_qr) .AND. PRESENT(qr3d)) THEN
193  IF (f_qr) THEN
194  do k=kts,kte
195  !NK=kme-1-K+kms
196  nk=ke-(k-ks)
197  qr1d(k)=qr3d(i,nk,j)
198  qr1d(k)=max(0.0_rp,qr1d(k))
199  enddo
200  ENDIF
201  ENDIF
202 
203 !
204  IF ( PRESENT( f_qi ) ) THEN
205  predicate = f_qi
206  ELSE
207  predicate = .false.
208  ENDIF
209 
210  IF ( predicate .AND. PRESENT( qi3d ) ) THEN
211  do k=ks,ke
212  !NK=kme-1-K+kms
213  nk=ke-(k-ks)
214  qi1d(k)=qi3d(i,nk,j)
215  qi1d(k)=max(0.0_rp,qi1d(k))
216  enddo
217  ELSE
218  IF (.NOT. warm_rain) THEN
219  do k=ks,ke
220  IF(t1d(k) < 273.15) THEN
221  qi1d(k)=qc1d(k)
222  qc1d(k)=0.
223  qs1d(k)=qr1d(k)
224  qr1d(k)=0.
225  ENDIF
226  enddo
227  ENDIF
228  ENDIF
229 
230  IF (PRESENT(f_qs) .AND. PRESENT(qs3d)) THEN
231  IF (f_qs) THEN
232  do k=ks,ke
233  !NK=kme-1-K+kms
234  nk=ke-(k-ks)
235  qs1d(k)=qs3d(i,nk,j)
236  qs1d(k)=max(0.0_rp,qs1d(k))
237  enddo
238  ENDIF
239  ENDIF
240 
241  IF (PRESENT(f_qg) .AND. PRESENT(qg3d)) THEN
242  IF (f_qg) THEN
243  do k=ks,ke
244  !NK=kme-1-K+kms
245  nk=ke-(k-ks)
246  qg1d(k)=qg3d(i,nk,j)
247  qg1d(k)=max(0.0_rp,qg1d(k))
248  enddo
249  ENDIF
250  ENDIF
251 
252  xlat0=xlat(i,j)
253  xlong0=xlong(i,j)
254  solins0=solins(i,j)
255  cossza0=cossza(i,j)
256  alb0=albedo(i,j)
257 ! slope code removed - factor now done in surface driver
258  CALL swpara(tten1d,sdown1d,gsw0,alb0,cossza0, &
259  t1d,qv1d,qc1d,qr1d,qi1d,qs1d,qg1d,p1d, &
260  rho01d,dz, &
261  r,cp,g,solins0, &
262  xlat0,xlong0, &
263  icloud,aer_dry1,aer_water1, &
264  kts,kte )
265  gsw(i,j)=gsw0
266 
267 ! DO K=kts,kte
268 ! NK=kme-1-K+kms
269 ! NK=KA-1-K+1
270 ! RTHRATEN(I,K,J)=RTHRATEN(I,K,J)+TTEN1D(NK)/pi3D(I,K,J)
271 ! SDOWN3D(I,K,J)=SDOWN1D(NK)
272 ! if(k==kte)then
273 ! NK=KA-1-(K+1)+1
274 ! SDOWN3D(I,K+1,J)=SDOWN1D(NK)
275 ! endif
276 ! ENDDO
277 
278  do k=ks,ke
279  nk=ke-(k-ks)
280  rthraten(i,k,j)=rthraten(i,k,j)+tten1d(nk)/pi3d(i,k,j)
281  sdown3d(i,k,j)=sdown1d(nk)
282  !print *,KS,KE,K,"<=",NK
283  enddo
284  k=ks-1 ; nk=ke+1
285  sdown3d(i,ks-1,j)=sdown1d(ke+1)
286  !print *,KS,KE,K,"<=",NK
287 !
288  ENDDO i_loop
289  ENDDO j_loop
290 
291  return

References scale_atmos_solarins::atmos_solarins_constant, scale_const::const_grav, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_tracer::k, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, and swpara().

Here is the call graph for this function:

◆ swpara()

subroutine scale_atmos_phy_rd_mm5sw::swpara ( real(rp), dimension( kts:kte ), intent(inout)  TTEN,
real(rp), dimension( kts:kte+1 ), intent(out)  SDOWN,
real(rp), intent(inout)  GSW,
real(rp), intent(in)  ALBEDO,
real(rp), intent(in)  cosSZA,
real(rp), dimension( kts:kte ), intent(in)  T,
real(rp), dimension( kts:kte ), intent(in)  QV,
real(rp), dimension( kts:kte ), intent(in)  QC,
real(rp), dimension( kts:kte ), intent(in)  QR,
real(rp), dimension( kts:kte ), intent(in)  QI,
real(rp), dimension( kts:kte ), intent(in)  QS,
real(rp), dimension( kts:kte ), intent(in)  QG,
real(rp), dimension( kts:kte ), intent(in)  P,
real(rp), dimension( kts:kte ), intent(in)  RHO0,
real(rp), dimension( kts:kte ), intent(in)  DZ,
real(rp), intent(in)  R,
real(rp), intent(in)  CP,
real(rp), intent(in)  G,
real(rp), intent(in)  solins,
real(rp), intent(in)  XXLAT,
real(rp), intent(in)  XXLON,
integer, intent(in)  ICLOUD,
real(rp), dimension( kts:kte )  aer_dry1,
real(rp), dimension( kts:kte )  aer_water1,
integer, intent(in)  kts,
integer, intent(in)  kte 
)

Definition at line 302 of file scale_atmos_phy_rd_mm5sw.F90.

302 
303  use scale_time, only: &
304  nowdate => time_nowdate !< current time [yyyy mm dd hh mm ss]
305  use scale_const, only: &
306  d2r => const_d2r ! degree to radian
307 !------------------------------------------------------------------
308 ! TO CALCULATE SHORT-WAVE ABSORPTION AND SCATTERING IN CLEAR
309 ! AIR AND REFLECTION AND ABSORPTION IN CLOUD LAYERS (STEPHENS,
310 ! 1984)
311 ! CHANGES:
312 ! REDUCE EFFECTS OF ICE CLOUDS AND PRECIP ON LIQUID WATER PATH
313 ! ADD EFFECT OF GRAUPEL
314 !------------------------------------------------------------------
315  IMPLICIT NONE
316 
317  INTEGER, INTENT(IN ) :: kts,kte
318 
319  real(RP), DIMENSION( kts:kte ), INTENT(IN ) :: &
320  RHO0, &
321  T, &
322  P, &
323  DZ, &
324  QV, &
325  QC, &
326  QR, &
327  QI, &
328  QS, &
329  QG
330 
331  real(RP), DIMENSION( kts:kte ), INTENT(INOUT) :: TTEN
332  real(RP), DIMENSION( kts:kte+1 ), INTENT(OUT) :: SDOWN
333 !
334  real(RP), INTENT(IN ) :: R,CP,G
335  real(RP), INTENT(IN) :: ALBEDO,cosSZA,solins
336  integer , INTENT(IN) :: icloud
337  real(RP), INTENT(INOUT) :: GSW
338  real(RP), INTENT(IN ) :: XXLAT,XXLON
339 
340  !real, INTENT(IN ) :: GMT,RADFRQ,XTIME,XLAT,XLONG,DEGRAD
341 !
342 ! For slope-dependent radiation
343 ! INTEGER, OPTIONAL, INTENT(IN) :: slope_rad,shadow
344 ! real, OPTIONAL, INTENT(IN) :: slp_azi,slope
345 !
346 
347 ! LOCAL VARS
348 
349  real(RP) :: XLWP( kts:kte )
350  real(RP) :: XATP( kts:kte )
351  real(RP) :: XWVP( kts:kte )
352  real(RP) :: aer_dry1( kts:kte )
353  real(RP) :: aer_water1( kts:kte )
354  real(RP) :: RO( kts:kte )
355 !
356  real(RP) :: ALBTAB(4,5), ABSTAB(4,5)
357  real(RP) :: XMUVAL(4)
358 
359  real(RP) :: beta
360  real(RP) :: DayOfYear
361 !------------------------------------------------------------------
362 
363  DATA albtab/0.0,0.0,0.0,0.0, &
364  69.0,58.0,40.0,15.0, &
365  90.0,80.0,70.0,60.0, &
366  94.0,90.0,82.0,78.0, &
367  96.0,92.0,85.0,80.0/
368 
369  DATA abstab/0.0,0.0,0.0,0.0, &
370  0.0,2.5,4.0,5.0, &
371  0.0,2.6,7.0,10.0, &
372  0.0,3.3,10.0,14.0, &
373  0.0,3.7,10.0,15.0/
374 
375  DATA xmuval/0.0,0.2,0.5,1.0/
376 
377  real(RP) :: bext340, absc, alba, alw, csza,dabsa,dsca,dabs
378  real(RP) :: bexth2o, dscld, ff,oldalb,oldabs,oldabc
379  real(RP) :: soltop, totabs, ugcm, uv,xabs,xabsa,wv
380  real(RP) :: wgm, xalb, xi, xsca, xmu,xabsc,trans0,yj
381  real(RP) :: HRANG,XT24,TLOCTM,tloc,dsec,lon,lat
382  real(RP) :: DECLIN
383  real(RP) :: ww
384  INTEGER :: iil,ii,jjl,ju,k,iu
385 
386 ! For slope-dependent radiation
387 
388  real(RP) :: diffuse_frac, corr_fac, csza_slp
389 
390  gsw=0.0
391  sdown=0.0
392  bext340=5.e-6_rp
393  bexth2o=5.e-6_rp
394  soltop = solins
395 
396  !XT24=MOD(XTIME+RADFRQ*0.5,1440.)
397  !TLOCTM=GMT+XT24/60.+XLONG/15.
398  !HRANG=15.*(TLOCTM-12.)*DEGRAD
399  !XXLAT=XLAT*DEGRAD
400  !CSZA=SIN(XXLAT)*SIN(DECLIN)+COS(XXLAT)*COS(DECLIN)*COS(HRANG)
401 
402  ! local time
403  lat = xxlat / d2r
404  lon = xxlon / d2r
405 
406  tloc = mod((nowdate(4) + int(lon/15.0_rp)),24 )
407  dsec = real(nowdate(5)*60 + nowdate(6)) / 60.0_rp /60.0_rp
408  tloctm = real(nowdate(4)) + lon/15.0_rp + dsec
409  hrang = 15.*(tloctm-12.)*d2r
410 
411  csza=cossza
412 
413 ! RETURN IF NIGHT
414  IF(csza <= 1.0e-9_rp)GOTO 7
415 !
416  DO k=kts, kte
417 
418 ! P in the unit of 10mb
419  ro(k)=p(k)/(r*t(k))
420  xwvp(k)=ro(k)*qv(k)*dz(k)*1000.
421 ! KG/M**2
422  xatp(k)=ro(k)*dz(k)
423  ENDDO
424 !
425 ! G/M**2
426 ! REDUCE WEIGHT OF LIQUID AND ICE IN SHORT-WAVE SCHEME
427 ! ADD GRAUPEL EFFECT (ASSUMED SAME AS RAIN)
428 !
429  IF (icloud==0)THEN
430  DO k=kts, kte
431  xlwp(k)=0.
432  ENDDO
433  ELSE
434  DO k=kts, kte
435  xlwp(k)=ro(k)*1000.*dz(k)*(qc(k)+0.1*qi(k)+0.05* &
436  qr(k)+0.02*qs(k)+0.05*qg(k))
437  ENDDO
438  ENDIF
439 !
440  xmu=csza
441  ! SDOWN(1)=SOLTOP*XMU !adachi
442 
443  sdown(kts)=soltop * max(cossza,0.0_rp)
444 
445 ! SET WW (G/M**2) LIQUID WATER PATH INTEGRATED DOWN
446 ! SET UV (G/M**2) WATER VAPOR PATH INTEGRATED DOWN
447  ww=0.
448  uv=0.
449  oldalb=0.
450  oldabc=0.
451  totabs=0.
452 ! CONTRIBUTIONS DUE TO CLEAR AIR AND CLOUD
453  dsca=0.
454  dabs=0.
455  dscld=0.
456 !
457 ! CONTRIBUTION DUE TO AEROSOLS (FOR CHEMISTRY)
458  dabsa=0.
459 !
460  DO 200 k=kts,kte
461  ww=ww+xlwp(k)
462  uv=uv+xwvp(k)
463 ! WGM IS WW/COS(THETA) (G/M**2)
464 ! UGCM IS UV/COS(THETA) (G/CM**2)
465  wgm=ww/xmu
466  ugcm=uv*0.0001/xmu
467 !
468  oldabs=totabs
469 ! WATER VAPOR ABSORPTION AS IN LACIS AND HANSEN (1974)
470  totabs=2.9*ugcm/((1.+141.5*ugcm)**0.635+5.925*ugcm)
471 ! APPROXIMATE RAYLEIGH + AEROSOL SCATTERING
472 ! XSCA=1.E-5*XATP(K)/XMU
473 ! XSCA=(1.E-5*XATP(K)+aer_dry1(K)*bext340+aer_water1(K)*bexth2o)/XMU
474  beta=0.4*(1.0-xmu)+0.1
475 ! CSSCA - CLEAR-SKY SCATTERING SET FROM NAMELIST SWRAD_SCAT
476  xsca=(cssca*xatp(k)+beta*aer_dry1(k)*bext340*dz(k) &
477  +beta*aer_water1(k)*bexth2o*dz(k))/xmu
478 
479 ! LAYER VAPOR ABSORPTION DONE FIRST
480  ! XABS=(TOTABS-OLDABS)*(SDOWN(1)-DSCLD-DSCA-DABSA)/SDOWN(K)
481  xabs=(totabs-oldabs)*(sdown(kts)-dscld-dsca-dabsa)/sdown(k)
482 !rs AEROSOL ABSORB (would be elemental carbon). So far XABSA = 0.
483  xabsa=0.
484  IF(xabs<0.)xabs=0.
485 !
486  alw=log10(wgm+1.0_rp)
487  IF(alw>3.999)alw=3.999_rp
488 !
489  DO ii=1,3
490  IF(xmu>xmuval(ii))THEN
491  iil=ii
492  iu=ii+1
493  xi=(xmu-xmuval(ii))/(xmuval(ii+1)-xmuval(ii))+float(iil)
494  ENDIF
495  ENDDO
496 !
497  jjl = int(alw)+1
498  ju = jjl+1
499  yj = alw+1.
500 ! CLOUD ALBEDO
501  alba=(albtab(iu,ju)*(xi-iil)*(yj-jjl) &
502  +albtab(iil,ju)*(iu-xi)*(yj-jjl) &
503  +albtab(iu,jjl)*(xi-iil)*(ju-yj) &
504  +albtab(iil,jjl)*(iu-xi)*(ju-yj)) &
505  /((iu-iil)*(ju-jjl))
506 ! CLOUD ABSORPTION
507  absc=(abstab(iu,ju)*(xi-iil)*(yj-jjl) &
508  +abstab(iil,ju)*(iu-xi)*(yj-jjl) &
509  +abstab(iu,jjl)*(xi-iil)*(ju-yj) &
510  +abstab(iil,jjl)*(iu-xi)*(ju-yj)) &
511  /((iu-iil)*(ju-jjl))
512 ! LAYER ALBEDO AND ABSORPTION
513  !XALB=(ALBA-OLDALB)*(SDOWN(1)-DSCA-DABS)/SDOWN(K)
514  !XABSC=(ABSC-OLDABC)*(SDOWN(1)-DSCA-DABS)/SDOWN(K)
515  xalb=(alba-oldalb)*(sdown(kts)-dsca-dabs)/sdown(k)
516  xabsc=(absc-oldabc)*(sdown(kts)-dsca-dabs)/sdown(k)
517  IF(xalb<0.)xalb=0.
518  IF(xabsc<0.)xabsc=0.
519  dscld=dscld+(xalb+xabsc)*sdown(k)*0.01
520  dsca=dsca+xsca*sdown(k)
521  dabs=dabs+xabs*sdown(k)
522  dabsa=dabsa+xabsa*sdown(k)
523  oldalb=alba
524  oldabc=absc
525 ! LAYER TRANSMISSIVITY
526  trans0=100.0_rp-xalb-xabsc-xabs*100.0_rp-xsca*100.0_rp
527  IF(trans0<1.)THEN
528  ff=99.0_rp/(xalb+xabsc+xabs*100.0_rp+xsca*100.0_rp)
529  xalb=xalb*ff
530  xabsc=xabsc*ff
531  xabs=xabs*ff
532  xsca=xsca*ff
533  trans0=1.0_rp
534  ENDIF
535  sdown(k+1)=max(1.0e-9_rp,sdown(k)*trans0*0.01_rp)
536  tten(k)=sdown(k)*(xabsc+xabs*100._rp+xabsa*100._rp)*0.01_rp/(ro(k)*cp*dz(k))
537  200 CONTINUE
538 !
539  gsw=(1.-albedo)*sdown(kte+1)
540 
541 ! IF (PRESENT(slope_rad)) THEN
542 !! Slope-dependent solar radiation part
543 !
544 ! if (slope_rad==1) then
545 !
546 !! Parameterize diffuse fraction of global solar radiation as a function of the ratio between TOA radiation and surface global radiation
547 !
548 ! diffuse_frac = min(1.,1/(max(0.1,2.1-2.8*log(log(SDOWN(kts)/max(SDOWN(kte+1),1.e-3))))))
549 ! 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
550 ! corr_fac = 1
551 ! goto 140
552 ! endif
553 !
554 !! cosine of zenith angle over sloping topography
555 !
556 ! csza_slp = ((SIN(XXLAT)*COS(HRANG))* &
557 ! (-cos(slp_azi)*sin(slope))-SIN(HRANG)*(sin(slp_azi)*sin(slope))+ &
558 ! (COS(XXLAT)*COS(HRANG))*cos(slope))* &
559 ! COS(DECLIN)+(COS(XXLAT)*(cos(slp_azi)*sin(slope))+ &
560 ! SIN(XXLAT)*cos(slope))*SIN(DECLIN)
561 !
562 ! csza_slp = 0
563 ! IF(csza_slp<=1.E-4) csza_slp = 0
564 !
565 ! Topographic shading
566 !
567 ! if (shadow==1) csza_slp = 0
568 !
569 ! Correction factor for sloping topography; the diffuse fraction of solar radiation is assumed to be unaffected by the slope
570 ! corr_fac = diffuse_frac + (1-diffuse_frac)*csza_slp/csza
571 !
572 ! 140 continue
573 !
574 ! GSW=(1.-ALBEDO)*SDOWN(kte+1)*corr_fac
575 !
576 ! endif
577 ! ENDIF
578 
579  !print *,"A1",TTEN(KE),SDOWN(KE),GSW,ALBEDO,cosSZA
580  !print *,"A2",QV(KE),QC(KE),QR(KE),QI(KE),QS(KE),QG(KE)
581  !print *,"A3",T(KE),P(KE),RHO0(KE),DZ(KE)
582  !print *,"A4",R,CP,G,solins
583  !print *,"A5",XXLAT,XXLON,XXLAT/D2R,XXLON/D2R
584  !print *,"A6",ICLOUD,aer_dry1(KE),aer_water1(KE)
585  !print *,"A7",kts,kte
586 
587  7 CONTINUE
588  return

References scale_const::const_d2r, float(), and scale_time::time_nowdate.

Referenced by swrad().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ swinit()

subroutine, public scale_atmos_phy_rd_mm5sw::swinit

Definition at line 593 of file scale_atmos_phy_rd_mm5sw.F90.

593 !--------------------------------------------------------------------
594  IMPLICIT NONE
595 !--------------------------------------------------------------------
596 ! LOGICAL , INTENT(IN) :: allowed_to_read
597 ! INTEGER , INTENT(IN) :: ids, ide, jds, jde, kds, kde, &
598 ! ims, ime, jms, jme, kms, kme, &
599 ! its, ite, jts, jte, kts, kte
600 ! real , INTENT(IN) :: swrad_scat
601 
602  LOGICAL :: allowed_to_read = .true.
603  real :: swrad_scat = 1
604 
605 
606 ! CSSCA - CLEAR-SKY SCATTERING SET FROM NAMELIST SWRAD_SCAT
607  cssca = swrad_scat * 1.e-5
608 
scale_const::const_grav
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:49
scale_atmos_solarins::atmos_solarins_constant
real(rp), public atmos_solarins_constant
Definition: scale_atmos_solarins.F90:43
scale_atmos_solarins
module atmosphere / SOLARINS
Definition: scale_atmos_solarins.F90:14
float
typedef float(real32_t)
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_time
module TIME
Definition: scale_time.F90:11
scale_time::time_nowdate
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:68
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:33