Radiation Scheme in CCPP
grrad.f
Go to the documentation of this file.
1 
7 ! ========================================================== !!!!!
8 ! 'module_radiation_driver' descriptions !!!!!
9 ! ========================================================== !!!!!
10 ! !
11 ! this is the radiation driver module. it prepares atmospheric !
12 ! profiles and invokes main radiation calculations. !
13 ! !
14 ! in module 'module_radiation_driver' there are twe externally !
15 ! callable subroutine: !
16 ! !
17 ! 'radinit' -- initialization routine !
18 ! input: !
19 ! ( si, NLAY, me ) !
20 ! output: !
21 ! ( none ) !
22 ! !
23 ! 'radupdate' -- update time sensitive data used by radiations !
24 ! input: !
25 ! ( idate,jdate,deltsw,deltim,lsswr, me ) !
26 ! output: !
27 ! ( slag,sdec,cdec,solcon ) !
28 ! !
29 ! 'grrad' -- setup and invoke main radiation calls !
30 ! input: !
31 ! ( prsi,prsl,prslk,tgrs,qgrs,tracer,vvl,slmsk, !
32 ! xlon,xlat,tsfc,snowd,sncovr,snoalb,zorl,hprim, !
33 ! alvsf,alnsf,alvwf,alnwf,facsf,facwf,fice,tisfc, !
34 ! sinlat,coslat,solhr,jdate,solcon, !
35 ! cv,cvt,cvb,fcice,frain,rrime,flgmin, !
36 ! icsdsw,icsdlw, ntcw,ncld,ntoz, NTRAC,NFXR, !
37 ! dtlw,dtsw, lsswr,lslwr,lssav, !
38 ! IX, IM, LM, me, lprnt, ipt, kdt,deltaq,sup,cnvw,cnvc, !
39 ! output: !
40 ! htrsw,topfsw,sfcfsw,dswcmp,uswcmp,sfalb,coszen,coszdg, !
41 ! htrlw,topflw,sfcflw,tsflw,semis,cldcov, !
42 ! input/output: !
43 ! fluxr !
44 ! optional output: !
45 ! htrlw0,htrsw0,htrswb,htrlwb !
46 ! !
47 ! !
48 ! external modules referenced: !
49 ! 'module physparam' in 'physparam.f' !
50 ! 'module funcphys' in 'funcphys.f' !
51 ! 'module physcons' in 'physcons.f' !
52 ! !
53 ! 'module module_radiation_gases' in 'radiation_gases.f' !
54 ! 'module module_radiation_aerosols' in 'radiation_aerosols.f' !
55 ! 'module module_radiation_surface' in 'radiation_surface.f' !
56 ! 'module module_radiation_clouds' in 'radiation_clouds.f' !
57 ! !
58 ! 'module module_radsw_cntr_para' in 'radsw_xxxx_param.f' !
59 ! 'module module_radsw_parameters' in 'radsw_xxxx_param.f' !
60 ! 'module module_radsw_main' in 'radsw_xxxx_main.f' !
61 ! !
62 ! 'module module_radlw_cntr_para' in 'radlw_xxxx_param.f' !
63 ! 'module module_radlw_parameters' in 'radlw_xxxx_param.f' !
64 ! 'module module_radlw_main' in 'radlw_xxxx_main.f' !
65 ! !
66 ! where xxxx may vary according to different scheme selection !
67 ! !
68 ! !
69 ! program history log: !
70 ! mm-dd-yy ncep - created program grrad !
71 ! 08-12-03 yu-tai hou - re-written for modulized radiations !
72 ! 11-06-03 yu-tai hou - modified !
73 ! 01-18-05 s. moorthi - NOAH/ICE model changes added !
74 ! 05-10-05 yu-tai hou - modified module structure !
75 ! 12-xx-05 s. moorthi - sfc lw flux adj by mean temperature !
76 ! 02-20-06 yu-tai hou - add time variation for co2 data, and !
77 ! solar const. add sfc emiss change !
78 ! 03-21-06 s. Moorthi - added surface temp over ice !
79 ! 07-28-06 yu-tai hou - add stratospheric vocanic aerosols !
80 ! 03-14-07 yu-tai hou - add generalized spectral band interp !
81 ! for aerosol optical prop. (sw and lw) !
82 ! 04-10-07 yu-tai hou - spectral band sw/lw heating rates !
83 ! 05-04-07 yu-tai hou - make options for clim based and modis !
84 ! based (h. wei and c. marshall) albedo !
85 ! 09-05-08 yu-tai hou - add the initial date and time 'idate' !
86 ! and control param 'ICTM' to the passing param list!
87 ! to handel different time/date requirements for !
88 ! external data (co2, aeros, solcon, ...) !
89 ! 10-10-08 yu-tai hou - add the ICTM=-2 option for combining !
90 ! initial condition data with seasonal cycle from !
91 ! climatology. !
92 ! 03-12-09 yu-tai hou - use two time stamps to keep tracking !
93 ! dates for init cond and fcst time. remove volcanic!
94 ! aerosols data in climate hindcast (ICTM=-2). !
95 ! 03-16-09 yu-tai hou - included sub-column clouds approx. !
96 ! control flags isubcsw/isubclw in initialization !
97 ! subroutine. passed auxiliary cloud control arrays !
98 ! icsdsw/icsdlw (if isubcsw/isubclw =2, it will be !
99 ! the user provided permutation seeds) to the sw/lw !
100 ! radiation calculation programs. also moved cloud !
101 ! overlapping control flags iovrsw/iovrlw from main !
102 ! radiation routines to the initialization routines.!
103 ! 04-02-09 yu-tai hou - modified surface control flag iems to !
104 ! have additional function of if the surface-air !
105 ! interface have the same or different temperature !
106 ! for radiation calculations. !
107 ! 04-03-09 yu-tai hou - modified to add lw surface emissivity !
108 ! as output variable. changed the sign of sfcnsw to !
109 ! be positive value denote solar flux goes into the !
110 ! ground (this is needed to reduce sign confusion !
111 ! in other part of model) !
112 ! 09-09-09 fanglin yang (thru s.moorthi) added QME5 QME6 to E-20!
113 ! 01-09-10 sarah lu - added gocart option, revised grrad for!
114 ! gocart coupling. calling argument modifed: ldiag3 !
115 ! removed; cldcov/fluxr sequence changed; cldcov is !
116 ! changed from accumulative to instant field and !
117 ! from input/output to output field !
118 ! 01-24-10 sarah lu - added aod to fluxr, added prslk and !
119 ! oz to setaer input argument (for gocart coupling),!
120 ! added tau_gocart to setaer output argument (for, !
121 ! aerosol diag by index of nv_aod) !
122 ! 07-08-10 s.moorthi - updated the NEMS version for new physics !
123 ! 07-28-10 yu-tai hou - changed grrad interface to allow all !
124 ! components of sw/lw toa/sfc instantaneous values !
125 ! being passed to the calling program. moved the !
126 ! computaion of sfc net sw flux (sfcnsw) to the !
127 ! calling program. merged carlos' nmmb modification.!
128 ! 07-30-10 s. moorthi - corrected some errors associated with !
129 ! unit changes !
130 ! 12-02-10 s. moorthi/y. hou - removed the use of aerosol flags !
131 ! 'iaersw' 'iaerlw' from radiations and replaced !
132 ! them by using the runtime variable laswflg and !
133 ! lalwflg defined in module radiation_aerosols. !
134 ! also replaced param nspc in grrad with the use of !
135 ! max_num_gridcomp in module radiation_aerosols. !
136 ! jun 2012 yu-tai hou - added sea/land madk 'slmsk' to the !
137 ! argument list of subrotine setaer call for the !
138 ! newly modified horizontal bi-linear interpolation !
139 ! in climatological aerosols schem. also moved the !
140 ! virtual temperature calculations in subroutines !
141 ! 'radiation_clouds' and 'radiation_aerosols' to !
142 ! 'grrad' to reduce repeat comps. renamed var oz as !
143 ! tracer to reflect that it carries various prog !
144 ! tracer quantities. !
145 ! - modified to add 4 compontents of sw !
146 ! surface downward fluxes to the output. (vis/nir; !
147 ! direct/diffused). re-arranged part of the fluxr !
148 ! variable fields and filled the unused slots for !
149 ! the new components. added check print of select !
150 ! data (co2 value for now). !
151 ! - changed the initialization subrution !
152 ! 'radinit' into two parts: 'radinit' is called at !
153 ! the start of model run to set up radiation related!
154 ! fixed parameters; and 'radupdate' is called in !
155 ! the time-loop to update time-varying data sets !
156 ! and module variables. !
157 ! sep 2012 h-m lin/y-t hou added option of extra top layer for !
158 ! models with low toa ceiling. the extra layer will !
159 ! help ozone absorption at higher altitude. !
160 ! nov 2012 yu-tai hou - modified control parameters through !
161 ! module 'physparam'. !
162 ! jan 2013 yu-tai hou - updated subr radupdate for including !
163 ! options of annual/monthly solar constant table. !
164 ! mar 2013 h-m lin/y-t hou corrected a bug in extra top layer !
165 ! when using ferrier microphysics. !
166 ! may 2013 s. mooorthi - removed fpkapx !
167 ! jul 2013 r. sun - added pdf cld and convective cloud water and!
168 ! cover for radiation !
169 ! aug 2013 s. moorthi - port from gfs to nems !
170 ! 13Feb2014 sarah lu - add aerodp to fluxr !
171 ! Apr 2014 Xingren Wu - add sfc SW downward fluxes nir/vis and !
172 ! sfcalb to export for A/O/I coupling !
173 ! jun 2014 y-t hou - revised code to include surface up and !
174 ! down spectral components sw fluxes as output. !
175 ! !
176 !!!!! ========================================================== !!!!!
177 !!!!! end descriptions !!!!!
178 !!!!! ========================================================== !!!!!
179 
180 
181 
182 !========================================!
183 
185 !........................................!
186 !
187  use physparam
188  use physcons, only : eps => con_eps, &
189  & epsm1 => con_epsm1, &
190  & fvirt => con_fvirt &
191  &, rocp => con_rocp
192  use funcphys, only : fpvs
193 
198  & aer_init, aer_update, &
199  & nspc1
201  & setemis
204  & diagcld1
205 
208  use module_radsw_main, only : rswinit, swrad
209 
211  & proflw_type, nbdlw
212  use module_radlw_main, only : rlwinit, lwrad
213 !
214  implicit none
215 !
216  private
217 
218 ! --- version tag and last revision date
219  character(40), parameter :: &
220  & VTAGRAD='NCEP-Radiation_driver v5.2 Jan 2013 '
221 ! & VTAGRAD='NCEP-Radiation_driver v5.1 Nov 2012 '
222 ! & VTAGRAD='NCEP-Radiation_driver v5.0 Aug 2012 '
223 
224 ! --- constant values
225  real (kind=kind_phys) :: qmin, qme5, qme6, epsq
226 ! parameter (QMIN=1.0e-10, QME5=1.0e-5, QME6=1.0e-6, EPSQ=1.0e-12)
227  parameter(qmin=1.0e-10, qme5=1.0e-7, qme6=1.0e-7, epsq=1.0e-12)
228 ! parameter (QMIN=1.0e-10, QME5=1.0e-20, QME6=1.0e-20, EPSQ=1.0e-12)
229  real, parameter :: prsmin = 1.0e-6 ! toa pressure minimum value in mb (hpa)
230 
231 ! --- control flags set in subr radinit:
232  integer :: itsfc =0 ! flag for lw sfc air/ground interface temp setting
233 
234 ! --- data input control variables set in subr radupdate:
235  integer :: month0=0, iyear0=0, monthd=0
236  logical :: loz1st =.true. ! first-time clim ozone data read flag
237 
238 ! --- optional extra top layer on top of low ceiling models
239  integer, parameter :: ltp = 0 ! no extra top layer
240 ! integer, parameter :: LTP = 1 ! add an extra top layer
241  logical, parameter :: lextop = (ltp > 0)
242 
243 ! --- publicly accessible module programs:
244 
245  public radinit, radupdate, grrad
246 
247 
248 ! =================
249  contains
250 ! =================
251 
258 !-----------------------------------
259  subroutine radinit &
260 !...................................
262 ! --- inputs:
263  & ( si, nlay, me )
264 ! --- outputs:
265 ! ( none )
266 
267 ! ================= subprogram documentation block ================ !
268 ! !
269 ! subprogram: radinit initialization of radiation calculations !
270 ! !
271 ! usage: call radinit !
272 ! !
273 ! attributes: !
274 ! language: fortran 90 !
275 ! machine: wcoss !
276 ! !
277 ! ==================== definition of variables ==================== !
278 ! !
279 ! input parameters: !
280 ! si : model vertical sigma interface !
281 ! NLAY : number of model vertical layers !
282 ! me : print control flag !
283 ! !
284 ! outputs: (none) !
285 ! !
286 ! external module variables: (in module physparam) !
287 ! isolar : solar constant cntrol flag !
288 ! = 0: use the old fixed solar constant in "physcon" !
289 ! =10: use the new fixed solar constant in "physcon" !
290 ! = 1: use noaa ann-mean tsi tbl abs-scale with cycle apprx!
291 ! = 2: use noaa ann-mean tsi tbl tim-scale with cycle apprx!
292 ! = 3: use cmip5 ann-mean tsi tbl tim-scale with cycl apprx!
293 ! = 4: use cmip5 mon-mean tsi tbl tim-scale with cycl apprx!
294 ! iaerflg : 3-digit aerosol flag (abc for volc, lw, sw) !
295 ! a:=0 use background stratospheric aerosol !
296 ! =1 include stratospheric vocanic aeros !
297 ! b:=0 no topospheric aerosol in lw radiation !
298 ! =1 compute tropspheric aero in 1 broad band for lw !
299 ! =2 compute tropspheric aero in multi bands for lw !
300 ! c:=0 no topospheric aerosol in sw radiation !
301 ! =1 include tropspheric aerosols for sw !
302 ! ico2flg : co2 data source control flag !
303 ! =0: use prescribed global mean co2 (old oper) !
304 ! =1: use observed co2 annual mean value only !
305 ! =2: use obs co2 monthly data with 2-d variation !
306 ! ictmflg : =yyyy#, external data ic time/date control flag !
307 ! = -2: same as 0, but superimpose seasonal cycle !
308 ! from climatology data set. !
309 ! = -1: use user provided external data for the !
310 ! forecast time, no extrapolation. !
311 ! = 0: use data at initial cond time, if not !
312 ! available, use latest, no extrapolation. !
313 ! = 1: use data at the forecast time, if not !
314 ! available, use latest and extrapolation. !
315 ! =yyyy0: use yyyy data for the forecast time, !
316 ! no further data extrapolation. !
317 ! =yyyy1: use yyyy data for the fcst. if needed, do !
318 ! extrapolation to match the fcst time. !
319 ! ioznflg : ozone data source control flag !
320 ! =0: use climatological ozone profile !
321 ! =1: use interactive ozone profile !
322 ! ialbflg : albedo scheme control flag !
323 ! =0: climatology, based on surface veg types !
324 ! =1: modis retrieval based surface albedo scheme !
325 ! iemsflg : emissivity scheme cntrl flag (ab 2-digit integer) !
326 ! a:=0 set sfc air/ground t same for lw radiation !
327 ! =1 set sfc air/ground t diff for lw radiation !
328 ! b:=0 use fixed sfc emissivity=1.0 (black-body) !
329 ! =1 use varying climtology sfc emiss (veg based) !
330 ! =2 future development (not yet) !
331 ! icldflg : cloud optical property scheme control flag !
332 ! =0: use diagnostic cloud scheme !
333 ! =1: use prognostic cloud scheme (default) !
334 ! icmphys : cloud microphysics scheme control flag !
335 ! =1 zhao/carr/sundqvist microphysics scheme !
336 ! =2 brad ferrier microphysics scheme !
337 ! =3 zhao/carr/sundqvist microphysics+pdf cloud & cnvc,cnvw!
338 ! iovrsw : control flag for cloud overlap in sw radiation !
339 ! iovrlw : control flag for cloud overlap in lw radiation !
340 ! =0: random overlapping clouds !
341 ! =1: max/ran overlapping clouds !
342 ! isubcsw : sub-column cloud approx control flag in sw radiation !
343 ! isubclw : sub-column cloud approx control flag in lw radiation !
344 ! =0: with out sub-column cloud approximation !
345 ! =1: mcica sub-col approx. prescribed random seed !
346 ! =2: mcica sub-col approx. provided random seed !
347 ! lsashal : shallow convection scheme flag !
348 ! lcrick : control flag for eliminating CRICK !
349 ! =t: apply layer smoothing to eliminate CRICK !
350 ! =f: do not apply layer smoothing !
351 ! lcnorm : control flag for in-cld condensate !
352 ! =t: normalize cloud condensate !
353 ! =f: not normalize cloud condensate !
354 ! lnoprec : precip effect in radiation flag (ferrier microphysics) !
355 ! =t: snow/rain has no impact on radiation !
356 ! =f: snow/rain has impact on radiation !
357 ! ivflip : vertical index direction control flag !
358 ! =0: index from toa to surface !
359 ! =1: index from surface to toa !
360 ! !
361 ! subroutines called: sol_init, aer_init, gas_init, cld_init, !
362 ! sfc_init, rlwinit, rswinit !
363 ! !
364 ! usage: call radinit !
365 ! !
366 ! =================================================================== !
367 !
368  implicit none
369 
370 ! --- inputs:
371  integer, intent(in) :: NLAY, me
372 
373  real (kind=kind_phys), intent(in) :: si(:)
374 
375 ! --- outputs: (none, to module variables)
376 
377 ! --- locals:
378 
379 !
380 !===> ... begin here
381 !
382 ! --- set up control variables
384  itsfc = iemsflg / 10 ! sfc air/ground temp control
385  loz1st = (ioznflg == 0) ! first-time clim ozone data read flag
386  month0 = 0
387  iyear0 = 0
388  monthd = 0
389 
390  if (me == 0) then
391 ! print *,' NEW RADIATION PROGRAM STRUCTURES -- SEP 01 2004'
392  print *,' NEW RADIATION PROGRAM STRUCTURES BECAME OPER. ', &
393  & ' May 01 2007'
394  print *, vtagrad !print out version tag
395  print *,' - Selected Control Flag settings: ICTMflg=',ictmflg, &
396  & ' ISOLar =',isolar, ' ICO2flg=',ico2flg,' IAERflg=',iaerflg, &
397  & ' IALBflg=',ialbflg,' IEMSflg=',iemsflg,' ICLDflg=',icldflg, &
398  & ' ICMPHYS=',icmphys,' IOZNflg=',ioznflg
399  print *,' IVFLIP=',ivflip,' IOVRSW=',iovrsw,' IOVRLW=',iovrlw, &
400  & ' ISUBCSW=',isubcsw,' ISUBCLW=',isubclw
401 ! write(0,*)' IVFLIP=',ivflip,' IOVRSW=',iovrsw,' IOVRLW=',iovrlw,&
402 ! & ' ISUBCSW=',isubcsw,' ISUBCLW=',isubclw
403  print *,' LSASHAL=',lsashal,' LCRICK=',lcrick,' LCNORM=',lcnorm,&
404  & ' LNOPREC=',lnoprec
405  print *,' LTP =',ltp,', add extra top layer =',lextop
406 
407  if ( ictmflg==0 .or. ictmflg==-2 ) then
408  print *,' Data usage is limited by initial condition!'
409  print *,' No volcanic aerosols'
410  endif
411 
412  if ( isubclw == 0 ) then
413  print *,' - ISUBCLW=',isubclw,' No McICA, use grid ', &
414  & 'averaged cloud in LW radiation'
415  elseif ( isubclw == 1 ) then
416  print *,' - ISUBCLW=',isubclw,' Use McICA with fixed ', &
417  & 'permutation seeds for LW random number generator'
418  elseif ( isubclw == 2 ) then
419  print *,' - ISUBCLW=',isubclw,' Use McICA with random ', &
420  & 'permutation seeds for LW random number generator'
421  else
422  print *,' - ERROR!!! ISUBCLW=',isubclw,' is not a ', &
423  & 'valid option '
424  stop
425  endif
426 
427  if ( isubcsw == 0 ) then
428  print *,' - ISUBCSW=',isubcsw,' No McICA, use grid ', &
429  & 'averaged cloud in SW radiation'
430  elseif ( isubcsw == 1 ) then
431  print *,' - ISUBCSW=',isubcsw,' Use McICA with fixed ', &
432  & 'permutation seeds for SW random number generator'
433  elseif ( isubcsw == 2 ) then
434  print *,' - ISUBCSW=',isubcsw,' Use McICA with random ', &
435  & 'permutation seeds for SW random number generator'
436  else
437  print *,' - ERROR!!! ISUBCSW=',isubcsw,' is not a ', &
438  & 'valid option '
439  stop
440  endif
441 
442  if ( isubcsw /= isubclw ) then
443  print *,' - *** Notice *** ISUBCSW /= ISUBCLW !!!', &
444  & isubcsw, isubclw
445  endif
446  endif
447 
457 ! Initialization
458 
459  call sol_init ( me ) ! --- ... astronomy initialization routine
460 
461  call aer_init ( nlay, me ) ! --- ... aerosols initialization routine
462 
463  call gas_init ( me ) ! --- ... co2 and other gases initialization routine
464 
465  call sfc_init ( me ) ! --- ... surface initialization routine
466 
467  call cld_init ( si, nlay, me) ! --- ... cloud initialization routine
468 
469  call rlwinit ( me ) ! --- ... lw radiation initialization routine
470 
471  call rswinit ( me ) ! --- ... sw radiation initialization routine
472 !
473  return
474 !...................................
475  end subroutine radinit
476 !-----------------------------------
478 
492 !-----------------------------------
493  subroutine radupdate &
494 !...................................
496 ! --- inputs:
497  & ( idate,jdate,deltsw,deltim,lsswr, me, &
498 ! --- outputs:
499  & slag,sdec,cdec,solcon &
500  & )
501 
502 ! ================= subprogram documentation block ================ !
503 ! !
504 ! subprogram: radupdate calls many update subroutines to check and !
505 ! update radiation required but time varying data sets and module !
506 ! variables. !
507 ! !
508 ! usage: call radupdate !
509 ! !
510 ! attributes: !
511 ! language: fortran 90 !
512 ! machine: ibm sp !
513 ! !
514 ! ==================== definition of variables ==================== !
515 ! !
516 ! input parameters: !
517 ! idate(8) : ncep absolute date and time of initial condition !
518 ! (yr, mon, day, t-zone, hr, min, sec, mil-sec) !
519 ! jdate(8) : ncep absolute date and time at fcst time !
520 ! (yr, mon, day, t-zone, hr, min, sec, mil-sec) !
521 ! deltsw : sw radiation calling frequency in seconds !
522 ! deltim : model timestep in seconds !
523 ! lsswr : logical flags for sw radiation calculations !
524 ! me : print control flag !
525 ! !
526 ! outputs: !
527 ! slag : equation of time in radians !
528 ! sdec, cdec : sin and cos of the solar declination angle !
529 ! solcon : sun-earth distance adjusted solar constant (w/m2) !
530 ! !
531 ! external module variables: !
532 ! isolar : solar constant cntrl (in module physparam) !
533 ! = 0: use the old fixed solar constant in "physcon" !
534 ! =10: use the new fixed solar constant in "physcon" !
535 ! = 1: use noaa ann-mean tsi tbl abs-scale with cycle apprx!
536 ! = 2: use noaa ann-mean tsi tbl tim-scale with cycle apprx!
537 ! = 3: use cmip5 ann-mean tsi tbl tim-scale with cycl apprx!
538 ! = 4: use cmip5 mon-mean tsi tbl tim-scale with cycl apprx!
539 ! ictmflg : =yyyy#, external data ic time/date control flag !
540 ! = -2: same as 0, but superimpose seasonal cycle !
541 ! from climatology data set. !
542 ! = -1: use user provided external data for the !
543 ! forecast time, no extrapolation. !
544 ! = 0: use data at initial cond time, if not !
545 ! available, use latest, no extrapolation. !
546 ! = 1: use data at the forecast time, if not !
547 ! available, use latest and extrapolation. !
548 ! =yyyy0: use yyyy data for the forecast time, !
549 ! no further data extrapolation. !
550 ! =yyyy1: use yyyy data for the fcst. if needed, do !
551 ! extrapolation to match the fcst time. !
552 ! !
553 ! module variables: !
554 ! loz1st : first-time clim ozone data read flag !
555 ! !
556 ! subroutines called: sol_update, aer_update, gas_update !
557 ! !
558 ! =================================================================== !
559 !
560  implicit none
561 
562 ! --- inputs:
563  integer, intent(in) :: idate(:), jdate(:), me
564  logical, intent(in) :: lsswr
565 
566  real (kind=kind_phys), intent(in) :: deltsw, deltim
567 
568 ! --- outputs:
569  real (kind=kind_phys), intent(out) :: slag, sdec, cdec, solcon
570 
571 ! --- locals:
572  integer :: iyear, imon, iday, ihour
573  integer :: kyear, kmon, kday, khour
574 
575  logical :: lmon_chg ! month change flag
576  logical :: lco2_chg ! cntrl flag for updating co2 data
577  logical :: lsol_chg ! cntrl flag for updating solar constant
578 !
579 !===> ... begin here
580 !
582 ! --- ... time stamp at fcst time
583 
584  iyear = jdate(1)
585  imon = jdate(2)
586  iday = jdate(3)
587  ihour = jdate(5)
588 
589 ! --- ... set up time stamp used for green house gases (** currently co2 only)
590 
591  if ( ictmflg==0 .or. ictmflg==-2 ) then ! get external data at initial condition time
592  kyear = idate(1)
593  kmon = idate(2)
594  kday = idate(3)
595  khour = idate(5)
596  else ! get external data at fcst or specified time
597  kyear = iyear
598  kmon = imon
599  kday = iday
600  khour = ihour
601  endif ! end if_ictmflg_block
602 
603  if ( month0 /= imon ) then
604  lmon_chg = .true.
605  month0 = imon
606  else
607  lmon_chg = .false.
608  endif
611 ! --- ... call astronomy update routine, yearly update, no time interpolation
612 
613  if (lsswr) then
614 
615  if ( isolar == 0 .or. isolar == 10 ) then
616  lsol_chg = .false.
617  elseif ( iyear0 /= iyear ) then
618  lsol_chg = .true.
619  else
620  lsol_chg = ( isolar==4 .and. lmon_chg )
621  endif
622  iyear0 = iyear
623 
624  call sol_update &
625 ! --- inputs:
626  & ( jdate,kyear,deltsw,deltim,lsol_chg, me, &
627 ! --- outputs:
628  & slag,sdec,cdec,solcon &
629  & )
630 
631  endif ! end_if_lsswr_block
634 ! --- ... call aerosols update routine, monthly update, no time interpolation
635 
636  if ( lmon_chg ) then
637  call aer_update ( iyear, imon, me )
638  endif
639 
642 ! --- ... call co2 and other gases update routine
643 
644  if ( monthd /= kmon ) then
645  monthd = kmon
646  lco2_chg = .true.
647  else
648  lco2_chg = .false.
649  endif
650 
651  call gas_update ( kyear,kmon,kday,khour,loz1st,lco2_chg, me )
652 
653  if ( loz1st ) loz1st = .false.
654 
656 ! --- ... call surface update routine (currently not needed)
657 ! call sfc_update ( iyear, imon, me )
658 
660 ! --- ... call clouds update routine (currently not needed)
661 ! call cld_update ( iyear, imon, me )
662 !
663  return
664 !...................................
665  end subroutine radupdate
666 !-----------------------------------
668 
808 !-----------------------------------
809  subroutine grrad &
810 !...................................
811 ! --- inputs:
812  & ( prsi,prsl,prslk,tgrs,qgrs,tracer,vvl,slmsk, &
813  & xlon,xlat,tsfc,snowd,sncovr,snoalb,zorl,hprim, &
814  & alvsf,alnsf,alvwf,alnwf,facsf,facwf,fice,tisfc, &
815  & sinlat,coslat,solhr,jdate,solcon, &
816  & cv,cvt,cvb,fcice,frain,rrime,flgmin, &
817  & icsdsw,icsdlw, ntcw,ncld,ntoz, ntrac,nfxr, &
818  & dtlw,dtsw, lsswr,lslwr,lssav, shoc_cld, &
819  & ix, im, lm, me, lprnt, ipt, kdt, deltaq,sup,cnvw,cnvc, &
820 ! --- outputs:
821  & htrsw,topfsw,sfcfsw,dswcmp,uswcmp,sfalb,coszen,coszdg, &
822  & htrlw,topflw,sfcflw,tsflw,semis,cldcov, &
823 ! --- input/output:
824  & fluxr &
825 !! --- optional outputs:
826  &, htrlw0,htrsw0,htrswb,htrlwb &
827  & )
828 
829 ! ================= subprogram documentation block ================ !
830 ! !
831 ! this program is the driver of radiation calculation subroutines. * !
832 ! It sets up profile variables for radiation input, including * !
833 ! clouds, surface albedos, atmospheric aerosols, ozone, etc. * !
834 ! * !
835 ! usage: call grrad * !
836 ! * !
837 ! subprograms called: * !
838 ! setalb, setemis, setaer, getozn, getgases, * !
839 ! progcld1, progcld2, diagcds, * !
840 ! swrad, lwrad, fpvs * !
841 ! * !
842 ! attributes: * !
843 ! language: fortran 90 * !
844 ! machine: ibm-sp, sgi * !
845 ! * !
846 ! * !
847 ! ==================== definition of variables ==================== !
848 ! !
849 ! input variables: !
850 ! prsi (IX,LM+1) : model level pressure in Pa !
851 ! prsl (IX,LM) : model layer mean pressure Pa !
852 ! prslk (IX,LM) : exner function = (p/p0)**rocp !
853 ! tgrs (IX,LM) : model layer mean temperature in k !
854 ! qgrs (IX,LM) : layer specific humidity in gm/gm !
855 ! tracer(IX,LM,NTRAC):layer prognostic tracer amount/mixing-ratio !
856 ! incl: oz, cwc, aeros, etc. !
857 ! vvl (IX,LM) : layer mean vertical velocity in pa/sec !
858 ! slmsk (IM) : sea/land mask array (sea:0,land:1,sea-ice:2) !
859 ! xlon (IM) : grid longitude in radians, ok for both 0->2pi !
860 ! or -pi -> +pi ranges !
861 ! xlat (IM) : grid latitude in radians, default to pi/2 -> !
862 ! -pi/2 range, otherwise adj in subr called !
863 ! tsfc (IM) : surface temperature in k !
864 ! snowd (IM) : snow depth water equivalent in mm !
865 ! sncovr(IM) : snow cover in fraction !
866 ! snoalb(IM) : maximum snow albedo in fraction !
867 ! zorl (IM) : surface roughness in cm !
868 ! hprim (IM) : topographic standard deviation in m !
869 ! alvsf (IM) : mean vis albedo with strong cosz dependency !
870 ! alnsf (IM) : mean nir albedo with strong cosz dependency !
871 ! alvwf (IM) : mean vis albedo with weak cosz dependency !
872 ! alnwf (IM) : mean nir albedo with weak cosz dependency !
873 ! facsf (IM) : fractional coverage with strong cosz dependen !
874 ! facwf (IM) : fractional coverage with weak cosz dependency !
875 ! fice (IM) : ice fraction over open water grid !
876 ! tisfc (IM) : surface temperature over ice fraction !
877 ! sinlat(IM) : sine of the grids' corresponding latitudes !
878 ! coslat(IM) : cosine of the grids' corresponding latitudes !
879 ! solhr : hour time after 00z at the t-stepe !
880 ! jdate (8) : current forecast date and time !
881 ! (yr, mon, day, t-zone, hr, min, sec, mil-sec) !
882 ! solcon : solar constant (sun-earth distant adjusted) !
883 ! cv (IM) : fraction of convective cloud !
884 ! cvt, cvb (IM) : convective cloud top/bottom pressure in pa !
885 ! fcice : fraction of cloud ice (in ferrier scheme) !
886 ! frain : fraction of rain water (in ferrier scheme) !
887 ! rrime : mass ratio of total to unrimed ice ( >= 1 ) !
888 ! flgmin : minimim large ice fraction !
889 ! icsdsw/icsdlw : auxiliary cloud control arrays passed to main !
890 ! (IM) radiations. if isubcsw/isubclw (input to init) !
891 ! are set to 2, the arrays contains provided !
892 ! random seeds for sub-column clouds generators !
893 ! ntcw : =0 no cloud condensate calculated !
894 ! >0 array index location for cloud condensate !
895 ! ncld : only used when ntcw .gt. 0 !
896 ! ntoz : =0 climatological ozone profile !
897 ! >0 interactive ozone profile !
898 ! NTRAC : dimension veriable for array oz !
899 ! NFXR : second dimension of input/output array fluxr !
900 ! dtlw, dtsw : time duration for lw/sw radiation call in sec !
901 ! lsswr, lslwr : logical flags for sw/lw radiation calls !
902 ! lssav : logical flag for store 3-d cloud field !
903 ! IX,IM : horizontal dimention and num of used points !
904 ! LM : vertical layer dimension !
905 ! me : control flag for parallel process !
906 ! lprnt : control flag for diagnostic print out !
907 ! ipt : index for diagnostic printout point !
908 ! kdt : time-step number !
909 ! deltaq : half width of uniform total water distribution !
910 ! sup : supersaturation in pdf cloud when t is very low!
911 ! cnvw : layer convective cloud water !
912 ! cnvc : layer convective cloud cover !
913 ! !
914 ! output variables: !
915 ! htrsw (IX,LM) : total sky sw heating rate in k/sec !
916 ! topfsw(IM) : sw radiation fluxes at toa, components: !
917 ! (check module_radsw_parameters for definition) !
918 ! %upfxc - total sky upward sw flux at toa (w/m**2) !
919 ! %dnflx - total sky downward sw flux at toa (w/m**2) !
920 ! %upfx0 - clear sky upward sw flux at toa (w/m**2) !
921 ! sfcfsw(IM) : sw radiation fluxes at sfc, components: !
922 ! (check module_radsw_parameters for definition) !
923 ! %upfxc - total sky upward sw flux at sfc (w/m**2) !
924 ! %dnfxc - total sky downward sw flux at sfc (w/m**2) !
925 ! %upfx0 - clear sky upward sw flux at sfc (w/m**2) !
926 ! %dnfx0 - clear sky downward sw flux at sfc (w/m**2) !
927 ! dswcmp(IX,4) : dn sfc sw spectral components: !
928 ! ( :, 1) - total sky sfc downward nir direct flux !
929 ! ( :, 2) - total sky sfc downward nir diffused flux !
930 ! ( :, 3) - total sky sfc downward uv+vis direct flux !
931 ! ( :, 4) - total sky sfc downward uv+vis diff flux !
932 ! uswcmp(IX,4) : up sfc sw spectral components: !
933 ! ( :, 1) - total sky sfc upward nir direct flux !
934 ! ( :, 2) - total sky sfc upward nir diffused flux !
935 ! ( :, 3) - total sky sfc upward uv+vis direct flux !
936 ! ( :, 4) - total sky sfc upward uv+vis diff flux !
937 ! sfalb (IM) : mean surface diffused sw albedo !
938 ! coszen(IM) : mean cos of zenith angle over rad call period !
939 ! coszdg(IM) : daytime mean cosz over rad call period !
940 ! htrlw (IX,LM) : total sky lw heating rate in k/sec !
941 ! topflw(IM) : lw radiation fluxes at top, component: !
942 ! (check module_radlw_paramters for definition) !
943 ! %upfxc - total sky upward lw flux at toa (w/m**2) !
944 ! %upfx0 - clear sky upward lw flux at toa (w/m**2) !
945 ! sfcflw(IM) : lw radiation fluxes at sfc, component: !
946 ! (check module_radlw_paramters for definition) !
947 ! %upfxc - total sky upward lw flux at sfc (w/m**2) !
948 ! %upfx0 - clear sky upward lw flux at sfc (w/m**2) !
949 ! %dnfxc - total sky downward lw flux at sfc (w/m**2) !
950 ! %dnfx0 - clear sky downward lw flux at sfc (w/m**2) !
951 ! semis (IM) : surface lw emissivity in fraction !
952 ! cldcov(IX,LM) : 3-d cloud fraction !
953 ! tsflw (IM) : surface air temp during lw calculation in k !
954 ! !
955 ! input and output variables: !
956 ! fluxr (IX,NFXR) : to save time accumulated 2-d fields defined as:!
957 ! 1 - toa total sky upwd lw radiation flux !
958 ! 2 - toa total sky upwd sw radiation flux !
959 ! 3 - sfc total sky upwd sw radiation flux !
960 ! 4 - sfc total sky dnwd sw radiation flux !
961 ! 5 - high domain cloud fraction !
962 ! 6 - mid domain cloud fraction !
963 ! 7 - low domain cloud fraction !
964 ! 8 - high domain mean cloud top pressure !
965 ! 9 - mid domain mean cloud top pressure !
966 ! 10 - low domain mean cloud top pressure !
967 ! 11 - high domain mean cloud base pressure !
968 ! 12 - mid domain mean cloud base pressure !
969 ! 13 - low domain mean cloud base pressure !
970 ! 14 - high domain mean cloud top temperature !
971 ! 15 - mid domain mean cloud top temperature !
972 ! 16 - low domain mean cloud top temperature !
973 ! 17 - total cloud fraction !
974 ! 18 - boundary layer domain cloud fraction !
975 ! 19 - sfc total sky dnwd lw radiation flux !
976 ! 20 - sfc total sky upwd lw radiation flux !
977 ! 21 - sfc total sky dnwd sw uv-b radiation flux !
978 ! 22 - sfc clear sky dnwd sw uv-b radiation flux !
979 ! 23 - toa incoming solar radiation flux !
980 ! 24 - sfc vis beam dnwd sw radiation flux !
981 ! 25 - sfc vis diff dnwd sw radiation flux !
982 ! 26 - sfc nir beam dnwd sw radiation flux !
983 ! 27 - sfc nir diff dnwd sw radiation flux !
984 ! 28 - toa clear sky upwd lw radiation flux !
985 ! 29 - toa clear sky upwd sw radiation flux !
986 ! 30 - sfc clear sky dnwd lw radiation flux !
987 ! 31 - sfc clear sky upwd sw radiation flux !
988 ! 32 - sfc clear sky dnwd sw radiation flux !
989 ! 33 - sfc clear sky upwd lw radiation flux !
990 !optional 34 - aeros opt depth at 550nm (all components) !
991 ! 35 - aeros opt depth at 550nm for du component !
992 ! 36 - aeros opt depth at 550nm for bc component !
993 ! 37 - aeros opt depth at 550nm for oc component !
994 ! 38 - aeros opt depth at 550nm for su component !
995 ! 39 - aeros opt depth at 550nm for ss component !
996 ! !
997 ! optional output variables: !
998 ! htrswb(IX,LM,NBDSW) : spectral band total sky sw heating rate !
999 ! htrlwb(IX,LM,NBDLW) : spectral band total sky lw heating rate !
1000 ! !
1001 ! !
1002 ! definitions of internal variable arrays: !
1003 ! !
1004 ! 1. fixed gases: (defined in 'module_radiation_gases') !
1005 ! gasvmr(:,:,1) - co2 volume mixing ratio !
1006 ! gasvmr(:,:,2) - n2o volume mixing ratio !
1007 ! gasvmr(:,:,3) - ch4 volume mixing ratio !
1008 ! gasvmr(:,:,4) - o2 volume mixing ratio !
1009 ! gasvmr(:,:,5) - co volume mixing ratio !
1010 ! gasvmr(:,:,6) - cf11 volume mixing ratio !
1011 ! gasvmr(:,:,7) - cf12 volume mixing ratio !
1012 ! gasvmr(:,:,8) - cf22 volume mixing ratio !
1013 ! gasvmr(:,:,9) - ccl4 volume mixing ratio !
1014 ! !
1015 ! 2. cloud profiles: (defined in 'module_radiation_clouds') !
1016 ! --- for prognostic cloud --- !
1017 ! clouds(:,:,1) - layer total cloud fraction !
1018 ! clouds(:,:,2) - layer cloud liq water path !
1019 ! clouds(:,:,3) - mean effective radius for liquid cloud !
1020 ! clouds(:,:,4) - layer cloud ice water path !
1021 ! clouds(:,:,5) - mean effective radius for ice cloud !
1022 ! clouds(:,:,6) - layer rain drop water path !
1023 ! clouds(:,:,7) - mean effective radius for rain drop !
1024 ! clouds(:,:,8) - layer snow flake water path !
1025 ! clouds(:,:,9) - mean effective radius for snow flake !
1026 ! --- for diagnostic cloud --- !
1027 ! clouds(:,:,1) - layer total cloud fraction !
1028 ! clouds(:,:,2) - layer cloud optical depth !
1029 ! clouds(:,:,3) - layer cloud single scattering albedo !
1030 ! clouds(:,:,4) - layer cloud asymmetry factor !
1031 ! !
1032 ! 3. surface albedo: (defined in 'module_radiation_surface') !
1033 ! sfcalb( :,1 ) - near ir direct beam albedo !
1034 ! sfcalb( :,2 ) - near ir diffused albedo !
1035 ! sfcalb( :,3 ) - uv+vis direct beam albedo !
1036 ! sfcalb( :,4 ) - uv+vis diffused albedo !
1037 ! !
1038 ! 4. sw aerosol profiles: (defined in 'module_radiation_aerosols') !
1039 ! faersw(:,:,:,1)- sw aerosols optical depth !
1040 ! faersw(:,:,:,2)- sw aerosols single scattering albedo !
1041 ! faersw(:,:,:,3)- sw aerosols asymmetry parameter !
1042 ! !
1043 ! 5. lw aerosol profiles: (defined in 'module_radiation_aerosols') !
1044 ! faerlw(:,:,:,1)- lw aerosols optical depth !
1045 ! faerlw(:,:,:,2)- lw aerosols single scattering albedo !
1046 ! faerlw(:,:,:,3)- lw aerosols asymmetry parameter !
1047 ! !
1048 ! 6. sw fluxes at toa: (defined in 'module_radsw_main') !
1049 ! (topfsw_type -- derived data type for toa rad fluxes) !
1050 ! topfsw(:)%upfxc - total sky upward flux at toa !
1051 ! topfsw(:)%dnfxc - total sky downward flux at toa !
1052 ! topfsw(:)%upfx0 - clear sky upward flux at toa !
1053 ! !
1054 ! 7. lw fluxes at toa: (defined in 'module_radlw_main') !
1055 ! (topflw_type -- derived data type for toa rad fluxes) !
1056 ! topflw(:)%upfxc - total sky upward flux at toa !
1057 ! topflw(:)%upfx0 - clear sky upward flux at toa !
1058 ! !
1059 ! 8. sw fluxes at sfc: (defined in 'module_radsw_main') !
1060 ! (sfcfsw_type -- derived data type for sfc rad fluxes) !
1061 ! sfcfsw(:)%upfxc - total sky upward flux at sfc !
1062 ! sfcfsw(:)%dnfxc - total sky downward flux at sfc !
1063 ! sfcfsw(:)%upfx0 - clear sky upward flux at sfc !
1064 ! sfcfsw(:)%dnfx0 - clear sky downward flux at sfc !
1065 ! !
1066 ! 9. lw fluxes at sfc: (defined in 'module_radlw_main') !
1067 ! (sfcflw_type -- derived data type for sfc rad fluxes) !
1068 ! sfcflw(:)%upfxc - total sky upward flux at sfc !
1069 ! sfcflw(:)%dnfxc - total sky downward flux at sfc !
1070 ! sfcflw(:)%dnfx0 - clear sky downward flux at sfc !
1071 ! !
1072 !! optional radiation outputs: !
1073 !! 10. sw flux profiles: (defined in 'module_radsw_main') !
1074 !! (profsw_type -- derived data type for rad vertical profiles) !
1075 !! fswprf(:,:)%upfxc - total sky upward flux !
1076 !! fswprf(:,:)%dnfxc - total sky downward flux !
1077 !! fswprf(:,:)%upfx0 - clear sky upward flux !
1078 !! fswprf(:,:)%dnfx0 - clear sky downward flux !
1079 !! !
1080 !! 11. lw flux profiles: (defined in 'module_radlw_main') !
1081 !! (proflw_type -- derived data type for rad vertical profiles) !
1082 !! flwprf(:,:)%upfxc - total sky upward flux !
1083 !! flwprf(:,:)%dnfxc - total sky downward flux !
1084 !! flwprf(:,:)%upfx0 - clear sky upward flux !
1085 !! flwprf(:,:)%dnfx0 - clear sky downward flux !
1086 !! !
1087 !! 12. sw sfc components: (defined in 'module_radsw_main') !
1088 !! (cmpfsw_type -- derived data type for component sfc fluxes) !
1089 !! scmpsw(:)%uvbfc - total sky downward uv-b flux at sfc !
1090 !! scmpsw(:)%uvbf0 - clear sky downward uv-b flux at sfc !
1091 !! scmpsw(:)%nirbm - total sky sfc downward nir direct flux !
1092 !! scmpsw(:)%nirdf - total sky sfc downward nir diffused flux !
1093 !! scmpsw(:)%visbm - total sky sfc downward uv+vis direct flx !
1094 !! scmpsw(:)%visdf - total sky sfc downward uv+vis diff flux !
1095 ! !
1096 ! external module variables: !
1097 ! ivflip : control flag for in/out vertical indexing !
1098 ! =0 index from toa to surface !
1099 ! =1 index from surface to toa !
1100 ! icmphys : cloud microphysics scheme control flag !
1101 ! =1 zhao/carr/sundqvist microphysics scheme !
1102 ! =2 brad ferrier microphysics scheme !
1103 ! =3 zhao/carr/sundqvist microphysics +pdf cloud !
1104 ! !
1105 ! module variables: !
1106 ! itsfc : =0 use same sfc skin-air/ground temp !
1107 ! =1 use diff sfc skin-air/ground temp (not yet) !
1108 ! !
1109 ! ====================== end of definitions ======================= !
1110 !
1111  implicit none
1112 
1113 ! --- inputs: (for rank>1 arrays, horizontal dimensioned by IX)
1114  integer, intent(in) :: IX,IM, LM, NTRAC, NFXR, me, &
1115  & ntoz, ntcw, ncld, ipt, kdt
1116  integer, intent(in) :: icsdsw(im), icsdlw(im), jdate(8)
1117 
1118  logical, intent(in) :: lsswr, lslwr, lssav, lprnt, shoc_cld
1119 
1120  real (kind=kind_phys), dimension(IX,LM+1), intent(in) :: prsi
1121 
1122  real (kind=kind_phys), dimension(IX,LM), intent(in) :: prsl, &
1123  & prslk, tgrs, qgrs, vvl, fcice, frain, rrime, deltaq, cnvw, &
1124  & cnvc
1125  real (kind=kind_phys), dimension(IM), intent(in) :: flgmin
1126  real(kind=kind_phys), intent(in) :: sup
1127 
1128  real (kind=kind_phys), dimension(IM), intent(in) :: slmsk, &
1129  & xlon, xlat, tsfc, snowd, zorl, hprim, alvsf, alnsf, alvwf, &
1130  & alnwf, facsf, facwf, cv, cvt, cvb, fice, tisfc, &
1131  & sncovr, snoalb, sinlat, coslat
1132 
1133  real (kind=kind_phys), intent(in) :: solcon, dtlw, dtsw, solhr, &
1134  & tracer(IX,LM,NTRAC)
1135 
1136  real (kind=kind_phys), dimension(IX,LM),intent(inout):: cldcov
1137 
1138 ! --- outputs: (horizontal dimensioned by IX)
1139  real (kind=kind_phys), dimension(IX,LM),intent(out):: htrsw,htrlw
1140 
1141  real (kind=kind_phys), dimension(IX,4), intent(out) :: dswcmp, &
1142  & uswcmp
1143 
1144  real (kind=kind_phys), dimension(IM), intent(out):: tsflw, &
1145  & sfalb, semis, coszen, coszdg
1146 
1147  type(topfsw_type), dimension(IM), intent(out) :: topfsw
1148  type(sfcfsw_type), dimension(IM), intent(out) :: sfcfsw
1149 
1150  type(topflw_type), dimension(IM), intent(out) :: topflw
1151  type(sfcflw_type), dimension(IM), intent(out) :: sfcflw
1152 
1153 ! --- variables are for both input and output:
1154  real (kind=kind_phys), intent(inout) :: fluxr(ix,nfxr)
1155 
1156 !! --- optional outputs:
1157  real (kind=kind_phys), dimension(IX,LM,NBDSW), optional, &
1158  & intent(out) :: htrswb
1159  real (kind=kind_phys), dimension(IX,LM,NBDLW), optional, &
1160  & intent(out) :: htrlwb
1161  real (kind=kind_phys), dimension(ix,lm), optional, &
1162  & intent(out) :: htrlw0
1163  real (kind=kind_phys), dimension(ix,lm), optional, &
1164  & intent(out) :: htrsw0
1165 
1166 ! --- local variables: (horizontal dimensioned by IM)
1167  real (kind=kind_phys), dimension(IM,LM+1+LTP):: plvl, tlvl
1168 
1169  real (kind=kind_phys), dimension(IM,LM+LTP) :: plyr, tlyr, qlyr, &
1170  & olyr, rhly, qstl, vvel, clw, prslk1, tem2da, tem2db, tvly
1171 
1172  real (kind=kind_phys), dimension(IM) :: tsfa, cvt1, cvb1, tem1d, &
1173  & sfcemis, tsfg, tskn
1174 
1175  real (kind=kind_phys), dimension(IM,LM+LTP,NF_CLDS) :: clouds
1176  real (kind=kind_phys), dimension(IM,LM+LTP,NF_VGAS) :: gasvmr
1177  real (kind=kind_phys), dimension(IM, NF_ALBD) :: sfcalb
1178  real (kind=kind_phys), dimension(IM, NSPC1) :: aerodp
1179  real (kind=kind_phys), dimension(IM,LM+LTP,NTRAC) :: tracer1
1180 
1181  real (kind=kind_phys), dimension(IM,LM+LTP,NBDSW,NF_AESW)::faersw
1182  real (kind=kind_phys), dimension(IM,LM+LTP,NBDLW,NF_AELW)::faerlw
1183 
1184  real (kind=kind_phys), dimension(IM,LM+LTP) :: htswc
1185  real (kind=kind_phys), dimension(IM,LM+LTP) :: htlwc
1186 
1187  real (kind=kind_phys), dimension(IM,LM+LTP) :: gcice, grain, grime
1188 
1189 !! --- may be used for optional sw/lw outputs:
1190 !! take out "!!" as needed
1191  real (kind=kind_phys), dimension(IM,LM+LTP) :: htsw0
1192 !! type (profsw_type), dimension(IM,LM+1+LTP) :: fswprf
1193  type(cmpfsw_type), dimension(IM) :: scmpsw
1194  real (kind=kind_phys), dimension(IM,LM+LTP,NBDSW) :: htswb
1195 
1196  real (kind=kind_phys), dimension(IM,LM+LTP) :: htlw0
1197 !! type (proflw_type), dimension(IM,LM+1+LTP) :: flwprf
1198  real (kind=kind_phys), dimension(IM,LM+LTP,NBDLW) :: htlwb
1199 
1200  real (kind=kind_phys) :: raddt, es, qs, delt, tem0d, cldsa(im,5)
1201 
1202  integer :: i, j, k, k1, lv, icec, itop, ibtc, nday, idxday(im), &
1203  & mbota(IM,3), mtopa(IM,3), LP1, nb, LMK, LMP, kd, lla, llb, &
1204  & lya, lyb, kt, kb
1205 
1206 ! --- for debug test use
1207 ! real (kind=kind_phys) :: temlon, temlat, alon, alat
1208 ! integer :: ipt
1209 ! logical :: lprnt1
1210 
1211 !
1212 !===> ... begin here
1213 !
1214  lp1 = lm + 1 ! num of in/out levels
1215 
1216 ! --- ... set local /level/layer indexes corresponding to in/out variables
1217 
1218  lmk = lm + ltp ! num of local layers
1219  lmp = lmk + 1 ! num of local levels
1220 
1221  if ( lextop ) then
1222  if ( ivflip == 1 ) then ! vertical from sfc upward
1223  kd = 0 ! index diff between in/out and local
1224  kt = 1 ! index diff between lyr and upper bound
1225  kb = 0 ! index diff between lyr and lower bound
1226  lla = lmk ! local index at the 2nd level from top
1227  llb = lmp ! local index at toa level
1228  lya = lm ! local index for the 2nd layer from top
1229  lyb = lp1 ! local index for the top layer
1230  else ! vertical from toa downward
1231  kd = 1 ! index diff between in/out and local
1232  kt = 0 ! index diff between lyr and upper bound
1233  kb = 1 ! index diff between lyr and lower bound
1234  lla = 2 ! local index at the 2nd level from top
1235  llb = 1 ! local index at toa level
1236  lya = 2 ! local index for the 2nd layer from top
1237  lyb = 1 ! local index for the top layer
1238  endif ! end if_ivflip_block
1239  else
1240  kd = 0
1241  if ( ivflip == 1 ) then ! vertical from sfc upward
1242  kt = 1 ! index diff between lyr and upper bound
1243  kb = 0 ! index diff between lyr and lower bound
1244  else ! vertical from toa downward
1245  kt = 0 ! index diff between lyr and upper bound
1246  kb = 1 ! index diff between lyr and lower bound
1247  endif ! end if_ivflip_block
1248  endif ! end if_lextop_block
1249 
1250  raddt = min(dtsw, dtlw)
1251 
1252 ! --- ... for debug test
1253 ! alon = 120.0
1254 ! alat = 29.5
1255 ! ipt = 0
1256 ! do i = 1, IM
1257 ! temlon = xlon(i) * 57.29578
1258 ! if (temlon < 0.0) temlon = temlon + 360.0
1259 ! temlat = xlat(i) * 57.29578
1260 ! lprnt1 = abs(temlon-alon) < 1.1 .and. abs(temlat-alat) < 1.1
1261 ! if ( lprnt1 ) then
1262 ! ipt = i
1263 ! exit
1264 ! endif
1265 ! enddo
1266 
1267 ! print *,' in grrad : raddt=',raddt
1269 ! --- ... setup surface ground temp and ground/air skin temp if required
1270 
1271  if ( itsfc == 0 ) then ! use same sfc skin-air/ground temp
1272  do i = 1, im
1273  tskn(i) = tsfc(i)
1274  tsfg(i) = tsfc(i)
1275  enddo
1276  else ! use diff sfc skin-air/ground temp
1277  do i = 1, im
1278 !! tskn(i) = ta (i) ! not yet
1279 !! tsfg(i) = tg (i) ! not yet
1280  tskn(i) = tsfc(i)
1281  tsfg(i) = tsfc(i)
1282  enddo
1283  endif
1285 ! --- ... prepare atmospheric profiles for radiation input
1286 !
1287 ! if (im > ipt) then
1288 ! write(0,*)' prsi=',prsi(ipt,1:10)
1289 ! write(0,*)' prsi=',prsl(ipt,1:10)
1290 ! write(0,*)' tgrs=',tgrs(ipt,1:10)
1291 ! endif
1292 
1293 ! convert pressure unit from pa to mb
1294  do k = 1, lm
1295  k1 = k + kd
1296  do i = 1, im
1297 ! plvl(i,k1) = 10.0 * prsi(i,k) ! cb (kpa) to mb (hpa)
1298 ! plyr(i,k1) = 10.0 * prsl(i,k) ! cb (kpa) to mb (hpa)
1299  plvl(i,k1) = 0.01 * prsi(i,k) ! pa to mb (hpa)
1300  plyr(i,k1) = 0.01 * prsl(i,k) ! pa to mb (hpa)
1301  tlyr(i,k1) = tgrs(i,k)
1302  prslk1(i,k1) = prslk(i,k)
1303 
1305 ! --- ... compute relative humidity
1306 ! es = min( prsl(i,k), 0.001 * fpvs( tgrs(i,k) ) ) ! fpvs in pa
1307  es = min( prsl(i,k), fpvs( tgrs(i,k) ) ) ! fpvs and prsl in pa
1308  qs = max( qmin, eps * es / (prsl(i,k) + epsm1*es) )
1309  rhly(i,k1) = max( 0.0, min( 1.0, max(qmin, qgrs(i,k))/qs ) )
1310  qstl(i,k1) = qs
1311  enddo
1312  enddo
1313  do j = 1, ntrac
1314  do k = 1, lm
1315  k1 = k + kd
1316  do i = 1, im
1317  tracer1(i,k1,j) = tracer(i,k,j)
1318  enddo
1319  enddo
1320  enddo
1321 
1322  do i = 1, im
1323 ! plvl(i,LP1+kd) = 10.0 * prsi(i,LP1) ! cb (kpa) to mb (hpa
1324  plvl(i,lp1+kd) = 0.01 * prsi(i,lp1) ! pa to mb (hpa)
1325  enddo
1326 
1327  if ( lextop ) then ! values for extra top layer
1328  do i = 1, im
1329  plvl(i,llb) = prsmin
1330  if ( plvl(i,lla) <= prsmin ) plvl(i,lla) = 2.0*prsmin
1331  plyr(i,lyb) = 0.5 * plvl(i,lla)
1332  tlyr(i,lyb) = tlyr(i,lya)
1333 ! prslk1(i,lyb) = (plyr(i,lyb)*0.001) ** rocp ! plyr in hPa
1334  prslk1(i,lyb) = (plyr(i,lyb)*0.00001) ** rocp ! plyr in Pa
1335 
1336  rhly(i,lyb) = rhly(i,lya)
1337  qstl(i,lyb) = qstl(i,lya)
1338  enddo
1339 
1340  do j = 1, ntrac
1341  do i = 1, im
1342 ! --- note: may need to take care the top layer amount
1343  tracer1(i,lyb,j) = tracer1(i,lya,j)
1344  enddo
1345  enddo
1346  endif
1347 
1348 ! --- ... extra variables needed for ferrier's microphysics
1349 
1350  if (icmphys == 2) then
1351  do k = 1, lm
1352  k1 = k + kd
1353 
1354  do i = 1, im
1355  gcice(i,k1) = fcice(i,k)
1356  grain(i,k1) = frain(i,k)
1357  grime(i,k1) = rrime(i,k)
1358  enddo
1359  enddo
1360 
1361  if ( lextop ) then
1362  do i = 1, im
1363  gcice(i,lyb) = fcice(i,lya)
1364  grain(i,lyb) = frain(i,lya)
1365  grime(i,lyb) = rrime(i,lya)
1366  enddo
1367  endif
1368  endif ! if_icmphys
1369 
1371 ! --- ... get layer ozone mass mixing ratio
1372 
1373  if (ntoz > 0) then ! interactive ozone generation
1374 
1375  do k = 1, lmk
1376  do i = 1, im
1377  olyr(i,k) = max( qmin, tracer1(i,k,ntoz) )
1378  enddo
1379  enddo
1380 
1381  else ! climatological ozone
1382 
1383 ! print *,' in grrad : calling getozn'
1384  call getozn &
1385 ! --- inputs:
1386  & ( prslk1,xlat, &
1387  & im, lmk, &
1388 ! --- outputs:
1389  & olyr &
1390  & )
1391 
1392  endif ! end_if_ntoz
1394 ! --- ... compute cosin of zenith angle
1395 
1396  call coszmn &
1397 ! --- inputs:
1398  & ( xlon,sinlat,coslat,solhr, im, me, &
1399 ! --- outputs:
1400  & coszen, coszdg &
1401  & )
1412 
1413 ! --- ... set up non-prognostic gas volume mixing ratioes
1414 
1415  call getgases &
1416 ! --- inputs:
1417  & ( plvl, xlon, xlat, &
1418  & im, lmk, &
1419 ! --- outputs:
1420  & gasvmr &
1421  & )
1422 
1424 ! --- ... get temperature at layer interface, and layer moisture
1425 
1426  do k = 2, lmk
1427  do i = 1, im
1428  tem2da(i,k) = log( plyr(i,k) )
1429  tem2db(i,k) = log( plvl(i,k) )
1430  enddo
1431  enddo
1432 
1433  if (ivflip == 0) then ! input data from toa to sfc
1434 
1435  do i = 1, im
1436  tem1d(i) = qme6
1437  tem2da(i,1) = log( plyr(i,1) )
1438  tem2db(i,1) = 1.0
1439  tsfa(i) = tlyr(i,lmk) ! sfc layer air temp
1440  tlvl(i,1) = tlyr(i,1)
1441  tlvl(i,lmp) = tskn(i)
1442  enddo
1443 
1444  do k = 1, lm
1445  k1 = k + kd
1446 
1447  do i = 1, im
1448  qlyr(i,k1) = max( tem1d(i), qgrs(i,k) )
1449  tem1d(i) = min( qme5, qlyr(i,k1) )
1450  tvly(i,k1) = tgrs(i,k) * (1.0 + fvirt*qlyr(i,k1)) ! virtual T (K)
1451  enddo
1452  enddo
1453 
1454  if ( lextop ) then
1455  do i = 1, im
1456  qlyr(i,lyb) = qlyr(i,lya)
1457  tvly(i,lyb) = tvly(i,lya)
1458  enddo
1459  endif
1460 
1461  do k = 2, lmk
1462  do i = 1, im
1463  tlvl(i,k) = tlyr(i,k) + (tlyr(i,k-1) - tlyr(i,k)) &
1464  & * (tem2db(i,k) - tem2da(i,k)) &
1465  & / (tem2da(i,k-1) - tem2da(i,k))
1466  enddo
1467  enddo
1468 
1469  else ! input data from sfc to toa
1470 
1471  do i = 1, im
1472  tem1d(i) = qme6
1473  tem2da(i,1) = log( plyr(i,1) )
1474  tem2db(i,1) = log( plvl(i,1) )
1475  tsfa(i) = tlyr(i,1) ! sfc layer air temp
1476  tlvl(i,1) = tskn(i)
1477  tlvl(i,lmp) = tlyr(i,lmk)
1478  enddo
1479 
1480  do k = lm, 1, -1
1481  do i = 1, im
1482  qlyr(i,k) = max( tem1d(i), qgrs(i,k) )
1483  tem1d(i) = min( qme5, qlyr(i,k) )
1484  tvly(i,k) = tgrs(i,k) * (1.0 + fvirt*qlyr(i,k)) ! virtual T (K)
1485  enddo
1486  enddo
1487 
1488  if ( lextop ) then
1489  do i = 1, im
1490  qlyr(i,lyb) = qlyr(i,lya)
1491  tvly(i,lyb) = tvly(i,lya)
1492  enddo
1493  endif
1494 
1495  do k = 1, lmk-1
1496  do i = 1, im
1497  tlvl(i,k+1) = tlyr(i,k) + (tlyr(i,k+1) - tlyr(i,k)) &
1498  & * (tem2db(i,k+1) - tem2da(i,k)) &
1499  & / (tem2da(i,k+1) - tem2da(i,k))
1500  enddo
1501  enddo
1502 
1503  endif ! end_if_ivflip
1505 ! --- ... check for daytime points
1506 
1507  nday = 0
1508  do i = 1, im
1509  if (coszen(i) >= 0.0001) then
1510  nday = nday + 1
1511  idxday(nday) = i
1512  endif
1513  enddo
1514 
1515 ! write(0,*)' plvl=',plvl(ipt,1:65)
1516 ! write(0,*)' plyr=',plyr(ipt,1:64)
1517 ! write(0,*)' tlyr=',tlyr(ipt,1:64)
1518 ! write(0,*)' tlvl=',tlvl(ipt,1:65)
1519 ! write(0,*)' qlyr=',qlyr(ipt,1:10)*1000
1520 
1523 ! --- ... setup aerosols property profile for radiation
1524 
1525 !check print *,' in grrad : calling setaer '
1526 
1527  call setaer &
1528 ! --- inputs:
1529  & ( plvl,plyr,prslk1,tvly,rhly,slmsk,tracer1,xlon,xlat, &
1530  & im,lmk,lmp, lsswr,lslwr, &
1531 ! --- outputs:
1532  & faersw,faerlw,aerodp &
1533  & )
1534 
1543 ! --- ... obtain cloud information for radiation calculations
1544 
1545  if (ntcw > 0) then ! prognostic cloud scheme
1546 
1547  do k = 1, lmk
1548  do i = 1, im
1549  clw(i,k) = 0.0
1550  enddo
1551 
1552  do j = 1, ncld
1553  lv = ntcw + j - 1
1554  do i = 1, im
1555  clw(i,k) = clw(i,k) + tracer1(i,k,lv) ! cloud condensate amount
1556  enddo
1557  enddo
1558  enddo
1559 
1560  do k = 1, lmk
1561  do i = 1, im
1562  if ( clw(i,k) < epsq ) clw(i,k) = 0.0
1563  enddo
1564  enddo
1565 
1566  if (icmphys == 1) then ! zhao/moorthi's prognostic cloud scheme
1567 
1568  call progcld1 &
1569 ! --- inputs:
1570  & ( plyr,plvl,tlyr,tvly,qlyr,qstl,rhly,clw, &
1571  & xlat,xlon,slmsk, &
1572  & im, lmk, lmp, shoc_cld, cldcov(1:im,1:lm), &
1573 ! --- outputs:
1574  & clouds,cldsa,mtopa,mbota &
1575  & )
1576 
1577  elseif (icmphys == 2) then ! ferrier's microphysics
1578 
1579 ! print *,' in grrad : calling progcld2'
1580  call progcld2 &
1581 ! --- inputs:
1582  & ( plyr,plvl,tlyr,tvly,qlyr,qstl,rhly,clw, &
1583  & xlat,xlon,slmsk, gcice,grain,grime,flgmin, &
1584  & im, lmk, lmp, &
1585 ! --- outputs:
1586  & clouds,cldsa,mtopa,mbota &
1587  & )
1588 !
1589  elseif(icmphys == 3) then ! zhao/moorthi's prognostic cloud+pdfcld
1590 
1591  call progcld3 &
1592 ! --- inputs:
1593  & ( plyr,plvl,tlyr,tvly,qlyr,qstl,rhly,clw,cnvw,cnvc, &
1594  & xlat,xlon,slmsk, &
1595  & im, lmk, lmp, &
1596  & deltaq, sup,kdt,me, &
1597 ! --- outputs:
1598  & clouds,cldsa,mtopa,mbota &
1599  & )
1600 
1601  endif ! end if_icmphys
1602 
1603  else ! diagnostic cloud scheme
1604 
1605  do i = 1, im
1606 ! cvt1(i) = 10.0 * cvt(i)
1607 ! cvb1(i) = 10.0 * cvb(i)
1608  cvt1(i) = 0.01 * cvt(i)
1609  cvb1(i) = 0.01 * cvb(i)
1610 
1611  enddo
1612 
1613  do k = 1, lm
1614  k1 = k + kd
1615 
1616  do i = 1, im
1617 ! vvel(i,k1) = 10.0 * vvl(i,k)
1618  vvel(i,k1) = 0.01 * vvl(i,k)
1619  enddo
1620  enddo
1621 
1622  if ( lextop ) then
1623  do i = 1, im
1624  vvel(i,lyb) = vvel(i,lya)
1625  enddo
1626  endif
1627 
1628 ! --- compute diagnostic cloud related quantities
1629 
1630  call diagcld1 &
1631 ! --- inputs:
1632  & ( plyr,plvl,tlyr,rhly,vvel,cv,cvt1,cvb1, &
1633  & xlat,xlon,slmsk, &
1634  & im, lmk, lmp, &
1635 ! --- outputs:
1636  & clouds,cldsa,mtopa,mbota &
1637  & )
1638 
1639  endif ! end_if_ntcw
1640 
1641 ! --- ... start radiation calculations
1642 ! remember to set heating rate unit to k/sec!
1643 
1644  if (lsswr) then
1645 
1648 
1649 ! --- setup surface albedo for sw radiation, incl xw (nov04) sea-ice
1650 
1651  call setalb &
1652 ! --- inputs:
1653  & ( slmsk,snowd,sncovr,snoalb,zorl,coszen,tsfg,tsfa,hprim, &
1654  & alvsf,alnsf,alvwf,alnwf,facsf,facwf,fice,tisfc, &
1655  & im, &
1656 ! --- outputs:
1657  & sfcalb &
1658  & )
1659 
1661 
1662 ! --- approximate mean surface albedo from vis- and nir- diffuse values
1663 
1664  do i = 1, im
1665  sfalb(i) = max(0.01, 0.5 * (sfcalb(i,2) + sfcalb(i,4)))
1666  enddo
1667 
1668  if (nday > 0) then
1669 
1671 ! print *,' in grrad : calling swrad'
1672 
1673  if ( present(htrswb) .and. present(htrsw0)) then
1674 
1675  call swrad &
1676 ! --- inputs:
1677  & ( plyr,plvl,tlyr,tlvl,qlyr,olyr,gasvmr, &
1678  & clouds,icsdsw,faersw,sfcalb, &
1679  & coszen,solcon, nday,idxday, &
1680  & im, lmk, lmp, lprnt, &
1681 ! --- outputs:
1682  & htswc,topfsw,sfcfsw &
1683 !! --- optional:
1684 !! &, HSW0=htsw0,FLXPRF=fswprf &
1685  &, hsw0=htsw0,hswb=htswb,fdncmp=scmpsw &
1686  & )
1687 
1688  do k = 1, lm
1689  k1 = k + kd
1690 
1691  do j = 1, nbdsw
1692  do i = 1, im
1693  htrswb(i,k,j) = htswb(i,k1,j)
1694  enddo
1695  enddo
1696  enddo
1697 
1698  else if ( present(htrswb) .and. .not. present(htrsw0) ) then
1699 
1700  call swrad &
1701 ! --- inputs:
1702  & ( plyr,plvl,tlyr,tlvl,qlyr,olyr,gasvmr, &
1703  & clouds,icsdsw,faersw,sfcalb, &
1704  & coszen,solcon, nday,idxday, &
1705  & im, lmk, lmp, lprnt, &
1706 ! --- outputs:
1707  & htswc,topfsw,sfcfsw &
1708 !! --- optional:
1709 !! &, hsw0=htsw0,flxprf=fswprf &
1710  &, hswb=htswb,fdncmp=scmpsw &
1711  & )
1712 
1713  do k = 1, lm
1714  k1 = k + kd
1715  do j = 1, nbdsw
1716  do i = 1, im
1717  htrswb(i,k,j) = htswb(i,k1,j)
1718  enddo
1719  enddo
1720  enddo
1721 
1722  else if ( present(htrsw0) .and. .not. present(htrswb) ) then
1723 
1724  call swrad &
1725 ! --- inputs:
1726  & ( plyr,plvl,tlyr,tlvl,qlyr,olyr,gasvmr, &
1727  & clouds,icsdsw,faersw,sfcalb, &
1728  & coszen,solcon, nday,idxday, &
1729  & im, lmk, lmp, lprnt, &
1730 ! --- outputs:
1731  & htswc,topfsw,sfcfsw &
1732 !! --- optional:
1733 !! &, hsw0=htsw0,flxprf=fswprf &
1734  &, hsw0=htsw0,fdncmp=scmpsw &
1735  & )
1736 
1737  else
1738 
1739  call swrad &
1740 ! --- inputs:
1741  & ( plyr,plvl,tlyr,tlvl,qlyr,olyr,gasvmr, &
1742  & clouds,icsdsw,faersw,sfcalb, &
1743  & coszen,solcon, nday,idxday, &
1744  & im, lmk, lmp, lprnt, &
1745 ! --- outputs:
1746  & htswc,topfsw,sfcfsw &
1747 !! --- optional:
1748 !! &, HSW0=htsw0,FLXPRF=fswprf,HSWB=htswb &
1749  &, fdncmp=scmpsw &
1750  & )
1751 
1752  endif
1753 
1754  do k = 1, lm
1755  k1 = k + kd
1756 
1757  do i = 1, im
1758  htrsw(i,k) = htswc(i,k1)
1759  enddo
1760  enddo
1761  if (present(htrsw0)) then
1762  do k = 1, lm
1763  k1 = k + kd
1764  do i = 1, im
1765  htrsw0(i,k) = htsw0(i,k1)
1766  enddo
1767  enddo
1768  endif
1769 
1770 ! --- surface down and up spectral component fluxes
1771 
1772  do i = 1, im
1773  dswcmp(i,1) = scmpsw(i)%nirbm
1774  dswcmp(i,2) = scmpsw(i)%nirdf
1775  dswcmp(i,3) = scmpsw(i)%visbm
1776  dswcmp(i,4) = scmpsw(i)%visdf
1777 
1778  uswcmp(i,1) = scmpsw(i)%nirbm * sfcalb(i,1)
1779  uswcmp(i,2) = scmpsw(i)%nirdf * sfcalb(i,2)
1780  uswcmp(i,3) = scmpsw(i)%visbm * sfcalb(i,3)
1781  uswcmp(i,4) = scmpsw(i)%visdf * sfcalb(i,4)
1782  enddo
1783 
1784  else ! if_nday_block
1785 
1786  do k = 1, lm
1787  do i = 1, im
1788  htrsw(i,k) = 0.0
1789  enddo
1790  enddo
1791 
1792  sfcfsw = sfcfsw_type( 0.0, 0.0, 0.0, 0.0 )
1793  topfsw = topfsw_type( 0.0, 0.0, 0.0 )
1794  scmpsw = cmpfsw_type( 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 )
1795 
1796  do k = 1, 4
1797  do i = 1, im
1798  dswcmp(i,k) = 0.0
1799  uswcmp(i,k) = 0.0
1800  enddo
1801  enddo
1802 
1803 !! --- optional:
1804 !! fswprf= profsw_type( 0.0, 0.0, 0.0, 0.0 )
1805 
1806  if ( present(htrswb) ) then
1807  do j = 1, nbdsw
1808  do k = 1, lm
1809  do i = 1, im
1810  htrswb(i,k,j) = 0.0
1811  enddo
1812  enddo
1813  enddo
1814  endif
1815  if ( present(htrsw0) ) then
1816  do k = 1, lm
1817  do i = 1, im
1818  htrsw0(i,k) = 0.0
1819  enddo
1820  enddo
1821  endif
1822 
1823  endif ! end_if_nday
1824 
1825  endif ! end_if_lsswr
1826 
1827 ! write(0,*)' htrsw=',htrsw(ipt,1:64)*86400
1828  if (lslwr) then
1829 
1832 ! --- setup surface emissivity for lw radiation
1833 
1834  call setemis &
1835 ! --- inputs:
1836  & ( xlon,xlat,slmsk,snowd,sncovr,zorl,tsfg,tsfa,hprim, &
1837  & im, &
1838 ! --- outputs:
1839  & sfcemis &
1840  & )
1842 ! print *,' in grrad : calling lwrad'
1843 
1844  if ( present(htrlwb) .and. present(htrlw0) ) then
1845 
1846  call lwrad &
1847 ! --- inputs:
1848  & ( plyr,plvl,tlyr,tlvl,qlyr,olyr,gasvmr, &
1849  & clouds,icsdlw,faerlw,sfcemis,tsfg, &
1850  & im, lmk, lmp, lprnt, &
1851 ! --- outputs:
1852  & htlwc,topflw,sfcflw &
1853 !! --- optional:
1854 !! &, HLW0=htlw0,FLXPRF=flwprf &
1855  &, hlw0=htlw0 &
1856  &, hlwb=htlwb &
1857  & )
1858 
1859  do k = 1, lm
1860  k1 = k + kd
1861 
1862  do j = 1, nbdlw
1863  do i = 1, im
1864  htrlwb(i,k,j) = htlwb(i,k1,j)
1865  enddo
1866  enddo
1867  enddo
1868 
1869  else if ( present(htrlwb) .and. .not. present(htrlw0) ) then
1870 
1871  call lwrad &
1872 ! --- inputs:
1873  & ( plyr,plvl,tlyr,tlvl,qlyr,olyr,gasvmr, &
1874  & clouds,icsdlw,faerlw,sfcemis,tsfg, &
1875  & im, lmk, lmp, lprnt, &
1876 ! --- outputs:
1877  & htlwc,topflw,sfcflw &
1878 !! --- optional:
1879 !! &, hlw0=htlw0,flxprf=flwprf &
1880  &, hlwb=htlwb &
1881  & )
1882 
1883  do k = 1, lm
1884  k1 = k + kd
1885 
1886  do j = 1, nbdlw
1887  do i = 1, im
1888  htrlwb(i,k,j) = htlwb(i,k1,j)
1889  enddo
1890  enddo
1891  enddo
1892  else if ( present(htrlw0) .and. .not. present(htrlwb) ) then
1893 
1894  !print *,'call lwrad saving clear sky component'
1895  call lwrad &
1896 ! --- inputs:
1897  & ( plyr,plvl,tlyr,tlvl,qlyr,olyr,gasvmr, &
1898  & clouds,icsdlw,faerlw,sfcemis,tsfg, &
1899  & im, lmk, lmp, lprnt, &
1900 ! --- outputs:
1901  & htlwc,topflw,sfcflw &
1902 !! --- optional:
1903 !! &, hlw0=htlw0,flxprf=flwprf &
1904  &, hlw0=htlw0 &
1905  & )
1906 
1907  else
1908 
1909  call lwrad &
1910 ! --- inputs:
1911  & ( plyr,plvl,tlyr,tlvl,qlyr,olyr,gasvmr, &
1912  & clouds,icsdlw,faerlw,sfcemis,tsfg, &
1913  & im, lmk, lmp, lprnt, &
1914 ! --- outputs:
1915  & htlwc,topflw,sfcflw &
1916 !! --- optional:
1917 !! &, HLW0=htlw0,FLXPRF=flwprf,HLWB=htlwb &
1918  & )
1919 
1920  endif
1921 
1922  do i = 1, im
1923  semis(i) = sfcemis(i)
1924 ! --- save surface air temp for diurnal adjustment at model t-steps
1925  tsflw(i) = tsfa(i)
1926  enddo
1927 
1928  do k = 1, lm
1929  k1 = k + kd
1930 
1931  do i = 1, im
1932  htrlw(i,k) = htlwc(i,k1)
1933  enddo
1934  enddo
1935 
1936  if (present(htrlw0)) then
1937  do k = 1, lm
1938  k1 = k + kd
1939  do i = 1, im
1940  htrlw0(i,k) = htlw0(i,k1)
1941  enddo
1942  enddo
1943  endif
1944 
1945  endif ! end_if_lslwr
1946 
1948 ! --- ... collect the fluxr data for wrtsfc
1949 
1950  if (lssav) then
1951 
1952  if ( lsswr ) then
1953  do i = 1, im
1954  fluxr(i,34) = fluxr(i,34) + dtsw*aerodp(i,1) ! total aod at 550nm
1955  fluxr(i,35) = fluxr(i,35) + dtsw*aerodp(i,2) ! DU aod at 550nm
1956  fluxr(i,36) = fluxr(i,36) + dtsw*aerodp(i,3) ! BC aod at 550nm
1957  fluxr(i,37) = fluxr(i,37) + dtsw*aerodp(i,4) ! OC aod at 550nm
1958  fluxr(i,38) = fluxr(i,38) + dtsw*aerodp(i,5) ! SU aod at 550nm
1959  fluxr(i,39) = fluxr(i,39) + dtsw*aerodp(i,6) ! SS aod at 550nm
1960  enddo
1961  endif
1962 ! --- save lw toa and sfc fluxes
1963 
1964  if (lslwr) then
1965  do i = 1, im
1966 ! --- lw total-sky fluxes
1967  fluxr(i,1 ) = fluxr(i,1 ) + dtlw * topflw(i)%upfxc ! total sky top lw up
1968  fluxr(i,19) = fluxr(i,19) + dtlw * sfcflw(i)%dnfxc ! total sky sfc lw dn
1969  fluxr(i,20) = fluxr(i,20) + dtlw * sfcflw(i)%upfxc ! total sky sfc lw up
1970 ! --- lw clear-sky fluxes
1971  fluxr(i,28) = fluxr(i,28) + dtlw * topflw(i)%upfx0 ! clear sky top lw up
1972  fluxr(i,30) = fluxr(i,30) + dtlw * sfcflw(i)%dnfx0 ! clear sky sfc lw dn
1973  fluxr(i,33) = fluxr(i,33) + dtlw * sfcflw(i)%upfx0 ! clear sky sfc lw up
1974  enddo
1975  endif
1976 
1977 ! --- save sw toa and sfc fluxes with proper diurnal sw wgt. coszen=mean cosz over daylight
1978 ! part of sw calling interval, while coszdg= mean cosz over entire interval
1979 
1980  if (lsswr) then
1981  do i = 1, im
1982  if (coszen(i) > 0.) then
1983 ! --- sw total-sky fluxes
1984 ! -------------------
1985  tem0d = dtsw * coszdg(i) / coszen(i)
1986  fluxr(i,2 ) = fluxr(i,2) + topfsw(i)%upfxc * tem0d ! total sky top sw up
1987  fluxr(i,3 ) = fluxr(i,3) + sfcfsw(i)%upfxc * tem0d ! total sky sfc sw up
1988  fluxr(i,4 ) = fluxr(i,4) + sfcfsw(i)%dnfxc * tem0d ! total sky sfc sw dn
1989 ! --- sw uv-b fluxes
1990 ! --------------
1991  fluxr(i,21) = fluxr(i,21) + scmpsw(i)%uvbfc * tem0d ! total sky uv-b sw dn
1992  fluxr(i,22) = fluxr(i,22) + scmpsw(i)%uvbf0 * tem0d ! clear sky uv-b sw dn
1993 ! --- sw toa incoming fluxes
1994 ! ----------------------
1995  fluxr(i,23) = fluxr(i,23) + topfsw(i)%dnfxc * tem0d ! top sw dn
1996 ! --- sw sfc flux components
1997 ! ----------------------
1998  fluxr(i,24) = fluxr(i,24) + scmpsw(i)%visbm * tem0d ! uv/vis beam sw dn
1999  fluxr(i,25) = fluxr(i,25) + scmpsw(i)%visdf * tem0d ! uv/vis diff sw dn
2000  fluxr(i,26) = fluxr(i,26) + scmpsw(i)%nirbm * tem0d ! nir beam sw dn
2001  fluxr(i,27) = fluxr(i,27) + scmpsw(i)%nirdf * tem0d ! nir diff sw dn
2002 ! --- sw clear-sky fluxes
2003 ! -------------------
2004  fluxr(i,29) = fluxr(i,29) + topfsw(i)%upfx0 * tem0d ! clear sky top sw up
2005  fluxr(i,31) = fluxr(i,31) + sfcfsw(i)%upfx0 * tem0d ! clear sky sfc sw up
2006  fluxr(i,32) = fluxr(i,32) + sfcfsw(i)%dnfx0 * tem0d ! clear sky sfc sw dn
2007  endif
2008  enddo
2009  endif
2010 
2011 ! --- save total and boundary layer clouds
2012 
2013  if (lsswr .or. lslwr) then
2014  do i = 1, im
2015  fluxr(i,17) = fluxr(i,17) + raddt * cldsa(i,4)
2016  fluxr(i,18) = fluxr(i,18) + raddt * cldsa(i,5)
2017  enddo
2018 
2019 ! --- save cld frac,toplyr,botlyr and top temp, note that the order
2020 ! of h,m,l cloud is reversed for the fluxr output.
2021 ! --- save interface pressure (pa) of top/bot
2022 
2023  do j = 1, 3
2024  do i = 1, im
2025  tem0d = raddt * cldsa(i,j)
2026  itop = mtopa(i,j) - kd
2027  ibtc = mbota(i,j) - kd
2028  fluxr(i, 8-j) = fluxr(i, 8-j) + tem0d
2029  fluxr(i,11-j) = fluxr(i,11-j) + tem0d * prsi(i,itop+kt)
2030  fluxr(i,14-j) = fluxr(i,14-j) + tem0d * prsi(i,ibtc+kb)
2031  fluxr(i,17-j) = fluxr(i,17-j) + tem0d * tgrs(i,itop)
2032  enddo
2033  enddo
2034  endif
2035 
2036  if (.not. shoc_cld) then
2037  do k = 1, lm
2038  k1 = k + kd
2039  do i = 1, im
2040  cldcov(i,k) = clouds(i,k1,1)
2041  enddo
2042  enddo
2043  endif
2044 
2045 ! --- save optional vertically integrated aerosol optical depth at
2046 ! wavelenth of 550nm aerodp(:,1), and other optional aod for
2047 ! individual species aerodp(:,2:NSPC1)
2048 
2049 ! if ( laswflg ) then
2050 ! if ( NFXR > 33 ) then
2051 ! do i = 1, IM
2052 ! fluxr(i,34) = fluxr(i,34) + dtsw*aerodp(i,1) ! total aod at 550nm (all species)
2053 ! enddo
2054 
2055 ! if ( lspcodp ) then
2056 ! do j = 2, NSPC1
2057 ! k = 33 + j
2058 
2059 ! do i = 1, IM
2060 ! fluxr(i,k) = fluxr(i,k) + dtsw*aerodp(i,j) ! aod at 550nm for indiv species
2061 ! enddo
2062 ! enddo
2063 ! endif ! end_if_lspcodp
2064 ! else
2065 ! print *,' !Error! Need to increase array fluxr size NFXR ',&
2066 ! & ' to be able to output aerosol optical depth'
2067 ! stop
2068 ! endif ! end_if_nfxr
2069 ! endif ! end_if_laswflg
2070 
2071  endif ! end_if_lssav
2072 !
2073  return
2074 !...................................
2075  end subroutine grrad
2076 !-----------------------------------
2078 
2079 !
2080 !........................................!
2081  end module module_radiation_driver !
2082 !========================================!
subroutine, public progcld1
This subroutine computes cloud related quantities using zhao/moorthi&#39;s prognostic cloud microphysics ...
integer, save ialbflg
surface albedo scheme control flag
Definition: physparam.f:213
real(kind=kind_phys) qmin
Definition: grrad.f:225
subroutine, public gas_init
This subroutine sets up ozone, co2, etc. parameters. if climatology ozone then read in monthly ozone ...
subroutine, public swrad
This subroutine is the main sw radiation routine.
Definition: radsw_main.f:440
integer, save ioznflg
ozone dta source control flag
Definition: physparam.f:175
integer, save iaerflg
aerosol effect control flag
Definition: physparam.f:155
define type construct for optional component downward fluxes at surface
Definition: radsw_param.f:110
integer, save isubcsw
sub-column cloud approx flag in sw radiation
Definition: physparam.f:226
integer, save iemsflg
surface emissivity scheme control flag
Definition: physparam.f:215
subroutine, public coszmn
This subroutine computes mean cos solar zenith angle over SW calling interval.
subroutine, public sol_update
This subroutine computes solar parameters at forecast time.
define type construct for optional radiation flux profiles
Definition: radlw_param.f:94
define type construct for radiation fluxes at toa
Definition: radsw_param.f:70
integer, parameter, public nf_vgas
number of gas species
This module computes cloud related quantities for radiation computations.
integer, save iovrsw
cloud overlapping control flag for sw
Definition: physparam.f:197
subroutine, public gas_update
This subroutine reads in 2-d monthly co2 data set for a specified year. data are in a 15 degree lat/l...
subroutine, public sfc_init
This subroutine is the initialization program for surface radiation related quantities (albedo...
This module contains some the most frequently used math and physics constants for gcm models...
Definition: physcons.f:29
define type construct for radiation fluxes at surface
Definition: radsw_param.f:82
subroutine, public radupdate
This subroutine calls many update subroutines to check and update radiation required but time varying...
Definition: grrad.f:495
integer, parameter nbdlw
Definition: radlw_param.f:122
integer, parameter nbdsw
Definition: radsw_param.f:144
subroutine, public grrad
This subroutine is the driver of radiation calculation subroutines. It sets up profile variables for ...
Definition: grrad.f:811
integer, parameter, public nf_aelw
num of output fields for lw rad
integer, save iovrlw
cloud overlapping control flag for lw
Definition: physparam.f:199
subroutine, public getgases
This subroutine sets up global distribution of radiation absorbing gases in volume mixing ratio...
This module defines commonly used control variables/parameters in physics related programs...
Definition: physparam.f:27
integer, save isolar
solar constant scheme control flag
Definition: physparam.f:143
This module contains SW band parameters set up.
Definition: radsw_param.f:58
character(40), parameter vtagrad
Definition: grrad.f:219
integer, save ico2flg
co2 data source control flag
Definition: physparam.f:171
integer, save isubclw
sub-column cloud approx flag in lw radiation
Definition: physparam.f:228
real(kind=kind_phys) qme5
Definition: grrad.f:225
integer, save icmphys
cloud microphysics scheme control flag
Definition: physparam.f:195
logical, save lnoprec
precip effect on radiation flag (ferrier microphysics)
Definition: physparam.f:205
This module contains climatological atmospheric aerosol schemes for radiation computations.
real(kind=kind_phys), parameter con_rocp
Definition: physcons.f:104
logical, save lcrick
eliminating CRICK control flag
Definition: physparam.f:201
integer, save ivflip
vertical profile indexing flag
Definition: physparam.f:224
real(kind=kind_phys) epsq
Definition: grrad.f:225
define type construct for optional radiation flux profiles
Definition: radsw_param.f:96
define type construct for radiation fluxes at surface
Definition: radlw_param.f:80
subroutine, public sol_init
This subroutine initializes astronomy process, and set up module constants.
This module sets up astronomy quantities for solar radiation calculations.
subroutine, public setemis
This subroutine computes surface emissivity for LW radiation.
subroutine, public rswinit
This subroutine initializes non-varying module variables, conversion factors, and look-up tables...
Definition: radsw_main.f:1272
subroutine, public aer_init
The initialization program to set up necessary parameters and working arrays.
This module contains LW band parameters set up.
Definition: radlw_param.f:58
This module includes ncep&#39;s modifications of the rrtm-lw radiation ! code from aer inc...
Definition: radlw_main.f:236
subroutine, public setalb
This subroutine computes four components of surface albedos (i.e., vis-nir, direct-diffused) accordin...
logical, parameter lextop
Definition: grrad.f:241
This module includes ncep&#39;s modifications of the rrtm-sw radiation code from aer inc.
Definition: radsw_main.f:260
integer, parameter, public nf_albd
num of sfc albedo components
subroutine, public lwrad
This subroutine is the main lw radiation routine.
Definition: radlw_main.f:415
logical, save lsashal
shallow convection flag
Definition: physparam.f:207
define type construct for radiation fluxes at toa
Definition: radlw_param.f:70
integer, parameter, public nf_clds
number of fields in cloud array
real(kind=kind_phys) qme6
Definition: grrad.f:225
subroutine, public getozn
This subroutine sets up climatological ozone profile for radiation calculation this code is originall...
integer, parameter, public nf_aesw
num of output fields for sw rad
subroutine, public cld_init
This subroutine is an initialization program for cloud-radiation calculations. It sets up boundary la...
This module sets up ozone climatological profiles and other constant gas profiles, such as co2, ch4, n2o, o2, and those of cfc gases. All data are entered as mixing ratio by volume, except ozone which is mass mixing ratio (g/g).
subroutine, public diagcld1
This subroutine computes cloud fractions for radiation calculations.
real(kind=kind_phys), parameter con_fvirt
Definition: physcons.f:107
subroutine, public radinit
This subroutine is the initialization of radiation calculations.
Definition: grrad.f:261
subroutine, public progcld2
This subroutine computes cloud related quantities using ferrier&#39;s prognostic cloud microphysics schem...
subroutine, public progcld3
This subroutine computes cloud related quantities using zhao/moorthi&#39;s prognostic cloud microphysics ...
subroutine, public rlwinit
This subroutine performs calculations necessary for the initialization of the longwave model...
Definition: radlw_main.f:1249
real(kind=kind_phys), parameter con_epsm1
Definition: physcons.f:109
integer, parameter, public nspc1
total+species
integer, parameter ltp
Definition: grrad.f:239
real(kind=kind_phys), parameter con_eps
Definition: physcons.f:108
This module sets up surface albedo for sw radiation and surface emissivity for lw radiation...
logical, save lcnorm
in-cld condensate control flag
Definition: physparam.f:203
integer, save ictmflg
external data time/date control flag
Definition: physparam.f:173
subroutine, public setaer
This subroutine computes aerosols optical properties.
real, parameter prsmin
Definition: grrad.f:229
This file is the radiation driver module. it prepares atmospheric profiles and invokes main radiation...
Definition: grrad.f:184
integer, save icldflg
cloud optical property scheme control flag
Definition: physparam.f:193
subroutine, public aer_update
This subroutine checks and updates time varying climatology aerosol data sets.