148 & ( me, isolar, solar_file, con_solr, con_solr_old, con_pi )
183 integer,
intent(in) :: me, isolar
184 character(len=26),
intent(in) :: solar_file
185 real(kind=kind_phys),
intent(in) :: con_solr, con_solr_old, con_pi
189 logical :: file_exist
194 if ( me == 0 ) print *, vtagast
211 if ( isolar == 0 )
then
214 print *,
' - Using old fixed solar constant =',
solc0
216 elseif ( isolar == 10 )
then
218 print *,
' - Using new fixed solar constant =',
solc0
220 elseif ( isolar == 1 )
then
224 print *,
' - Using NOAA annual mean TSI table in ABS scale', &
225 &
' with cycle approximation (old values)!'
229 if ( .not. file_exist )
then
233 print *,
' Requested solar data file "',
solar_fname, &
235 print *,
' Using the default solar constant value =',
solc0,&
236 &
' reset control flag isolflg=',
isolflg
239 elseif ( isolar == 2 )
then
243 print *,
' - Using NOAA annual mean TSI table in TIM scale', &
244 &
' with cycle approximation (new values)!'
248 if ( .not. file_exist )
then
252 print *,
' Requested solar data file "',
solar_fname, &
254 print *,
' Using the default solar constant value =',
solc0,&
255 &
' reset control flag isolflg=',
isolflg
258 elseif ( isolar == 3 )
then
262 print *,
' - Using CMIP5 annual mean TSI table in TIM scale', &
263 &
' with cycle approximation'
267 if ( .not. file_exist )
then
271 print *,
' Requested solar data file "',
solar_fname, &
273 print *,
' Using the default solar constant value =',
solc0,&
274 &
' reset control flag isolflg=',
isolflg
277 elseif ( isolar == 4 )
then
281 print *,
' - Using CMIP5 monthly mean TSI table in TIM scale', &
282 &
' with cycle approximation'
286 if ( .not. file_exist )
then
290 print *,
' Requested solar data file "',
solar_fname, &
292 print *,
' Using the default solar constant value =',
solc0,&
293 &
' reset control flag isolflg=',
isolflg
300 print *,
' - !!! ERROR in selection of solar constant data', &
301 &
' source, ISOL =',isolar
302 print *,
' Using the default solar constant value =',
solc0, &
303 &
' reset control flag isolflg=',
isolflg
328 & ( jdate,kyear,deltsw,deltim,lsol_chg, me, &
329 & slag, sdec, cdec, solcon, con_pi, errmsg, errflg &
383 integer,
intent(in) :: jdate(:), kyear, me
384 logical,
intent(in) :: lsol_chg
386 real (kind=kind_phys),
intent(in) :: deltsw, deltim, con_pi
389 real (kind=kind_phys),
intent(out) :: slag, sdec, cdec, solcon
390 character(len=*),
intent(out) :: errmsg
391 integer,
intent(out) :: errflg
394 real (kind=kind_phys),
parameter :: hrday = 1.0/24.0
395 real (kind=kind_phys),
parameter :: minday= 1.0/1440.0
396 real (kind=kind_phys),
parameter :: secday= 1.0/86400.0
398 real (kind=kind_phys) :: smean, solc1, dtswh, smon(12)
399 real (kind=kind_phys) :: fjd, fjd1, dlt, r1, alp
401 integer :: jd, jd1, iyear, imon, iday, ihr, imin, isec
403 integer :: i, iyr, iyr1, iyr2, jyr, nn, nswr, icy1, icy2, icy
405 logical :: file_exist
406 character :: cline*60
434 if ( .not. file_exist )
then
436 errmsg =
"ERROR(radiation_astronomy): solar constant file"//&
447 read (niradsf, * ) iyr1,iyr2,icy1,icy2,smean,cline(1:60)
452 print *,
' Updating solar constant with cycle approx'
453 print *,
' Opened solar constant data file: ',
solar_fname
492 if ( iyr < iyr1 )
then
493 icy = icy1 - iyr1 + 1
494 lab_dowhile1 :
do while ( iyr < iyr1 )
499 print *,
' *** Year',iyear,
' out of table range!', &
501 print *,
' Using the closest-cycle year (',iyr,
')'
503 elseif ( iyr > iyr2 )
then
504 icy = iyr2 - icy2 + 1
505 lab_dowhile2 :
do while ( iyr > iyr2 )
510 print *,
' *** Year',iyear,
' out of table range!', &
512 print *,
' Using the closest-cycle year (',iyr,
')'
520 lab_dowhile3 :
do while ( i >= iyr1 )
523 read (niradsf,*) jyr, solc1
525 if ( i == iyr .and. iyr == jyr )
then
526 solc0 = smean + solc1
529 print *,
' CHECK: Solar constant data used for year',&
540 lab_dowhile4 :
do while ( i >= iyr1 )
543 read (niradsf,*) jyr, smon(1:12)
545 if ( i == iyr .and. iyr == jyr )
then
549 solc0 = smean + smon(imon)
552 print *,
' CHECK: Solar constant data used for year',&
553 & iyr,
' and month',imon
571 jd1 = iw3jdn(iyear,imon,iday)
578 fjd1= 0.5 + float(ihr)*hrday + float(imin)*minday &
579 & + float(isec)*secday
581 fjd1= float(ihr - 12)*hrday + float(imin)*minday &
582 & + float(isec)*secday
593 & ( jd, fjd, con_pi, &
599 solcon =
solc0 / (r1*r1)
612 & ( jd, fjd, dlt, alp, r1, solcon )
619 nswr = max(1, nint(deltsw/deltim))
620 dtswh = deltsw /
f3600
636 print *,
' for cosz calculations: nswr,deltim,deltsw,dtswh =', &
637 & nswr,deltim,deltsw,dtswh,
' anginc,nstp =',
anginc,
nstp
656 & ( jd, fjd, con_pi, &
688 real (kind=kind_phys),
intent(in) :: fjd, con_pi
689 integer,
intent(in) :: jd
692 real (kind=kind_phys),
intent(out) :: r1, dlt, alp
695 real (kind=kind_phys),
parameter :: cyear = 365.25
696 real (kind=kind_phys),
parameter :: ccr = 1.3e-6
697 real (kind=kind_phys),
parameter :: tpp = 1.55
699 real (kind=kind_phys),
parameter :: svt6 = 78.035
701 integer,
parameter :: jdor = 2415020
704 real (kind=kind_phys) :: dat, t1, year, tyear, ec, angin, ador, &
705 & deleqn, sni, tini, er, qq, e1, ep, cd, eq, date, em, &
708 integer :: jdoe, iter
714 t1 = float(jd - jdor) / 36525.0
718 year = 0.25964134e0 + 0.304e-5 * t1
719 tyear= 0.24219879e0 - 0.614e-5 * t1
723 ec = 0.01675104e0 - (0.418e-4 + 0.126e-6 * t1) * t1
724 angin= 23.452294e0 - (0.0130125e0 + 0.164e-5 * t1) * t1
727 jdoe = ador + (svt6 * cyear) / (year - tyear)
731 deleqn= float(jdoe - jd) * (year - tyear) / cyear
733 sni = sin( angin /
degrad )
734 tini = 1.0 / tan( angin /
degrad )
735 er = sqrt( (1.0 + ec) / (1.0 - ec) )
736 qq = deleqn *
tpi / year
744 lab_do_1 :
do while ( cd > ccr )
746 ep = e1 - (e1 - ec*sin(e1) - qq) / (1.0 - ec*cos(e1))
752 write(6,*)
' ITERATION COUNT FOR LOOP 32 =', iter
753 write(6,*)
' E, EP, CD =', e1, ep, cd
759 eq = 2.0 * atan( er * tan( 0.5*e1 ) )
763 dat = float(jd - jdor) - tpp + fjd
764 date = mod(dat, year)
768 em =
tpi * date / year
773 lab_do_2 :
do while ( cr > ccr )
775 ep = e1 - (e1 - ec*sin(e1) - em) / (1.0 - ec*cos(e1))
781 write(6,*)
' ITERATION COUNT FOR LOOP 31 =', iter
787 w1 = 2.0 * atan( er * tan( 0.5*e1 ) )
789 r1 = 1.0 - ec*cos(e1)
791 sindec = sni * sin(w1 - eq)
795 alp = asin( tan(dlt)*tini )
798 if (tst < 0.0) alp = con_pi - alp
799 if (alp < 0.0) alp = alp +
tpi
801 sun =
tpi * (date - deleqn) / year
802 if (sun < 0.0) sun = sun +
tpi
803 sollag = sun - alp - 0.03255e0
824 & ( xlon,sinlat,coslat,solhr, im, me, &
861 integer,
intent(in) :: im, me
863 real (kind=kind_phys),
intent(in) :: sinlat(:), coslat(:), &
867 real (kind=kind_phys),
intent(out) :: coszen(:), coszdg(:)
870 real (kind=kind_phys) :: coszn, cns, solang, rstp
872 integer :: istsun(im), i, it, j, lat
877 rstp = 1.0 / float(
nstp)
889 coszen(i) = coszen(i) + max(0.0, coszn)
890 if (coszn >
czlimt) istsun(i) = istsun(i) + 1
897 coszdg(i) = coszen(i) * rstp
898 if (istsun(i) > 0 .and. coszen(i) /= 0.0_kind_phys)
then
899 coszen(i) = coszen(i) / istsun(i)
918 & ( jd, fjd, dlt, alp, r1, solc &
947 integer,
intent(in) :: jd
949 real (kind=kind_phys),
intent(in) :: fjd, dlt, alp, r1, solc
954 real (kind=kind_phys),
parameter :: sixty = 60.0
956 character(LEN=1),
parameter :: sign =
'-'
957 character(LEN=1),
parameter :: sigb =
' '
959 character(LEN=1) :: dsig
960 character(LEN=4) :: month(12)
962 data month /
'JAN.',
'FEB.',
'MAR.',
'APR.',
'MAY ',
'JUNE', &
963 &
'JULY',
'AUG.',
'SEP.',
'OCT.',
'NOV ',
'DEC.' /
965 integer :: iday, imon, iyear, ihr, ltd, ltm, &
966 & ihalp, iyy, jda, mfjd, idaywk, idayyr
967 real (kind=kind_phys) :: xmin, dltd, dltm, dlts, halp, ymin, &
976 mfjd= nint( fjd*1440.0 )
978 xmin= float(mfjd) - (ihr + 12)*sixty
981 mfjd= nint( fjd*1440.0 )
983 xmin= float(mfjd) - (ihr - 12)*sixty
988 call w3fs26(jda, iyear,imon,iday, idaywk,idayyr)
994 dltm = sixty * (abs(dltd) - abs(float(ltd)))
996 dlts = sixty * (dltm - float(ltm))
998 if ((dltd < 0.0) .and. (ltd == 0.0))
then
1004 halp = 6.0 * alp /
hpi
1006 ymin = abs(halp - float(ihalp)) * sixty
1008 asec = (ymin - float(iyy)) * sixty
1013 print 101, iday, month(imon), iyear, ihr, xmin, jd, fjd
1014 101
format(
'0 FORECAST DATE',9x,i3,a5,i6,
' AT',i3,
' HRS',f6.2,
' MINS'/&
1015 &
' JULIAN DAY',12x,i8,2x,
'PLUS',f11.6)
1017 print 102, r1, halp, ihalp, iyy, asec
1018 102
format(
' RADIUS VECTOR',9x,f10.7/
' RIGHT ASCENSION OF SUN', &
1019 &
f12.7,
' HRS, OR',i4,
' HRS',i4,
' MINS',f6.1,
' SECS')
1021 print 103, dltd, dsig, ltd, ltm, dlts, eqt, eqsec,
sollag, solc
1022 103
format(
' DECLINATION OF THE SUN',
f12.7,
' DEGS, OR ',a1,i3, &
1023 &
' DEGS',i4,
' MINS',f6.1,
' SECS'/
' EQUATION OF TIME',6x, &
1024 &
f12.7,
' MINS, OR',f10.2,
' SECS, OR',f9.6,
' RADIANS'/ &
1025 &
' SOLAR CONSTANT',8x,
f12.7,
' (DISTANCE AJUSTED)'//)
subroutine, public sol_init(me, isolar, solar_file, con_solr, con_solr_old, con_pi)
This subroutine initializes astronomy process, and set up module constants.
subroutine, public sol_update(jdate, kyear, deltsw, deltim, lsol_chg, me, slag, sdec, cdec, solcon, con_pi, errmsg, errflg)
This subroutine computes solar parameters at forecast time.
subroutine solar(jd, fjd, con_pi, r1, dlt, alp)
This subroutine computes radius vector, declination and right ascension of sun, and equation of time.
subroutine, public coszmn(xlon, sinlat, coslat, solhr, im, me, coszen, coszdg)
This subroutine computes mean cos solar zenith angle over SW calling interval.
subroutine prtime(jd, fjd, dlt, alp, r1, solc)
This subroutine prints out forecast date, time, and astronomy quantities.