80 use module_iounitdef
, only : niradsf
87 character(40),
parameter :: &
88 & VTAGAST=
'NCEP-Radiation_astronomy v5.2 Jan 2013 ' 93 real (kind=kind_phys),
parameter ::
tpi = 2.0 *
con_pi 94 real (kind=kind_phys),
parameter ::
hpi = 0.5 *
con_pi 95 real (kind=kind_phys),
parameter ::
f12 = 12.0
96 real (kind=kind_phys),
parameter ::
f3600 = 3600.0
97 real (kind=kind_phys),
parameter ::
czlimt = 0.0001
99 real (kind=kind_phys),
parameter ::
pid12 = (2.0*asin(1.0))/
f12 177 integer,
intent(in) :: me
182 logical :: file_exist
186 if ( me == 0 ) print *,
vtagast 198 print *,
' - Using old fixed solar constant =',
solc0 200 elseif (
isolar == 10 )
then 202 print *,
' - Using new fixed solar constant =',
solc0 204 elseif (
isolar == 1 )
then 208 print *,
' - Using NOAA annual mean TSI table in ABS scale', &
209 &
' with cycle approximation (old values)!' 213 if ( .not. file_exist )
then 217 print *,
' Requested solar data file "',
solar_fname, &
219 print *,
' Using the default solar constant value =',
solc0,&
220 &
' reset control flag isolflg=',
isolflg 223 elseif (
isolar == 2 )
then 227 print *,
' - Using NOAA annual mean TSI table in TIM scale', &
228 &
' with cycle approximation (new values)!' 232 if ( .not. file_exist )
then 236 print *,
' Requested solar data file "',
solar_fname, &
238 print *,
' Using the default solar constant value =',
solc0,&
239 &
' reset control flag isolflg=',
isolflg 242 elseif (
isolar == 3 )
then 246 print *,
' - Using CMIP5 annual mean TSI table in TIM scale', &
247 &
' with cycle approximation' 251 if ( .not. file_exist )
then 255 print *,
' Requested solar data file "',
solar_fname, &
257 print *,
' Using the default solar constant value =',
solc0,&
258 &
' reset control flag isolflg=',
isolflg 261 elseif (
isolar == 4 )
then 265 print *,
' - Using CMIP5 monthly mean TSI table in TIM scale', &
266 &
' with cycle approximation' 270 if ( .not. file_exist )
then 274 print *,
' Requested solar data file "',
solar_fname, &
276 print *,
' Using the default solar constant value =',
solc0,&
277 &
' reset control flag isolflg=',
isolflg 284 print *,
' - !!! ERROR in selection of solar constant data', &
285 &
' source, ISOL =',
isolar 286 print *,
' Using the default solar constant value =',
solc0, &
287 &
' reset control flag isolflg=',
isolflg 311 & ( jdate,kyear,deltsw,deltim,lsol_chg, me, &
313 & slag, sdec, cdec, solcon &
365 integer,
intent(in) :: jdate(:), kyear, me
366 logical,
intent(in) :: lsol_chg
368 real (kind=kind_phys),
intent(in) :: deltsw, deltim
371 real (kind=kind_phys),
intent(out) :: slag, sdec, cdec, solcon
374 real (kind=kind_phys),
parameter :: hrday = 1.0/24.0
375 real (kind=kind_phys),
parameter :: minday= 1.0/1440.0
376 real (kind=kind_phys),
parameter :: secday= 1.0/86400.0
378 real (kind=kind_phys) :: smean, solc1, dtswh, smon(12)
379 real (kind=kind_phys) :: fjd, fjd1, dlt, r1, alp
381 integer :: jd, jd1, iyear, imon, iday, ihr, imin, isec
383 integer :: i, iyr, iyr1, iyr2, jyr, nn, nswr, icy1, icy2, icy
385 logical :: file_exist
386 character :: cline*60
410 if ( .not. file_exist )
then 411 print *,
' !!! ERROR! Can not find solar constant file!!!' 421 read (niradsf, * ) iyr1,iyr2,icy1,icy2,smean,cline(1:60)
426 print *,
' Updating solar constant with cycle approx' 427 print *,
' Opened solar constant data file: ',
solar_fname 466 if ( iyr < iyr1 )
then 467 icy = icy1 - iyr1 + 1
468 lab_dowhile1 :
do while ( iyr < iyr1 )
473 print *,
' *** Year',iyear,
' out of table range!', &
475 print *,
' Using the closest-cycle year (',iyr,
')' 477 elseif ( iyr > iyr2 )
then 478 icy = iyr2 - icy2 + 1
479 lab_dowhile2 :
do while ( iyr > iyr2 )
484 print *,
' *** Year',iyear,
' out of table range!', &
486 print *,
' Using the closest-cycle year (',iyr,
')' 494 lab_dowhile3 :
do while ( i >= iyr1 )
497 read (niradsf,*) jyr, solc1
499 if ( i == iyr .and. iyr == jyr )
then 500 solc0 = smean + solc1
503 print *,
' CHECK: Solar constant data used for year',&
514 lab_dowhile4 :
do while ( i >= iyr1 )
517 read (niradsf,*) jyr, smon(1:12)
519 if ( i == iyr .and. iyr == jyr )
then 523 solc0 = smean + smon(imon)
526 print *,
' CHECK: Solar constant data used for year',&
527 & iyr,
' and month',imon
545 jd1 = iw3jdn(iyear,imon,iday)
552 fjd1= 0.5 + float(ihr)*hrday + float(imin)*minday &
553 & + float(isec)*secday
555 fjd1= float(ihr - 12)*hrday + float(imin)*minday &
556 & + float(isec)*secday
572 solcon =
solc0 / (r1*r1)
584 & ( jd, fjd, dlt, alp, r1, solcon )
591 nswr = nint(deltsw / deltim)
592 dtswh = deltsw /
f3600 594 if ( deltsw >=
f3600 )
then 595 nn = max(6, min(12, nint(
f3600/deltim) ))
596 nstp = nint(dtswh) * nn + 1
598 nstp = max(2, min(20, nswr)) + 1
605 print *,
' for cosz calculations: nswr,deltim,deltsw,dtswh =', &
606 & nswr,deltim,deltsw,dtswh,
' anginc,nstp =',
anginc,
nstp 654 real (kind=kind_phys),
intent(in) :: fjd
655 integer,
intent(in) :: jd
658 real (kind=kind_phys),
intent(out) :: r1, dlt, alp
661 real (kind=kind_phys),
parameter :: cyear = 365.25
662 real (kind=kind_phys),
parameter :: ccr = 1.3e-6
663 real (kind=kind_phys),
parameter :: tpp = 1.55
665 real (kind=kind_phys),
parameter :: svt6 = 78.035
667 integer,
parameter :: jdor = 2415020
670 real (kind=kind_phys) :: dat, t1, year, tyear, ec, angin, ador, &
671 & deleqn, sni, tini, er, qq, e1, ep, cd, eq, date, em, &
674 integer :: jdoe, iter
680 t1 = float(jd - jdor) / 36525.0
684 year = 0.25964134e0 + 0.304e-5 * t1
685 tyear= 0.24219879e0 - 0.614e-5 * t1
689 ec = 0.01675104e0 - (0.418e-4 + 0.126e-6 * t1) * t1
690 angin= 23.452294e0 - (0.0130125e0 + 0.164e-5 * t1) * t1
693 jdoe = ador + (svt6 * cyear) / (year - tyear)
697 deleqn= float(jdoe - jd) * (year - tyear) / cyear
699 sni = sin( angin /
degrad )
700 tini = 1.0 / tan( angin /
degrad )
701 er = sqrt( (1.0 + ec) / (1.0 - ec) )
702 qq = deleqn *
tpi / year
710 lab_do_1 :
do while ( cd > ccr )
712 ep = e1 - (e1 - ec*sin(e1) - qq) / (1.0 - ec*cos(e1))
718 write(6,*)
' ITERATION COUNT FOR LOOP 32 =', iter
719 write(6,*)
' E, EP, CD =', e1, ep, cd
725 eq = 2.0 * atan( er * tan( 0.5*e1 ) )
729 dat = float(jd - jdor) - tpp + fjd
730 date = mod(dat, year)
734 em =
tpi * date / year
739 lab_do_2 :
do while ( cr > ccr )
741 ep = e1 - (e1 - ec*sin(e1) - em) / (1.0 - ec*cos(e1))
747 write(6,*)
' ITERATION COUNT FOR LOOP 31 =', iter
753 w1 = 2.0 * atan( er * tan( 0.5*e1 ) )
755 r1 = 1.0 - ec*cos(e1)
757 sindec = sni * sin(w1 - eq)
761 alp = asin( tan(dlt)*tini )
764 if (tst < 0.0) alp =
con_pi - alp
765 if (alp < 0.0) alp = alp +
tpi 767 sun =
tpi * (date - deleqn) / year
768 if (sun < 0.0) sun = sun +
tpi 769 sollag = sun - alp - 0.03255e0
789 & ( xlon,sinlat,coslat,solhr, im, me, &
827 integer,
intent(in) :: IM, me
829 real (kind=kind_phys),
intent(in) :: sinlat(:), coslat(:), &
833 real (kind=kind_phys),
intent(out) :: coszen(:), coszdg(:)
836 real (kind=kind_phys) :: coszn, cns, ss, cc, solang, rstp
838 integer :: istsun(im), i, it, j, lat
843 rstp = 1.0 / float(
nstp)
857 coszn = ss + cc * cos(cns + xlon(i))
858 coszen(i) = coszen(i) + max(0.0, coszn)
859 if (coszn >
czlimt) istsun(i) = istsun(i) + 1
866 coszdg(i) = coszen(i) * rstp
867 if (istsun(i) > 0) coszen(i) = coszen(i) / istsun(i)
880 & ( jd, fjd, dlt, alp, r1, solc &
910 integer,
intent(in) :: jd
912 real (kind=kind_phys),
intent(in) :: fjd, dlt, alp, r1, solc
917 real (kind=kind_phys),
parameter :: sixty = 60.0
919 character(LEN=1),
parameter :: sign =
'-' 920 character(LEN=1),
parameter :: sigb =
' ' 922 character(LEN=1) :: dsig
923 character(LEN=4) :: month(12)
925 data month /
'JAN.',
'FEB.',
'MAR.',
'APR.',
'MAY ',
'JUNE', &
926 &
'JULY',
'AUG.',
'SEP.',
'OCT.',
'NOV ',
'DEC.' /
928 integer :: iday, imon, iyear, ihr, ltd, ltm, &
929 & ihalp, iyy, jda, mfjd, idaywk, idayyr
930 real (kind=kind_phys) :: xmin, dltd, dltm, dlts, halp, ymin, &
939 mfjd= nint( fjd*1440.0 )
941 xmin= float(mfjd) - (ihr + 12)*sixty
944 mfjd= nint( fjd*1440.0 )
946 xmin= float(mfjd) - (ihr - 12)*sixty
951 call w3fs26(jda, iyear,imon,iday, idaywk,idayyr)
957 dltm = sixty * (abs(dltd) - abs(float(ltd)))
959 dlts = sixty * (dltm - float(ltm))
961 if ((dltd < 0.0) .and. (ltd == 0.0))
then 967 halp = 6.0 * alp /
hpi 969 ymin = abs(halp - float(ihalp)) * sixty
971 asec = (ymin - float(iyy)) * sixty
976 print 101, iday, month(imon), iyear, ihr, xmin, jd, fjd
977 101
format(
'0 FORECAST DATE',9x,i3,a5,i6,
' AT',i3,
' HRS',f6.2,
' MINS'/&
978 &
' JULIAN DAY',12x,i8,2x,
'PLUS',f11.6)
980 print 102, r1, halp, ihalp, iyy, asec
981 102
format(
' RADIUS VECTOR',9x,f10.7/
' RIGHT ASCENSION OF SUN', &
982 &
f12.7,
' HRS, OR',i4,
' HRS',i4,
' MINS',f6.1,
' SECS')
984 print 103, dltd, dsig, ltd, ltm, dlts, eqt, eqsec,
sollag, solc
985 103
format(
' DECLINATION OF THE SUN',
f12.7,
' DEGS, OR ',a1,i3, &
986 &
' DEGS',i4,
' MINS',f6.1,
' SECS'/
' EQUATION OF TIME',6x, &
987 &
f12.7,
' MINS, OR',f10.2,
' SECS, OR',f9.6,
' RADIANS'/ &
988 &
' SOLAR CONSTANT',8x,
f12.7,
' (DISTANCE AJUSTED)'//)
real(kind=kind_phys) anginc
real(kind=kind_phys), parameter czlimt
real(kind=kind_phys) sollag
real(kind=kind_phys), parameter tpi
real(kind=kind_phys) cosdec
real(kind=kind_phys), parameter degrad
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.
real(kind=kind_phys), parameter con_solr_old
solar constant (W/m2)-liu(2002)
real(kind=kind_phys) sindec
subroutine solar
This subroutine computes radius vector, declination and right ascension of sun, and equation of time...
This module contains some the most frequently used math and physics constants for gcm models...
character(26) solar_fname
real(kind=kind_phys), parameter f12
real(kind=kind_phys), public solc0
character(40), parameter vtagast
This module defines commonly used control variables/parameters in physics related programs...
integer, save isolar
solar constant scheme control flag
subroutine prtime
This subroutine prints out forecast date, time, and astronomy quantities.
real(kind=kind_phys), parameter con_solr
solar constant (W/m2)-nasa-sorce tim(2008)
real(kind=kind_phys), parameter pid12
real(kind=kind_phys), parameter f3600
subroutine, public sol_init
This subroutine initializes astronomy process, and set up module constants.
This module sets up astronomy quantities for solar radiation calculations.
real(kind=kind_phys), parameter hpi
real(kind=kind_phys), parameter con_pi
real(kind=kind_phys), dimension(12) smon_sav
character, save solar_file
external solar constant data table