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
329 & ( jdate,kyear,deltsw,deltim,lsol_chg, me, &
330 & slag, sdec, cdec, solcon, con_pi, errmsg, errflg &
384 integer,
intent(in) :: jdate(:), kyear, me
385 logical,
intent(in) :: lsol_chg
387 real (kind=kind_phys),
intent(in) :: deltsw, deltim, con_pi
390 real (kind=kind_phys),
intent(out) :: slag, sdec, cdec, solcon
391 character(len=*),
intent(out) :: errmsg
392 integer,
intent(out) :: errflg
395 real (kind=kind_phys),
parameter :: hrday = 1.0/24.0
396 real (kind=kind_phys),
parameter :: minday= 1.0/1440.0
397 real (kind=kind_phys),
parameter :: secday= 1.0/86400.0
399 real (kind=kind_phys) :: smean, solc1, dtswh, smon(12)
400 real (kind=kind_phys) :: fjd, fjd1, dlt, r1, alp
402 integer :: jd, jd1, iyear, imon, iday, ihr, imin, isec
404 integer :: i, iyr, iyr1, iyr2, jyr, nn, nswr, icy1, icy2, icy
406 logical :: file_exist
407 character :: cline*60
435 if ( .not. file_exist )
then
436 print *,
' !!! ERROR! Can not find solar constant file!!!'
438 errmsg =
"ERROR(radiation_astronomy): solar constant file"//&
449 read (niradsf, * ) iyr1,iyr2,icy1,icy2,smean,cline(1:60)
454 print *,
' Updating solar constant with cycle approx'
455 print *,
' Opened solar constant data file: ',
solar_fname
494 if ( iyr < iyr1 )
then
495 icy = icy1 - iyr1 + 1
496 lab_dowhile1 :
do while ( iyr < iyr1 )
501 print *,
' *** Year',iyear,
' out of table range!', &
503 print *,
' Using the closest-cycle year (',iyr,
')'
505 elseif ( iyr > iyr2 )
then
506 icy = iyr2 - icy2 + 1
507 lab_dowhile2 :
do while ( iyr > iyr2 )
512 print *,
' *** Year',iyear,
' out of table range!', &
514 print *,
' Using the closest-cycle year (',iyr,
')'
522 lab_dowhile3 :
do while ( i >= iyr1 )
525 read (niradsf,*) jyr, solc1
527 if ( i == iyr .and. iyr == jyr )
then
528 solc0 = smean + solc1
531 print *,
' CHECK: Solar constant data used for year',&
542 lab_dowhile4 :
do while ( i >= iyr1 )
545 read (niradsf,*) jyr, smon(1:12)
547 if ( i == iyr .and. iyr == jyr )
then
551 solc0 = smean + smon(imon)
554 print *,
' CHECK: Solar constant data used for year',&
555 & iyr,
' and month',imon
573 jd1 = iw3jdn(iyear,imon,iday)
580 fjd1= 0.5 + float(ihr)*hrday + float(imin)*minday &
581 & + float(isec)*secday
583 fjd1= float(ihr - 12)*hrday + float(imin)*minday &
584 & + float(isec)*secday
595 & ( jd, fjd, con_pi, &
601 solcon =
solc0 / (r1*r1)
614 & ( jd, fjd, dlt, alp, r1, solcon )
621 nswr = max(1, nint(deltsw/deltim))
622 dtswh = deltsw /
f3600
638 print *,
' for cosz calculations: nswr,deltim,deltsw,dtswh =', &
639 & nswr,deltim,deltsw,dtswh,
' anginc,nstp =',
anginc,
nstp
659 & ( jd, fjd, con_pi, &
691 real (kind=kind_phys),
intent(in) :: fjd, con_pi
692 integer,
intent(in) :: jd
695 real (kind=kind_phys),
intent(out) :: r1, dlt, alp
698 real (kind=kind_phys),
parameter :: cyear = 365.25
699 real (kind=kind_phys),
parameter :: ccr = 1.3e-6
700 real (kind=kind_phys),
parameter :: tpp = 1.55
702 real (kind=kind_phys),
parameter :: svt6 = 78.035
704 integer,
parameter :: jdor = 2415020
707 real (kind=kind_phys) :: dat, t1, year, tyear, ec, angin, ador, &
708 & deleqn, sni, tini, er, qq, e1, ep, cd, eq, date, em, &
711 integer :: jdoe, iter
717 t1 = float(jd - jdor) / 36525.0
721 year = 0.25964134e0 + 0.304e-5 * t1
722 tyear= 0.24219879e0 - 0.614e-5 * t1
726 ec = 0.01675104e0 - (0.418e-4 + 0.126e-6 * t1) * t1
727 angin= 23.452294e0 - (0.0130125e0 + 0.164e-5 * t1) * t1
730 jdoe = ador + (svt6 * cyear) / (year - tyear)
734 deleqn= float(jdoe - jd) * (year - tyear) / cyear
736 sni = sin( angin /
degrad )
737 tini = 1.0 / tan( angin /
degrad )
738 er = sqrt( (1.0 + ec) / (1.0 - ec) )
739 qq = deleqn *
tpi / year
747 lab_do_1 :
do while ( cd > ccr )
749 ep = e1 - (e1 - ec*sin(e1) - qq) / (1.0 - ec*cos(e1))
755 write(6,*)
' ITERATION COUNT FOR LOOP 32 =', iter
756 write(6,*)
' E, EP, CD =', e1, ep, cd
762 eq = 2.0 * atan( er * tan( 0.5*e1 ) )
766 dat = float(jd - jdor) - tpp + fjd
767 date = mod(dat, year)
771 em =
tpi * date / year
776 lab_do_2 :
do while ( cr > ccr )
778 ep = e1 - (e1 - ec*sin(e1) - em) / (1.0 - ec*cos(e1))
784 write(6,*)
' ITERATION COUNT FOR LOOP 31 =', iter
790 w1 = 2.0 * atan( er * tan( 0.5*e1 ) )
792 r1 = 1.0 - ec*cos(e1)
794 sindec = sni * sin(w1 - eq)
798 alp = asin( tan(dlt)*tini )
801 if (tst < 0.0) alp = con_pi - alp
802 if (alp < 0.0) alp = alp +
tpi
804 sun =
tpi * (date - deleqn) / year
805 if (sun < 0.0) sun = sun +
tpi
806 sollag = sun - alp - 0.03255e0
828 & ( xlon,sinlat,coslat,solhr, im, me, &
865 integer,
intent(in) :: im, me
867 real (kind=kind_phys),
intent(in) :: sinlat(:), coslat(:), &
871 real (kind=kind_phys),
intent(out) :: coszen(:), coszdg(:)
874 real (kind=kind_phys) :: coszn, cns, solang, rstp
876 integer :: istsun(im), i, it, j, lat
881 rstp = 1.0 / float(
nstp)
893 coszen(i) = coszen(i) + max(0.0, coszn)
894 if (coszn >
czlimt) istsun(i) = istsun(i) + 1
901 coszdg(i) = coszen(i) * rstp
902 if (istsun(i) > 0 .and. coszen(i) /= 0.0_kind_phys)
then
903 coszen(i) = coszen(i) / istsun(i)
923 & ( jd, fjd, dlt, alp, r1, solc &
952 integer,
intent(in) :: jd
954 real (kind=kind_phys),
intent(in) :: fjd, dlt, alp, r1, solc
959 real (kind=kind_phys),
parameter :: sixty = 60.0
961 character(LEN=1),
parameter :: sign =
'-'
962 character(LEN=1),
parameter :: sigb =
' '
964 character(LEN=1) :: dsig
965 character(LEN=4) :: month(12)
967 data month /
'JAN.',
'FEB.',
'MAR.',
'APR.',
'MAY ',
'JUNE', &
968 &
'JULY',
'AUG.',
'SEP.',
'OCT.',
'NOV ',
'DEC.' /
970 integer :: iday, imon, iyear, ihr, ltd, ltm, &
971 & ihalp, iyy, jda, mfjd, idaywk, idayyr
972 real (kind=kind_phys) :: xmin, dltd, dltm, dlts, halp, ymin, &
981 mfjd= nint( fjd*1440.0 )
983 xmin= float(mfjd) - (ihr + 12)*sixty
986 mfjd= nint( fjd*1440.0 )
988 xmin= float(mfjd) - (ihr - 12)*sixty
993 call w3fs26(jda, iyear,imon,iday, idaywk,idayyr)
999 dltm = sixty * (abs(dltd) - abs(float(ltd)))
1001 dlts = sixty * (dltm - float(ltm))
1003 if ((dltd < 0.0) .and. (ltd == 0.0))
then
1009 halp = 6.0 * alp /
hpi
1011 ymin = abs(halp - float(ihalp)) * sixty
1013 asec = (ymin - float(iyy)) * sixty
1018 print 101, iday, month(imon), iyear, ihr, xmin, jd, fjd
1019 101
format(
'0 FORECAST DATE',9x,i3,a5,i6,
' AT',i3,
' HRS',f6.2,
' MINS'/&
1020 &
' JULIAN DAY',12x,i8,2x,
'PLUS',f11.6)
1022 print 102, r1, halp, ihalp, iyy, asec
1023 102
format(
' RADIUS VECTOR',9x,f10.7/
' RIGHT ASCENSION OF SUN', &
1024 &
f12.7,
' HRS, OR',i4,
' HRS',i4,
' MINS',f6.1,
' SECS')
1026 print 103, dltd, dsig, ltd, ltm, dlts, eqt, eqsec,
sollag, solc
1027 103
format(
' DECLINATION OF THE SUN',
f12.7,
' DEGS, OR ',a1,i3, &
1028 &
' DEGS',i4,
' MINS',f6.1,
' SECS'/
' EQUATION OF TIME',6x, &
1029 &
f12.7,
' MINS, OR',f10.2,
' SECS, OR',f9.6,
' RADIANS'/ &
1030 &
' 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.