CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
radiation_astronomy.f
1
4
5! ========================================================== !!!!!
6! 'module_radiation_astronomy' description !!!!!
7! ========================================================== !!!!!
8! !
9! set up astronomy quantities for solar radiation calculations. !
10! !
11! in module 'module_radiation_astronomy', externally accessable !
12! subroutines are listed below: !
13! !
14! 'sol_init' -- initialization !
15! input: !
16! ( me ) !
17! output: !
18! ( none ) !
19! !
20! 'sol_update' -- update astronomy related quantities !
21! input: !
22! ( jdate,kyear,deltsw,deltim,lsol_chg, me ) !
23! output: !
24! ( slag,sdec,cdec,solcon,errmsg,errflg) !
25! !
26! 'coszmn' -- compute cosin of zenith angles !
27! input: !
28! ( xlon,sinlat,coslat,solhr,IM, me ) !
29! output: !
30! ( coszen,coszdg ) !
31! !
32! program history log: !
33! - a collection of programs to track solar-earth position !
34! may 1977 --- ray orzol (gfdl) created program compjd to !
35! computes julian day and fraction from year,month,dayand,time!
36! jun 1977 --- robert white (gfdl) created program cdate to !
37! computes calendar month, day, year from julian day value. !
38! jul 1977 --- robert white (gfdl) created program solar to !
39! computes radius vector, declination and right ascension of !
40! sun, equation of time, hour angle, fractional daylight, and !
41! latitudinal mean zenith angle. !
42! fall 1988 --- hualu pan, updated to limit the iterations in !
43! newton method and also ccr reduced to avoid non-convergence.!
44! jul 1989 --- kenneth campana modified subr solar and created !
45! subr zenith for computations of effective mean cosz and !
46! daylight fraction. !
47! oct 1990 --- yu-tai hou created subr coszmn to replace !
48! the latitudinal mean cosz by time mean cosz at grid points. !
49! may 1998 --- mark iredell y2k compliance !
50! dec 2003 --- yu-tai hou combined compjd and fcstim and !
51! rewritten programs in fortran 90 compatable modular form. !
52! feb 2006 --- yu-tai hou add 11-yr solar constant cycle !
53! mar 2009 --- yu-tai hou modified solinit for climate !
54! hindcast situation responding to ic time. !
55! aug 2012 --- yu-tai hou modified coszmn to allows sw !
56! radiation calling interval less than 1 hr limit and linked !
57! model time step with numb of cosz evaluations. also changed !
58! the initialization subroutine 'solinit' into two parts: !
59! 'sol_init' is called at the start of run to set up module !
60! parameters; and 'sol_update' is called within the time !
61! loop to check and update data sets. !
62! nov 2012 --- yu-tai hou modified control parameters thru !
63! model 'physparam'. !
64! jan 2013 --- yu-tai hou modified to include new solar !
65! constant tables (noaa_a0, noaa_an, cmip_an, cmip_mn) !
66! !
67!!!!! ========================================================== !!!!!
68!!!!! end descriptions !!!!!
69!!!!! ========================================================== !!!!!
70
71
72
87
90!
91 use machine, only : kind_phys
92 use module_iounitdef, only : niradsf
93!
94 implicit none
95!
96 private
97
98! --- version tag and last revision date
99 character(40), parameter :: &
100 & VTAGAST='NCEP-Radiation_astronomy v5.2 Jan 2013 '
101! & VTAGAST='NCEP-Radiation_astronomy v5.1 Nov 2012 '
102
103! Parameter constants
104 real (kind=kind_phys) :: degrad
105 real (kind=kind_phys) :: tpi
106 real (kind=kind_phys) :: hpi
107 real (kind=kind_phys) :: pid12
108 real (kind=kind_phys), parameter :: f12 = 12.0
109 real (kind=kind_phys), parameter :: f3600 = 3600.0
110 real (kind=kind_phys), parameter :: czlimt = 0.0001 ! ~ cos(89.99427)
111! real (kind=kind_phys), parameter :: pid12 = (2.0*asin(1.0))/f12
112
113! Module variable (to be set in module_radiation_astronomy::sol_init):
114 real (kind=kind_phys), public :: solc0
115 integer :: isolflg = 10
116 character(26) :: solar_fname = ' '
117
118! Module variables (to be set in module_radiation_astronomy::sol_update)
119
120! equation of time
121 real (kind=kind_phys) :: sollag=0.0
122! sine of the solar declination angle
123 real (kind=kind_phys) :: sindec=0.0
124! cosine of the solar declination angle
125 real (kind=kind_phys) :: cosdec=0.0
126! solar angle increment per interation of cosz calc
127 real (kind=kind_phys) :: anginc=0.0
128! saved monthly solar constants (isolflg=4 only)
129 real (kind=kind_phys) :: smon_sav(12)
130
131! saved year of data used
132 integer :: iyr_sav =0
133! total number of zenith angle iterations
134 integer :: nstp =6
135
136 public sol_init, sol_update, coszmn
137
138
139! =================
140 contains
141! =================
142
147 subroutine sol_init &
148 & ( me, isolar, solar_file, con_solr, con_solr_old, con_pi ) ! --- inputs
149! --- outputs: ( none )
150
151! =================================================================== !
152! !
153! initialize astronomy process, set up module constants. !
154! !
155! inputs: !
156! me - print message control flag !
157! isolar - = 0: use the old fixed solar constant in "GFS_typedefs" !
158! =10: use the new fixed solar constant in "GFS_typedefs" !
159! = 1: use noaa ann-mean tsi tbl abs-scale with cyc apprx !
160! = 2: use noaa ann-mean tsi tbl tim-scale with cyc apprx !
161! = 3: use cmip5 ann-mean tsi tbl tim-scale with cyc apprx!
162! = 4: use cmip5 mon-mean tsi tbl tim-scale with cyc apprx!
163! solar_file - external solar constant data table !
164! !
165! outputs: (to module variable) !
166! ( none ) !
167! !
168! internal module variable: !
169! isolflg - internal solar constant scheme control flag !
170! solc0 - solar constant (w/m**2) !
171! solar_fname-file name for solar constant table assigned based on !
172! the scheme control flag, isolflg. !
173! !
174! usage: call sol_init !
175! !
176! subprograms called: none !
177! !
178! =================================================================== !
179!
180 implicit none
181
182! --- input:
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
186! --- output: ( none )
187
188! --- local:
189 logical :: file_exist
190 integer :: imonth
191!
192!===> ... begin here
193!
194 if ( me == 0 ) print *, vtagast !print out version tag
195
196 degrad = 180.0/con_pi
197 tpi = 2.0 * con_pi
198 hpi = 0.5 * con_pi
199 pid12 = con_pi/f12
200
201! --- initialization
202 isolflg = isolar
203 solc0 = con_solr
204 solar_fname = solar_file
205 iyr_sav = 0
206 nstp = 6
207 do imonth = 1,12
208 smon_sav(imonth) = con_solr
209 enddo
210
211 if ( isolar == 0 ) then
212 solc0 = con_solr_old
213 if ( me == 0 ) then
214 print *,' - Using old fixed solar constant =', solc0
215 endif
216 elseif ( isolar == 10 ) then
217 if ( me == 0 ) then
218 print *,' - Using new fixed solar constant =', solc0
219 endif
220 elseif ( isolar == 1 ) then ! noaa ann-mean tsi in absolute scale
221 solar_fname(15:26) = 'noaa_a0.txt'
222
223 if ( me == 0 ) then
224 print *,' - Using NOAA annual mean TSI table in ABS scale', &
225 & ' with cycle approximation (old values)!'
226 endif
227
228 inquire (file=solar_fname, exist=file_exist)
229 if ( .not. file_exist ) then
230 isolflg = 10
231
232 if ( me == 0 ) then
233 print *,' Requested solar data file "',solar_fname, &
234 & '" not found!'
235 print *,' Using the default solar constant value =',solc0,&
236 & ' reset control flag isolflg=',isolflg
237 endif
238 endif
239 elseif ( isolar == 2 ) then ! noaa ann-mean tsi in tim scale
240 solar_fname(15:26) = 'noaa_an.txt'
241
242 if ( me == 0 ) then
243 print *,' - Using NOAA annual mean TSI table in TIM scale', &
244 & ' with cycle approximation (new values)!'
245 endif
246
247 inquire (file=solar_fname, exist=file_exist)
248 if ( .not. file_exist ) then
249 isolflg = 10
250
251 if ( me == 0 ) then
252 print *,' Requested solar data file "',solar_fname, &
253 & '" not found!'
254 print *,' Using the default solar constant value =',solc0,&
255 & ' reset control flag isolflg=',isolflg
256 endif
257 endif
258 elseif ( isolar == 3 ) then ! cmip5 ann-mean tsi in tim scale
259 solar_fname(15:26) = 'cmip_an.txt'
260
261 if ( me == 0 ) then
262 print *,' - Using CMIP5 annual mean TSI table in TIM scale', &
263 & ' with cycle approximation'
264 endif
265
266 inquire (file=solar_fname, exist=file_exist)
267 if ( .not. file_exist ) then
268 isolflg = 10
269
270 if ( me == 0 ) then
271 print *,' Requested solar data file "',solar_fname, &
272 & '" not found!'
273 print *,' Using the default solar constant value =',solc0,&
274 & ' reset control flag isolflg=',isolflg
275 endif
276 endif
277 elseif ( isolar == 4 ) then ! cmip5 mon-mean tsi in tim scale
278 solar_fname(15:26) = 'cmip_mn.txt'
279
280 if ( me == 0 ) then
281 print *,' - Using CMIP5 monthly mean TSI table in TIM scale', &
282 & ' with cycle approximation'
283 endif
284
285 inquire (file=solar_fname, exist=file_exist)
286 if ( .not. file_exist ) then
287 isolflg = 10
288
289 if ( me == 0 ) then
290 print *,' Requested solar data file "',solar_fname, &
291 & '" not found!'
292 print *,' Using the default solar constant value =',solc0,&
293 & ' reset control flag isolflg=',isolflg
294 endif
295 endif
296 else ! selection error
297 isolflg = 10
298
299 if ( me == 0 ) then
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
304 endif
305 endif ! end if_isolar_block
306!
307 return
308!...................................
309 end subroutine sol_init
310!-----------------------------------
311
327!-----------------------------------
328 subroutine sol_update &
329 & ( jdate,kyear,deltsw,deltim,lsol_chg, me, & ! --- inputs
330 & slag, sdec, cdec, solcon, con_pi, errmsg, errflg & ! --- outputs
331 & )
332
333! =================================================================== !
334! !
335! sol_update computes solar parameters at forecast time !
336! !
337! inputs: !
338! jdate(8)- ncep absolute date and time at fcst time !
339! (yr, mon, day, t-zone, hr, min, sec, mil-sec) !
340! kyear - usually kyear=jdate(1). if not, it is for hindcast mode,!
341! and it is usually the init cond time and serves as the !
342! upper limit of data can be used. !
343! deltsw - time duration in seconds per sw calculation !
344! deltim - timestep in seconds !
345! lsol_chg- logical flags for change solar constant !
346! me - print message control flag !
347! !
348! outputs: !
349! slag - equation of time in radians !
350! sdec, cdec - sin and cos of the solar declination angle !
351! solcon - sun-earth distance adjusted solar constant (w/m2) !
352! errmsg - CCPP error message !
353! errflg - CCPP error flag !
354! !
355! !
356! module variable: !
357! solc0 - solar constant (w/m**2) not adjusted by earth-sun dist !
358! isolflg - solar constant control flag !
359! = 0: use the old fixed solar constant !
360! =10: use the new fixed solar constant !
361! = 1: use noaa ann-mean tsi tbl abs-scale with cycle apprx !
362! = 2: use noaa ann-mean tsi tbl tim-scale with cycle apprx !
363! = 3: use cmip5 ann-mean tsi tbl tim-scale with cycle apprx!
364! = 4: use cmip5 mon-mean tsi tbl tim-scale with cycle apprx!
365! solar_fname-external solar constant data table !
366! sindec - sine of the solar declination angle !
367! cosdec - cosine of the solar declination angle !
368! anginc - solar angle increment per iteration for cosz calc !
369! nstp - total number of zenith angle iterations !
370! smon_sav- saved monthly solar constants (isolflg=4 only) !
371! iyr_sav - saved year of data previously used !
372! !
373! usage: call sol_update !
374! !
375! subprograms called: solar, prtime !
376! !
377! external functions called: iw3jdn !
378! !
379! =================================================================== !
380!
381 implicit none
382
383! --- input:
384 integer, intent(in) :: jdate(:), kyear, me
385 logical, intent(in) :: lsol_chg
386
387 real (kind=kind_phys), intent(in) :: deltsw, deltim, con_pi
388
389! --- output:
390 real (kind=kind_phys), intent(out) :: slag, sdec, cdec, solcon
391 character(len=*), intent(out) :: errmsg
392 integer, intent(out) :: errflg
393
394! --- locals:
395 real (kind=kind_phys), parameter :: hrday = 1.0/24.0 ! frc day/hour
396 real (kind=kind_phys), parameter :: minday= 1.0/1440.0 ! frc day/minute
397 real (kind=kind_phys), parameter :: secday= 1.0/86400.0 ! frc day/second
398
399 real (kind=kind_phys) :: smean, solc1, dtswh, smon(12)
400 real (kind=kind_phys) :: fjd, fjd1, dlt, r1, alp
401
402 integer :: jd, jd1, iyear, imon, iday, ihr, imin, isec
403 integer :: iw3jdn
404 integer :: i, iyr, iyr1, iyr2, jyr, nn, nswr, icy1, icy2, icy
405
406 logical :: file_exist
407 character :: cline*60
408!
409!===> ... begin here
410!
411! Initialize the CCPP error handling variables
412 errmsg = ''
413 errflg = 0
414
415! --- ... forecast time
416 iyear = jdate(1)
417 imon = jdate(2)
418 iday = jdate(3)
419 ihr = jdate(5)
420 imin = jdate(6)
421 isec = jdate(7)
422
423 if ( lsol_chg ) then ! get solar constant from data table
424
425 if ( iyr_sav == iyear ) then ! same year, no new reading necessary
426 if ( isolflg==4 ) then
427 solc0 = smon_sav(imon)
428 endif
429 else ! need to read in new data
430 iyr_sav = iyear
431
432! --- ... check to see if the solar constant data file existed
433
434 inquire (file=solar_fname, exist=file_exist)
435 if ( .not. file_exist ) then
436 print *,' !!! ERROR! Can not find solar constant file!!!'
437 errflg = 1
438 errmsg = "ERROR(radiation_astronomy): solar constant file"//&
439 & " not found"
440 return
441 else
442 iyr = iyear
443
444 close(niradsf)
445 open (niradsf,file=solar_fname,form='formatted', &
446 & status='old')
447 rewind niradsf
448
449 read (niradsf, * ) iyr1,iyr2,icy1,icy2,smean,cline(1:60)
450! read (NIRADSF, 24) iyr1,iyr2,icy1,icy2,smean,cline
451! 24 format(4i5,f8.2,a60)
452
453 if ( me == 0 ) then
454 print *,' Updating solar constant with cycle approx'
455 print *,' Opened solar constant data file: ',solar_fname
456!check print *, iyr1, iyr2, icy1, icy2, smean, cline
457 endif
458
459! --- ... check if there is a upper year limit put on the data table
460
461! if ( iyear /= kyear ) then
462! icy = icy1 - iyr1 + 1 ! range of the earlest cycle in data table
463! if ( kyear-iyr1 < icy ) then ! need data range at least icy years
464 ! to perform cycle approximation
465! if ( me == 0 ) then
466! print *,' *** the requested year',iyear,' and upper',&
467! & 'limit',kyear,' do not fit the range of data ', &
468! & 'table of iyr1, iyr2 =',iyr1,iyr2
469! print *,' USE FIXED SOLAR CONSTANT=',con_solr
470! endif
471! solc0 = con_solr
472! isolflg = 10
473
474! elseif ( kyear < iyr2 ) then
475
476! --- ... because the usage limit put on the historical data table,
477! skip those unused data records at first
478
479! i = iyr2
480! Lab_dowhile0 : do while ( i > kyear )
481! read (NIRADSF,26) jyr, solc1
482! 26 format(i4,f10.4)
483! read (NIRADSF,*) jyr, solc1
484! i = i - 1
485! enddo Lab_dowhile0
486
487! iyr2 = kyear ! next record will serve the upper limit
488
489! endif ! end if_kyear_block
490! endif ! end if_iyear_block
491
492! --- ... checking the cycle range
493
494 if ( iyr < iyr1 ) then
495 icy = icy1 - iyr1 + 1 ! range of the earlest cycle in data table
496 lab_dowhile1 : do while ( iyr < iyr1 )
497 iyr = iyr + icy
498 enddo lab_dowhile1
499
500 if ( me == 0 ) then
501 print *,' *** Year',iyear,' out of table range!', &
502 & iyr1, iyr2
503 print *,' Using the closest-cycle year (',iyr,')'
504 endif
505 elseif ( iyr > iyr2 ) then
506 icy = iyr2 - icy2 + 1 ! range of the latest cycle in data table
507 lab_dowhile2 : do while ( iyr > iyr2 )
508 iyr = iyr - icy
509 enddo lab_dowhile2
510
511 if ( me == 0 ) then
512 print *,' *** Year',iyear,' out of table range!', &
513 & iyr1, iyr2
514 print *,' Using the closest-cycle year (',iyr,')'
515 endif
516 endif
517
518! --- ... locate the right record for the year of data
519
520 if ( isolflg < 4 ) then ! use annual mean data tables
521 i = iyr2
522 lab_dowhile3 : do while ( i >= iyr1 )
523! read (NIRADSF,26) jyr, solc1
524! 26 format(i4,f10.4)
525 read (niradsf,*) jyr, solc1
526
527 if ( i == iyr .and. iyr == jyr ) then
528 solc0 = smean + solc1
529
530 if (me == 0) then
531 print *,' CHECK: Solar constant data used for year',&
532 & iyr, solc1, solc0
533 endif
534 exit lab_dowhile3
535 else
536!check if(me == 0) print *,' Skip solar const data for yr',i
537 i = i - 1
538 endif
539 enddo lab_dowhile3
540 elseif ( isolflg == 4 ) then ! use monthly mean data tables
541 i = iyr2
542 lab_dowhile4 : do while ( i >= iyr1 )
543! read (NIRADSF,26) jyr, smon(:)
544! 26 format(i4,12f10.4)
545 read (niradsf,*) jyr, smon(1:12)
546
547 if ( i == iyr .and. iyr == jyr ) then
548 do nn = 1, 12
549 smon_sav(nn) = smean + smon(nn)
550 enddo
551 solc0 = smean + smon(imon)
552
553 if (me == 0) then
554 print *,' CHECK: Solar constant data used for year',&
555 & iyr,' and month',imon
556 endif
557 exit lab_dowhile4
558 else
559!check if(me == 0) print *,' Skip solar const data for yr',i
560 i = i - 1
561 endif
562 enddo lab_dowhile4
563 endif ! end if_isolflg_block
564
565 close ( niradsf )
566 endif ! end if_file_exist_block
567
568 endif ! end if_iyr_sav_block
569 endif ! end if_lsol_chg_block
570
571! --- ... calculate forecast julian day and fraction of julian day
572
573 jd1 = iw3jdn(iyear,imon,iday)
574
575! --- ... unlike in normal applications, where day starts from 0 hr,
576! in astronomy applications, day stats from noon.
577
578 if (ihr < 12) then
579 jd1 = jd1 - 1
580 fjd1= 0.5 + float(ihr)*hrday + float(imin)*minday &
581 & + float(isec)*secday
582 else
583 fjd1= float(ihr - 12)*hrday + float(imin)*minday &
584 & + float(isec)*secday
585 endif
586
587 fjd1 = fjd1 + jd1
588
589 jd = int(fjd1)
590 fjd = fjd1 - jd
591
593 call solar &
594! --- inputs:
595 & ( jd, fjd, con_pi, &
596! --- outputs:
597 & r1, dlt, alp &
598 & )
599
600! --- ... calculate sun-earth distance adjustment factor appropriate to date
601 solcon = solc0 / (r1*r1)
602
603 slag = sollag
604 sdec = sindec
605 cdec = cosdec
606
607! --- ... diagnostic print out
608
609 if (me == 0) then
610
612 call prtime &
613! --- inputs:
614 & ( jd, fjd, dlt, alp, r1, solcon )
615! --- outputs: ( none )
616
617 endif
618
619! --- ... setting up calculation parameters used by subr coszmn
620
621 nswr = max(1, nint(deltsw/deltim)) ! number of mdl t-step per sw call
622 dtswh = deltsw / f3600 ! time length in hours
623
624! if ( deltsw >= f3600 ) then ! for longer sw call interval
625! nn = max(6, min(12, nint(f3600/deltim) )) ! num of calc per hour
626! nstp = nint(dtswh) * nn + 1 ! num of calc per sw call
627! else ! for shorter sw sw call interval
628! nstp = max(2, min(20, nswr)) + 1
629!! nn = nint( float(nstp-1)/dtswh )
630! endif
631
632! anginc = pid12 * dtswh / float(nstp-1) ! solar angle inc during each calc step
633
634 nstp = max(6, nswr)
635 anginc = pid12 * dtswh / float(nstp)
636
637 if ( me == 0 ) then
638 print *,' for cosz calculations: nswr,deltim,deltsw,dtswh =', &
639 & nswr,deltim,deltsw,dtswh,' anginc,nstp =',anginc,nstp
640 endif
641
642! if (me == 0) print*,'in sol_update completed sr solar'
643!
644 return
645!...................................
646 end subroutine sol_update
647!-----------------------------------
648
657!-----------------------------------
658 subroutine solar &
659 & ( jd, fjd, con_pi, & ! --- inputs
660 & r1, dlt, alp & ! --- outputs
661 & )
662
663! =================================================================== !
664! !
665! solar computes radius vector, declination and right ascension of !
666! sun, and equation of time. !
667! !
668! inputs: !
669! jd - julian day !
670! fjd - fraction of the julian day !
671! !
672! outputs: !
673! r1 - earth-sun radius vector !
674! dlt - declination of sun in radians !
675! alp - right ascension of sun in radians !
676! !
677! module variables: !
678! sollag - equation of time in radians !
679! sindec - sine of declination angle !
680! cosdec - cosine of declination angle !
681! !
682! usage: call solar !
683! !
684! external subroutines called: none !
685! !
686! =================================================================== !
687!
688 implicit none
689
690! --- inputs:
691 real (kind=kind_phys), intent(in) :: fjd, con_pi
692 integer, intent(in) :: jd
693
694! --- outputs:
695 real (kind=kind_phys), intent(out) :: r1, dlt, alp
696
697! --- locals:
698 real (kind=kind_phys), parameter :: cyear = 365.25 ! days of year
699 real (kind=kind_phys), parameter :: ccr = 1.3e-6 ! iteration limit
700 real (kind=kind_phys), parameter :: tpp = 1.55 ! days between epoch and
701 ! perihelion passage of 1900
702 real (kind=kind_phys), parameter :: svt6 = 78.035 ! days between perihelion passage
703 ! and march equinox of 1900
704 integer, parameter :: jdor = 2415020 ! jd of epoch which is january
705 ! 0, 1900 at 12 hours ut
706
707 real (kind=kind_phys) :: dat, t1, year, tyear, ec, angin, ador, &
708 & deleqn, sni, tini, er, qq, e1, ep, cd, eq, date, em, &
709 & cr, w1, tst, sun
710
711 integer :: jdoe, iter
712
713!===> ... begin here
714
715! --- ... computes time in julian centuries after epoch
716
717 t1 = float(jd - jdor) / 36525.0
718
719! --- ... computes length of anomalistic and tropical years (minus 365 days)
720
721 year = 0.25964134e0 + 0.304e-5 * t1
722 tyear= 0.24219879e0 - 0.614e-5 * t1
723
724! --- ... computes orbit eccentricity and angle of earth's inclination from t
725
726 ec = 0.01675104e0 - (0.418e-4 + 0.126e-6 * t1) * t1
727 angin= 23.452294e0 - (0.0130125e0 + 0.164e-5 * t1) * t1
728
729 ador = jdor
730 jdoe = ador + (svt6 * cyear) / (year - tyear)
731
732! --- ... deleqn is updated svt6 for current date
733
734 deleqn= float(jdoe - jd) * (year - tyear) / cyear
735 year = year + 365.0
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
740
741! --- ... determine true anomaly at equinox
742
743 e1 = 1.0
744 cd = 1.0
745 iter = 0
746
747 lab_do_1 : do while ( cd > ccr )
748
749 ep = e1 - (e1 - ec*sin(e1) - qq) / (1.0 - ec*cos(e1))
750 cd = abs(e1 - ep)
751 e1 = ep
752 iter = iter + 1
753
754 if (iter > 10) then
755 write(6,*) ' ITERATION COUNT FOR LOOP 32 =', iter
756 write(6,*) ' E, EP, CD =', e1, ep, cd
757 exit lab_do_1
758 endif
759
760 enddo lab_do_1
761
762 eq = 2.0 * atan( er * tan( 0.5*e1 ) )
763
764! --- ... date is days since last perihelion passage
765
766 dat = float(jd - jdor) - tpp + fjd
767 date = mod(dat, year)
768
769! --- ... solve orbit equations by newton's method
770
771 em = tpi * date / year
772 e1 = 1.0
773 cr = 1.0
774 iter = 0
775
776 lab_do_2 : do while ( cr > ccr )
777
778 ep = e1 - (e1 - ec*sin(e1) - em) / (1.0 - ec*cos(e1))
779 cr = abs(e1 - ep)
780 e1 = ep
781 iter = iter + 1
782
783 if (iter > 10) then
784 write(6,*) ' ITERATION COUNT FOR LOOP 31 =', iter
785 exit lab_do_2
786 endif
787
788 enddo lab_do_2
789
790 w1 = 2.0 * atan( er * tan( 0.5*e1 ) )
791
792 r1 = 1.0 - ec*cos(e1)
793
794 sindec = sni * sin(w1 - eq)
795 cosdec = sqrt( 1.0 - sindec*sindec )
796
797 dlt = asin( sindec )
798 alp = asin( tan(dlt)*tini )
799
800 tst = cos( w1 - eq )
801 if (tst < 0.0) alp = con_pi - alp
802 if (alp < 0.0) alp = alp + tpi
803
804 sun = tpi * (date - deleqn) / year
805 if (sun < 0.0) sun = sun + tpi
806 sollag = sun - alp - 0.03255e0
807!
808 return
809!...................................
810 end subroutine solar
811!-----------------------------------
812
826!-----------------------------------
827 subroutine coszmn &
828 & ( xlon,sinlat,coslat,solhr, im, me, & ! --- inputs
829 & coszen, coszdg & ! --- outputs
830 & )
831
832! =================================================================== !
833! !
834! coszmn computes mean cos solar zenith angle over sw calling interval !
835! !
836! inputs: !
837! xlon (IM) - grids' longitudes in radians, work both on zonal !
838! 0->2pi and -pi->+pi arrangements !
839! sinlat(IM) - sine of the corresponding latitudes !
840! coslat(IM) - cosine of the corresponding latitudes !
841! solhr - time after 00z in hours !
842! IM - num of grids in horizontal dimension !
843! me - print message control flag !
844! !
845! outputs: !
846! coszen(IM) - average of cosz for daytime only in sw call interval
847! coszdg(IM) - average of cosz over entire sw call interval !
848! !
849! module variables: !
850! sollag - equation of time !
851! sindec - sine of the solar declination angle !
852! cosdec - cosine of the solar declination angle !
853! anginc - solar angle increment per iteration for cosz calc !
854! nstp - total number of zenith angle iterations !
855! !
856! usage: call comzmn !
857! !
858! external subroutines called: none !
859! !
860! =================================================================== !
861!
862 implicit none
863
864! --- inputs:
865 integer, intent(in) :: im, me
866
867 real (kind=kind_phys), intent(in) :: sinlat(:), coslat(:), &
868 & xlon(:), solhr
869
870! --- outputs:
871 real (kind=kind_phys), intent(out) :: coszen(:), coszdg(:)
872
873! --- locals:
874 real (kind=kind_phys) :: coszn, cns, solang, rstp
875
876 integer :: istsun(im), i, it, j, lat
877
878!===> ... begin here
879
880 solang = pid12 * (solhr - f12) ! solar angle at present time
881 rstp = 1.0 / float(nstp)
882
883 do i = 1, im
884 coszen(i) = 0.0
885 istsun(i) = 0
886 enddo
887
888 do it = 1, nstp
889 cns = solang + (float(it)-0.5)*anginc + sollag
890 do i = 1, im
891 coszn = sindec * sinlat(i) + cosdec * coslat(i) &
892 & * cos(cns+xlon(i))
893 coszen(i) = coszen(i) + max(0.0, coszn)
894 if (coszn > czlimt) istsun(i) = istsun(i) + 1
895 enddo
896 enddo
897
898! --- ... compute time averages
899
900 do i = 1, im
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)
904 endif
905 enddo
906!
907 return
908!...................................
909 end subroutine coszmn
910!-----------------------------------
911
921!-----------------------------------
922 subroutine prtime &
923 & ( jd, fjd, dlt, alp, r1, solc & ! --- inputs
924 & ) ! --- outputs: ( none )
925
926! =================================================================== !
927! !
928! prtime prints out forecast date, time, and astronomy quantities. !
929! !
930! inputs: !
931! jd - forecast julian day !
932! fjd - forecast fraction of julian day !
933! dlt - declination angle of sun in radians !
934! alp - right ascension of sun in radians !
935! r1 - earth-sun radius vector in meter !
936! solc - solar constant in w/m^2 !
937! !
938! outputs: ( none ) !
939! !
940! module variables: !
941! sollag - equation of time in radians !
942! !
943! usage: call prtime !
944! !
945! external subroutines called: w3fs26 !
946! !
947! =================================================================== !
948!
949 implicit none
950
951! --- inputs:
952 integer, intent(in) :: jd
953
954 real (kind=kind_phys), intent(in) :: fjd, dlt, alp, r1, solc
955
956! --- outputs: ( none )
957
958! --- locals:
959 real (kind=kind_phys), parameter :: sixty = 60.0
960
961 character(LEN=1), parameter :: sign = '-'
962 character(LEN=1), parameter :: sigb = ' '
963
964 character(LEN=1) :: dsig
965 character(LEN=4) :: month(12)
966
967 data month / 'JAN.','FEB.','MAR.','APR.','MAY ','JUNE', &
968 & 'JULY','AUG.','SEP.','OCT.','NOV ','DEC.' /
969
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, &
973 & asec, eqt, eqsec
974
975!===> ... begin here
976
977! --- ... get forecast hour and minute from fraction of julian day
978
979 if (fjd >= 0.5) then
980 jda = jd + 1
981 mfjd= nint( fjd*1440.0 )
982 ihr = mfjd / 60 - 12
983 xmin= float(mfjd) - (ihr + 12)*sixty
984 else
985 jda = jd
986 mfjd= nint( fjd*1440.0 )
987 ihr = mfjd / 60 + 12
988 xmin= float(mfjd) - (ihr - 12)*sixty
989 endif
990
991! --- ... get forecast year, month, and day from julian day
992
993 call w3fs26(jda, iyear,imon,iday, idaywk,idayyr)
994
995! -- ... compute solar parameters
996
997 dltd = degrad * dlt
998 ltd = dltd
999 dltm = sixty * (abs(dltd) - abs(float(ltd)))
1000 ltm = dltm
1001 dlts = sixty * (dltm - float(ltm))
1002
1003 if ((dltd < 0.0) .and. (ltd == 0.0)) then
1004 dsig = sign
1005 else
1006 dsig = sigb
1007 endif
1008
1009 halp = 6.0 * alp / hpi
1010 ihalp= halp
1011 ymin = abs(halp - float(ihalp)) * sixty
1012 iyy = ymin
1013 asec = (ymin - float(iyy)) * sixty
1014
1015 eqt = 228.55735 * sollag
1016 eqsec= sixty * eqt
1017
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)
1021
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')
1025
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)'//)
1031
1032!
1033 return
1034!...................................
1035 end subroutine prtime
1036!-----------------------------------
1037
1038!
1039!...........................................!
1040 end module module_radiation_astronomy !
1042!===========================================!
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.
real(kind=kind_phys), parameter f3600
real(kind=kind_phys), public solc0
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.
real(kind=kind_phys), dimension(12) smon_sav
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.
real(kind=kind_phys), parameter f12
subroutine prtime(jd, fjd, dlt, alp, r1, solc)
This subroutine prints out forecast date, time, and astronomy quantities.
real(kind=kind_phys), parameter czlimt
this module defines fortran unit numbers for input/output data files for the ncep gfs model.
Definition iounitdef.f:53
This module sets up astronomy quantities for solar radiation calculations.