CCPP SciDoc v7.0.x  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!...................................
308 end subroutine sol_init
309!-----------------------------------
310
326!-----------------------------------
327 subroutine sol_update &
328 & ( jdate,kyear,deltsw,deltim,lsol_chg, me, & ! --- inputs
329 & slag, sdec, cdec, solcon, con_pi, errmsg, errflg & ! --- outputs
330 & )
331
332! =================================================================== !
333! !
334! sol_update computes solar parameters at forecast time !
335! !
336! inputs: !
337! jdate(8)- ncep absolute date and time at fcst time !
338! (yr, mon, day, t-zone, hr, min, sec, mil-sec) !
339! kyear - usually kyear=jdate(1). if not, it is for hindcast mode,!
340! and it is usually the init cond time and serves as the !
341! upper limit of data can be used. !
342! deltsw - time duration in seconds per sw calculation !
343! deltim - timestep in seconds !
344! lsol_chg- logical flags for change solar constant !
345! me - print message control flag !
346! !
347! outputs: !
348! slag - equation of time in radians !
349! sdec, cdec - sin and cos of the solar declination angle !
350! solcon - sun-earth distance adjusted solar constant (w/m2) !
351! errmsg - CCPP error message !
352! errflg - CCPP error flag !
353! !
354! !
355! module variable: !
356! solc0 - solar constant (w/m**2) not adjusted by earth-sun dist !
357! isolflg - solar constant control flag !
358! = 0: use the old fixed solar constant !
359! =10: use the new fixed solar constant !
360! = 1: use noaa ann-mean tsi tbl abs-scale with cycle apprx !
361! = 2: use noaa ann-mean tsi tbl tim-scale with cycle apprx !
362! = 3: use cmip5 ann-mean tsi tbl tim-scale with cycle apprx!
363! = 4: use cmip5 mon-mean tsi tbl tim-scale with cycle apprx!
364! solar_fname-external solar constant data table !
365! sindec - sine of the solar declination angle !
366! cosdec - cosine of the solar declination angle !
367! anginc - solar angle increment per iteration for cosz calc !
368! nstp - total number of zenith angle iterations !
369! smon_sav- saved monthly solar constants (isolflg=4 only) !
370! iyr_sav - saved year of data previously used !
371! !
372! usage: call sol_update !
373! !
374! subprograms called: solar, prtime !
375! !
376! external functions called: iw3jdn !
377! !
378! =================================================================== !
379!
380 implicit none
381
382! --- input:
383 integer, intent(in) :: jdate(:), kyear, me
384 logical, intent(in) :: lsol_chg
385
386 real (kind=kind_phys), intent(in) :: deltsw, deltim, con_pi
387
388! --- output:
389 real (kind=kind_phys), intent(out) :: slag, sdec, cdec, solcon
390 character(len=*), intent(out) :: errmsg
391 integer, intent(out) :: errflg
392
393! --- locals:
394 real (kind=kind_phys), parameter :: hrday = 1.0/24.0 ! frc day/hour
395 real (kind=kind_phys), parameter :: minday= 1.0/1440.0 ! frc day/minute
396 real (kind=kind_phys), parameter :: secday= 1.0/86400.0 ! frc day/second
397
398 real (kind=kind_phys) :: smean, solc1, dtswh, smon(12)
399 real (kind=kind_phys) :: fjd, fjd1, dlt, r1, alp
400
401 integer :: jd, jd1, iyear, imon, iday, ihr, imin, isec
402 integer :: iw3jdn
403 integer :: i, iyr, iyr1, iyr2, jyr, nn, nswr, icy1, icy2, icy
404
405 logical :: file_exist
406 character :: cline*60
407!
408!===> ... begin here
409!
410! Initialize the CCPP error handling variables
411 errmsg = ''
412 errflg = 0
413
414! --- ... forecast time
415 iyear = jdate(1)
416 imon = jdate(2)
417 iday = jdate(3)
418 ihr = jdate(5)
419 imin = jdate(6)
420 isec = jdate(7)
421
422 if ( lsol_chg ) then ! get solar constant from data table
423
424 if ( iyr_sav == iyear ) then ! same year, no new reading necessary
425 if ( isolflg==4 ) then
426 solc0 = smon_sav(imon)
427 endif
428 else ! need to read in new data
429 iyr_sav = iyear
430
431! --- ... check to see if the solar constant data file existed
432
433 inquire (file=solar_fname, exist=file_exist)
434 if ( .not. file_exist ) then
435 errflg = 1
436 errmsg = "ERROR(radiation_astronomy): solar constant file"//&
437 & " not found"
438 return
439 else
440 iyr = iyear
441
442 close(niradsf)
443 open (niradsf,file=solar_fname,form='formatted', &
444 & status='old')
445 rewind niradsf
446
447 read (niradsf, * ) iyr1,iyr2,icy1,icy2,smean,cline(1:60)
448! read (NIRADSF, 24) iyr1,iyr2,icy1,icy2,smean,cline
449! 24 format(4i5,f8.2,a60)
450
451 if ( me == 0 ) then
452 print *,' Updating solar constant with cycle approx'
453 print *,' Opened solar constant data file: ',solar_fname
454!check print *, iyr1, iyr2, icy1, icy2, smean, cline
455 endif
456
457! --- ... check if there is a upper year limit put on the data table
458
459! if ( iyear /= kyear ) then
460! icy = icy1 - iyr1 + 1 ! range of the earlest cycle in data table
461! if ( kyear-iyr1 < icy ) then ! need data range at least icy years
462 ! to perform cycle approximation
463! if ( me == 0 ) then
464! print *,' *** the requested year',iyear,' and upper',&
465! & 'limit',kyear,' do not fit the range of data ', &
466! & 'table of iyr1, iyr2 =',iyr1,iyr2
467! print *,' USE FIXED SOLAR CONSTANT=',con_solr
468! endif
469! solc0 = con_solr
470! isolflg = 10
471
472! elseif ( kyear < iyr2 ) then
473
474! --- ... because the usage limit put on the historical data table,
475! skip those unused data records at first
476
477! i = iyr2
478! Lab_dowhile0 : do while ( i > kyear )
479! read (NIRADSF,26) jyr, solc1
480! 26 format(i4,f10.4)
481! read (NIRADSF,*) jyr, solc1
482! i = i - 1
483! enddo Lab_dowhile0
484
485! iyr2 = kyear ! next record will serve the upper limit
486
487! endif ! end if_kyear_block
488! endif ! end if_iyear_block
489
490! --- ... checking the cycle range
491
492 if ( iyr < iyr1 ) then
493 icy = icy1 - iyr1 + 1 ! range of the earlest cycle in data table
494 lab_dowhile1 : do while ( iyr < iyr1 )
495 iyr = iyr + icy
496 enddo lab_dowhile1
497
498 if ( me == 0 ) then
499 print *,' *** Year',iyear,' out of table range!', &
500 & iyr1, iyr2
501 print *,' Using the closest-cycle year (',iyr,')'
502 endif
503 elseif ( iyr > iyr2 ) then
504 icy = iyr2 - icy2 + 1 ! range of the latest cycle in data table
505 lab_dowhile2 : do while ( iyr > iyr2 )
506 iyr = iyr - icy
507 enddo lab_dowhile2
508
509 if ( me == 0 ) then
510 print *,' *** Year',iyear,' out of table range!', &
511 & iyr1, iyr2
512 print *,' Using the closest-cycle year (',iyr,')'
513 endif
514 endif
515
516! --- ... locate the right record for the year of data
517
518 if ( isolflg < 4 ) then ! use annual mean data tables
519 i = iyr2
520 lab_dowhile3 : do while ( i >= iyr1 )
521! read (NIRADSF,26) jyr, solc1
522! 26 format(i4,f10.4)
523 read (niradsf,*) jyr, solc1
524
525 if ( i == iyr .and. iyr == jyr ) then
526 solc0 = smean + solc1
527
528 if (me == 0) then
529 print *,' CHECK: Solar constant data used for year',&
530 & iyr, solc1, solc0
531 endif
532 exit lab_dowhile3
533 else
534!check if(me == 0) print *,' Skip solar const data for yr',i
535 i = i - 1
536 endif
537 enddo lab_dowhile3
538 elseif ( isolflg == 4 ) then ! use monthly mean data tables
539 i = iyr2
540 lab_dowhile4 : do while ( i >= iyr1 )
541! read (NIRADSF,26) jyr, smon(:)
542! 26 format(i4,12f10.4)
543 read (niradsf,*) jyr, smon(1:12)
544
545 if ( i == iyr .and. iyr == jyr ) then
546 do nn = 1, 12
547 smon_sav(nn) = smean + smon(nn)
548 enddo
549 solc0 = smean + smon(imon)
550
551 if (me == 0) then
552 print *,' CHECK: Solar constant data used for year',&
553 & iyr,' and month',imon
554 endif
555 exit lab_dowhile4
556 else
557!check if(me == 0) print *,' Skip solar const data for yr',i
558 i = i - 1
559 endif
560 enddo lab_dowhile4
561 endif ! end if_isolflg_block
562
563 close ( niradsf )
564 endif ! end if_file_exist_block
565
566 endif ! end if_iyr_sav_block
567 endif ! end if_lsol_chg_block
568
569! --- ... calculate forecast julian day and fraction of julian day
570
571 jd1 = iw3jdn(iyear,imon,iday)
572
573! --- ... unlike in normal applications, where day starts from 0 hr,
574! in astronomy applications, day stats from noon.
575
576 if (ihr < 12) then
577 jd1 = jd1 - 1
578 fjd1= 0.5 + float(ihr)*hrday + float(imin)*minday &
579 & + float(isec)*secday
580 else
581 fjd1= float(ihr - 12)*hrday + float(imin)*minday &
582 & + float(isec)*secday
583 endif
584
585 fjd1 = fjd1 + jd1
586
587 jd = int(fjd1)
588 fjd = fjd1 - jd
589
591 call solar &
592! --- inputs:
593 & ( jd, fjd, con_pi, &
594! --- outputs:
595 & r1, dlt, alp &
596 & )
597
598! --- ... calculate sun-earth distance adjustment factor appropriate to date
599 solcon = solc0 / (r1*r1)
600
601 slag = sollag
602 sdec = sindec
603 cdec = cosdec
604
605! --- ... diagnostic print out
606
607 if (me == 0) then
608
610 call prtime &
611! --- inputs:
612 & ( jd, fjd, dlt, alp, r1, solcon )
613! --- outputs: ( none )
614
615 endif
616
617! --- ... setting up calculation parameters used by subr coszmn
618
619 nswr = max(1, nint(deltsw/deltim)) ! number of mdl t-step per sw call
620 dtswh = deltsw / f3600 ! time length in hours
621
622! if ( deltsw >= f3600 ) then ! for longer sw call interval
623! nn = max(6, min(12, nint(f3600/deltim) )) ! num of calc per hour
624! nstp = nint(dtswh) * nn + 1 ! num of calc per sw call
625! else ! for shorter sw sw call interval
626! nstp = max(2, min(20, nswr)) + 1
627!! nn = nint( float(nstp-1)/dtswh )
628! endif
629
630! anginc = pid12 * dtswh / float(nstp-1) ! solar angle inc during each calc step
631
632 nstp = max(6, nswr)
633 anginc = pid12 * dtswh / float(nstp)
634
635 if ( me == 0 ) then
636 print *,' for cosz calculations: nswr,deltim,deltsw,dtswh =', &
637 & nswr,deltim,deltsw,dtswh,' anginc,nstp =',anginc,nstp
638 endif
639
640! if (me == 0) print*,'in sol_update completed sr solar'
641!
642!...................................
643 end subroutine sol_update
644!-----------------------------------
645
654!-----------------------------------
655 subroutine solar &
656 & ( jd, fjd, con_pi, & ! --- inputs
657 & r1, dlt, alp & ! --- outputs
658 & )
659
660! =================================================================== !
661! !
662! solar computes radius vector, declination and right ascension of !
663! sun, and equation of time. !
664! !
665! inputs: !
666! jd - julian day !
667! fjd - fraction of the julian day !
668! !
669! outputs: !
670! r1 - earth-sun radius vector !
671! dlt - declination of sun in radians !
672! alp - right ascension of sun in radians !
673! !
674! module variables: !
675! sollag - equation of time in radians !
676! sindec - sine of declination angle !
677! cosdec - cosine of declination angle !
678! !
679! usage: call solar !
680! !
681! external subroutines called: none !
682! !
683! =================================================================== !
684!
685 implicit none
686
687! --- inputs:
688 real (kind=kind_phys), intent(in) :: fjd, con_pi
689 integer, intent(in) :: jd
690
691! --- outputs:
692 real (kind=kind_phys), intent(out) :: r1, dlt, alp
693
694! --- locals:
695 real (kind=kind_phys), parameter :: cyear = 365.25 ! days of year
696 real (kind=kind_phys), parameter :: ccr = 1.3e-6 ! iteration limit
697 real (kind=kind_phys), parameter :: tpp = 1.55 ! days between epoch and
698 ! perihelion passage of 1900
699 real (kind=kind_phys), parameter :: svt6 = 78.035 ! days between perihelion passage
700 ! and march equinox of 1900
701 integer, parameter :: jdor = 2415020 ! jd of epoch which is january
702 ! 0, 1900 at 12 hours ut
703
704 real (kind=kind_phys) :: dat, t1, year, tyear, ec, angin, ador, &
705 & deleqn, sni, tini, er, qq, e1, ep, cd, eq, date, em, &
706 & cr, w1, tst, sun
707
708 integer :: jdoe, iter
709
710!===> ... begin here
711
712! --- ... computes time in julian centuries after epoch
713
714 t1 = float(jd - jdor) / 36525.0
715
716! --- ... computes length of anomalistic and tropical years (minus 365 days)
717
718 year = 0.25964134e0 + 0.304e-5 * t1
719 tyear= 0.24219879e0 - 0.614e-5 * t1
720
721! --- ... computes orbit eccentricity and angle of earth's inclination from t
722
723 ec = 0.01675104e0 - (0.418e-4 + 0.126e-6 * t1) * t1
724 angin= 23.452294e0 - (0.0130125e0 + 0.164e-5 * t1) * t1
725
726 ador = jdor
727 jdoe = ador + (svt6 * cyear) / (year - tyear)
728
729! --- ... deleqn is updated svt6 for current date
730
731 deleqn= float(jdoe - jd) * (year - tyear) / cyear
732 year = year + 365.0
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
737
738! --- ... determine true anomaly at equinox
739
740 e1 = 1.0
741 cd = 1.0
742 iter = 0
743
744 lab_do_1 : do while ( cd > ccr )
745
746 ep = e1 - (e1 - ec*sin(e1) - qq) / (1.0 - ec*cos(e1))
747 cd = abs(e1 - ep)
748 e1 = ep
749 iter = iter + 1
750
751 if (iter > 10) then
752 write(6,*) ' ITERATION COUNT FOR LOOP 32 =', iter
753 write(6,*) ' E, EP, CD =', e1, ep, cd
754 exit lab_do_1
755 endif
756
757 enddo lab_do_1
758
759 eq = 2.0 * atan( er * tan( 0.5*e1 ) )
760
761! --- ... date is days since last perihelion passage
762
763 dat = float(jd - jdor) - tpp + fjd
764 date = mod(dat, year)
765
766! --- ... solve orbit equations by newton's method
767
768 em = tpi * date / year
769 e1 = 1.0
770 cr = 1.0
771 iter = 0
772
773 lab_do_2 : do while ( cr > ccr )
774
775 ep = e1 - (e1 - ec*sin(e1) - em) / (1.0 - ec*cos(e1))
776 cr = abs(e1 - ep)
777 e1 = ep
778 iter = iter + 1
779
780 if (iter > 10) then
781 write(6,*) ' ITERATION COUNT FOR LOOP 31 =', iter
782 exit lab_do_2
783 endif
784
785 enddo lab_do_2
786
787 w1 = 2.0 * atan( er * tan( 0.5*e1 ) )
788
789 r1 = 1.0 - ec*cos(e1)
790
791 sindec = sni * sin(w1 - eq)
792 cosdec = sqrt( 1.0 - sindec*sindec )
793
794 dlt = asin( sindec )
795 alp = asin( tan(dlt)*tini )
796
797 tst = cos( w1 - eq )
798 if (tst < 0.0) alp = con_pi - alp
799 if (alp < 0.0) alp = alp + tpi
800
801 sun = tpi * (date - deleqn) / year
802 if (sun < 0.0) sun = sun + tpi
803 sollag = sun - alp - 0.03255e0
804!
805!...................................
806 end subroutine solar
807!-----------------------------------
808
822!-----------------------------------
823 subroutine coszmn &
824 & ( xlon,sinlat,coslat,solhr, im, me, & ! --- inputs
825 & coszen, coszdg & ! --- outputs
826 & )
827
828! =================================================================== !
829! !
830! coszmn computes mean cos solar zenith angle over sw calling interval !
831! !
832! inputs: !
833! xlon (IM) - grids' longitudes in radians, work both on zonal !
834! 0->2pi and -pi->+pi arrangements !
835! sinlat(IM) - sine of the corresponding latitudes !
836! coslat(IM) - cosine of the corresponding latitudes !
837! solhr - time after 00z in hours !
838! IM - num of grids in horizontal dimension !
839! me - print message control flag !
840! !
841! outputs: !
842! coszen(IM) - average of cosz for daytime only in sw call interval
843! coszdg(IM) - average of cosz over entire sw call interval !
844! !
845! module variables: !
846! sollag - equation of time !
847! sindec - sine of the solar declination angle !
848! cosdec - cosine of the solar declination angle !
849! anginc - solar angle increment per iteration for cosz calc !
850! nstp - total number of zenith angle iterations !
851! !
852! usage: call comzmn !
853! !
854! external subroutines called: none !
855! !
856! =================================================================== !
857!
858 implicit none
859
860! --- inputs:
861 integer, intent(in) :: im, me
862
863 real (kind=kind_phys), intent(in) :: sinlat(:), coslat(:), &
864 & xlon(:), solhr
865
866! --- outputs:
867 real (kind=kind_phys), intent(out) :: coszen(:), coszdg(:)
868
869! --- locals:
870 real (kind=kind_phys) :: coszn, cns, solang, rstp
871
872 integer :: istsun(im), i, it, j, lat
873
874!===> ... begin here
875
876 solang = pid12 * (solhr - f12) ! solar angle at present time
877 rstp = 1.0 / float(nstp)
878
879 do i = 1, im
880 coszen(i) = 0.0
881 istsun(i) = 0
882 enddo
883
884 do it = 1, nstp
885 cns = solang + (float(it)-0.5)*anginc + sollag
886 do i = 1, im
887 coszn = sindec * sinlat(i) + cosdec * coslat(i) &
888 & * cos(cns+xlon(i))
889 coszen(i) = coszen(i) + max(0.0, coszn)
890 if (coszn > czlimt) istsun(i) = istsun(i) + 1
891 enddo
892 enddo
893
894! --- ... compute time averages
895
896 do i = 1, im
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)
900 endif
901 enddo
902!
903!...................................
904 end subroutine coszmn
905!-----------------------------------
906
916!-----------------------------------
917 subroutine prtime &
918 & ( jd, fjd, dlt, alp, r1, solc & ! --- inputs
919 & ) ! --- outputs: ( none )
920
921! =================================================================== !
922! !
923! prtime prints out forecast date, time, and astronomy quantities. !
924! !
925! inputs: !
926! jd - forecast julian day !
927! fjd - forecast fraction of julian day !
928! dlt - declination angle of sun in radians !
929! alp - right ascension of sun in radians !
930! r1 - earth-sun radius vector in meter !
931! solc - solar constant in w/m^2 !
932! !
933! outputs: ( none ) !
934! !
935! module variables: !
936! sollag - equation of time in radians !
937! !
938! usage: call prtime !
939! !
940! external subroutines called: w3fs26 !
941! !
942! =================================================================== !
943!
944 implicit none
945
946! --- inputs:
947 integer, intent(in) :: jd
948
949 real (kind=kind_phys), intent(in) :: fjd, dlt, alp, r1, solc
950
951! --- outputs: ( none )
952
953! --- locals:
954 real (kind=kind_phys), parameter :: sixty = 60.0
955
956 character(LEN=1), parameter :: sign = '-'
957 character(LEN=1), parameter :: sigb = ' '
958
959 character(LEN=1) :: dsig
960 character(LEN=4) :: month(12)
961
962 data month / 'JAN.','FEB.','MAR.','APR.','MAY ','JUNE', &
963 & 'JULY','AUG.','SEP.','OCT.','NOV ','DEC.' /
964
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, &
968 & asec, eqt, eqsec
969
970!===> ... begin here
971
972! --- ... get forecast hour and minute from fraction of julian day
973
974 if (fjd >= 0.5) then
975 jda = jd + 1
976 mfjd= nint( fjd*1440.0 )
977 ihr = mfjd / 60 - 12
978 xmin= float(mfjd) - (ihr + 12)*sixty
979 else
980 jda = jd
981 mfjd= nint( fjd*1440.0 )
982 ihr = mfjd / 60 + 12
983 xmin= float(mfjd) - (ihr - 12)*sixty
984 endif
985
986! --- ... get forecast year, month, and day from julian day
987
988 call w3fs26(jda, iyear,imon,iday, idaywk,idayyr)
989
990! -- ... compute solar parameters
991
992 dltd = degrad * dlt
993 ltd = dltd
994 dltm = sixty * (abs(dltd) - abs(float(ltd)))
995 ltm = dltm
996 dlts = sixty * (dltm - float(ltm))
997
998 if ((dltd < 0.0) .and. (ltd == 0.0)) then
999 dsig = sign
1000 else
1001 dsig = sigb
1002 endif
1003
1004 halp = 6.0 * alp / hpi
1005 ihalp= halp
1006 ymin = abs(halp - float(ihalp)) * sixty
1007 iyy = ymin
1008 asec = (ymin - float(iyy)) * sixty
1009
1010 eqt = 228.55735 * sollag
1011 eqsec= sixty * eqt
1012
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)
1016
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')
1020
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)'//)
1026
1027!
1028!...................................
1029 end subroutine prtime
1030!-----------------------------------
1031
1032!
1033!...........................................!
1034 end module module_radiation_astronomy !
1036!===========================================!
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.