FV3 Bundle
astronomy.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 
20 !> \author Fei Liu <Fei.Liu@noaa.gov>
21 !! \link http://www.gfdl.noaa.gov/fms-cgi-bin/cvsweb.cgi/FMS/ \endlink
22 !!
23 !! \brief astronomy_mod provides astronomical variables for use
24 !! by other modules within fms. The only currently used interface is
25 !! for determination of astronomical values needed by the shortwave
26 !! radiation packages.
27 !!
28 !! Modules Included:
29 !!
30 !! <table>
31 !! <tr>
32 !! <th>Module Name</th>
33 !! <th>Functions Included</th>
34 !! </tr>
35 !! <tr>
36 !! <td>fms_mod</td>
37 !! <td>open_namelist_file, fms_init, mpp_pe, mpp_root_pe, stdlog,
38 !! file_exist, write_version_number, check_nml_error, error_mesg,
39 !! FATAL, NOTE, WARNING, close_file</td>
40 !! </tr>
41 !! <tr>
42 !! <td>time_manager_mod</td>
43 !! <td>time_type, set_time, get_time, get_date_julian, set_date_julian,
44 !! set_date, length_of_year, time_manager_init, operator(-),
45 !! operator(+), operator( // ), operator(<)</td>
46 !! </tr>
47 !! <tr>
48 !! <td>constants_mod</td>
49 !! <td>constants_init, PI</td>
50 !! </tr>
51 !! <tr>
52 !! <td>mpp_mod</td>
53 !! <td>input_nml_file</td>
54 !! </tr>
55 !! </table>
57 
58 
59 use fms_mod, only: open_namelist_file, fms_init, &
60  mpp_pe, mpp_root_pe, stdlog, &
61  file_exist, write_version_number, &
63  fatal, note, warning, close_file
68  operator(-), operator(+), &
69  operator( // ), operator(<)
71 use mpp_mod, only: input_nml_file
72 
73 !--------------------------------------------------------------------
74 
75 implicit none
76 private
77 
78 !---------------------------------------------------------------------
79 !----------- version number for this module --------------------------
80 
81 ! Include variable "version" to be written to log file.
82 #include<file_version.h>
83 
84 
85 !---------------------------------------------------------------------
86 !------- interfaces --------
87 
88 public &
94 
95 !> \page diurnal_solar diurnal_solar Interface
96 !!
97 !! ~~~~~~~~~~{.f90}
98 !! call diurnal_solar (lat, lon, time, cosz, fracday, rrsun, dt_time)
99 !! call diurnal_solar (lat, lon, gmt, time_since_ae, cosz, fracday, rrsun, dt)
100 !! ~~~~~~~~~~
101 !!
102 !! The first option (used in conjunction with time_manager_mod)
103 !! generates the real variables gmt and time_since_ae from the
104 !! time_type input, and then calls diurnal_solar with these real inputs.
105 !!
106 !! The time of day is set by
107 !! ~~~~~~~~~~{.f90}
108 !! real, intent(in) :: gmt
109 !! ~~~~~~~~~~
110 !! The time of year is set by
111 !! ~~~~~~~~~~{.f90}
112 !! real, intent(in) :: time_since_ae
113 !! ~~~~~~~~~~
114 !! with time_type input, both of these are extracted from
115 !! ~~~~~~~~~~{.f90}
116 !! type(time_type), intent(in) :: time
117 !! ~~~~~~~~~~
118 !!
119 !! Separate routines exist within this interface for scalar,
120 !! 1D or 2D input and output fields:
121 !!
122 !! ~~~~~~~~~~{.f90}
123 !! real, intent(in), dimension(:,:) :: lat, lon
124 !! real, intent(in), dimension(:) :: lat, lon
125 !! real, intent(in) :: lat, lon
126 !!
127 !! real, intent(out), dimension(:,:) :: cosz, fracday
128 !! real, intent(out), dimension(:) :: cosz, fracday
129 !! real, intent(out) :: cosz, fracday
130 !! ~~~~~~~~~~
131 !!
132 !! One may also average the output fields over the time interval
133 !! between gmt and gmt + dt by including the optional argument dt (or
134 !! dt_time). dt is measured in radians and must be less than pi
135 !! (1/2 day). This average is computed analytically, and should be
136 !! exact except for the fact that changes in earth-sun distance over
137 !! the time interval dt are ignored. In the context of a diurnal GCM,
138 !! this option should always be employed to insure that the total flux
139 !! at the top of the atmosphere is not modified by time truncation error.
140 !!
141 !! ~~~~~~~~~~{.f90}
142 !! real, intent(in), optional :: dt
143 !! type(time_type), optional :: dt_time
144 !! ~~~~~~~~~~
145 !! (see test.90 for examples of the use of these types)
146 !!
147 !! \param [in] <lat> Latitudes of model grid points [radians]
148 !! \param [in] <lon> Longitudes of model grid points [radians]
149 !! \param [in] <gmt> Time of day at longitude 0.0; midnight = 0.0, one day = 2 * pi [radians]
150 !! \param [in] <time_since_ae> Time of year; autumnal equinox = 0.0, one year = 2 * pi [radians]
151 !! \param [in] <time> Time at which astronomical values are desired (time_type variable) [seconds, days]
152 !! \param [out] <cosz> Cosine of solar zenith angle, set to zero when entire period is in darkness [dimensionless]
153 !! \param [out] <fracday> Daylight fraction of time interval [dimensionless]
154 !! \param [out] <rrsun> Earth-Sun distance (r) relative to semi-major axis of orbital ellipse (a):(a/r)**2 [dimensionless]
155 !! \param [in] <dt> OPTIONAL: Time interval after gmt over which the astronomical variables are to be
156 !! averaged. this produces averaged output rather than instantaneous. [radians], (1 day = 2 * pi)
157 !! \param [in] <dt_time> OPTIONAL: Time interval after gmt over which the astronomical variables are to be
158 !! averaged. this produces averaged output rather than instantaneous. time_type, [days, seconds]
159 !! \param [in] <allow_negative_cosz> Allow negative values for cosz?
160 !! \param [out] <half_day_out> half_day_out
161 interface diurnal_solar
162  module procedure diurnal_solar_2d
163  module procedure diurnal_solar_1d
164  module procedure diurnal_solar_0d
165  module procedure diurnal_solar_cal_2d
166  module procedure diurnal_solar_cal_1d
167  module procedure diurnal_solar_cal_0d
168 end interface
169 
170 !> \page daily_mean_solar daily_mean_solar Interface
171 !!
172 !! ~~~~~~~~~~{.f90}
173 !! call daily_mean_solar (lat, time, cosz, fracday, rrsun)
174 !! call daily_mean_solar (lat, time_since_ae, cosz, fracday, rrsun)
175 !! call daily_mean_solar (lat, time, cosz, solar)
176 !! call daily_mean_solar (lat, time_since_ae, cosz, solar)
177 !! ~~~~~~~~~~
178 !!
179 !! The first option (used in conjunction with time_manager_mod)
180 !! generates the real variable time_since_ae from the time_type
181 !! input time, and then calls daily_mean_solar with this real input
182 !! (option 2). The third and fourth options correspond to the first
183 !! and second and are used with then spectral 2-layer model, where
184 !! only cosz and solar are desired as output. These routines generate
185 !! dummy arguments and then call option 2, where the calculation is done.
186 !!
187 !! The time of year is set by
188 !! ~~~~~~~~~~{.f90}
189 !! real, intent(in) :: time_since_ae
190 !! ~~~~~~~~~~
191 !! With time_type input, the time of year is extracted from
192 !! ~~~~~~~~~~{.f90}
193 !! type(time_type), intent(in) :: time
194 !! ~~~~~~~~~~
195 !!
196 !! Separate routines exist within this interface for scalar,
197 !! 1D or 2D input and output fields:
198 !!
199 !! ~~~~~~~~~~{.f90}
200 !! real, intent(in), dimension(:,:) :: lat
201 !! real, intent(in), dimension(:) :: lat
202 !! real, intent(in) :: lat
203 !!
204 !! real, intent(out), dimension(:,:) :: cosz, fracday
205 !! real, intent(out), dimension(:) :: cosz, fracday
206 !! real, intent(out) :: cosz, fracday
207 !! ~~~~~~~~~~
208 !!
209 !! \param [in] <lat> Latitudes of model grid points [radians]
210 !! \param [in] <time_since_ae> Time of year; autumnal equinox = 0.0, one year = 2 * pi [radians]
211 !! \param [in] <time> Time at which astronomical values are desired (time_type variable) [seconds, days]
212 !! \param [out] <cosz> Cosine of solar zenith angle, set to zero when entire period is in darkness [dimensionless]
213 !! \param [out] <fracday> Daylight fraction of time interval [dimensionless]
214 !! \param [out] <rrsun> Earth-Sun distance (r) relative to semi-major axis of orbital ellipse (a):(a/r)**2 [dimensionless]
215 !! \param [out] <solar> shortwave flux factor: cosine of zenith angle * daylight fraction / (earth-sun distance squared) [dimensionless]
217  module procedure daily_mean_solar_2d
218  module procedure daily_mean_solar_1d
219  module procedure daily_mean_solar_2level
220  module procedure daily_mean_solar_0d
221  module procedure daily_mean_solar_cal_2d
222  module procedure daily_mean_solar_cal_1d
223  module procedure daily_mean_solar_cal_2level
224  module procedure daily_mean_solar_cal_0d
225 end interface
226 
227 !! \page annual_mean_solar annual_mean_solar Interface
228 !!
229 !! ~~~~~~~~~~{.f90}
230 !! call annual_mean_solar (js, je, lat, cosz, solar, fracday, rrsun)
231 !! call annual_mean_solar (lat, cosz, solar)
232 !! ~~~~~~~~~~
233 !!
234 !! The second interface above is used by the spectral 2-layer model,
235 !! which requires only cosz and solar as output arguments, and which
236 !! makes this call during the initialization phase of the model.
237 !! Separate routines exist within this interface for 1D or 2D input
238 !! and output fields:
239 !!
240 !! ~~~~~~~~~~{.f90}
241 !! real, intent(in), dimension(:,:) :: lat
242 !! real, intent(in), dimension(:) :: lat
243 !!
244 !! real, intent(out), dimension(:,:) :: cosz, solar, fracday
245 !! real, intent(out), dimension(:) :: cosz, solar, fracday
246 !! ~~~~~~~~~~
247 !!
248 !! \param [in] <jst> Starting subdomain j indices of data in the physics wiondow being integrated
249 !! \param [in] <jnd> Ending subdomain j indices of data in the physics wiondow being integrated
250 !! \param [in] <lat> Latitudes of model grid points [radians]
251 !! \param [out] <cosz> cosz is the average over the year of the cosine of an effective zenith angle
252 !! that would produce the correct daily solar flux if the sun were fixed at that
253 !! single position for the period of daylight on the given day. in this average,
254 !! the daily mean effective cosz is weighted by the daily mean solar flux. [dimensionless]
255 !! \param [out] <solar> Normalized solar flux, averaged over the year, equal to the product of
256 !! fracday*cosz*rrsun [dimensionless]
257 !! \param [out] <fracday> Daylight fraction calculated so as to make the average flux (solar) equal to the
258 !! product of the flux-weighted avg cosz * this fracday * assumed annual mean avg
259 !! Earth-Sun distance of 1.0. [dimensionless]
260 !! \param [out] <rrsun> Annual mean Earth-Sun distance (r) relative to semi-major axis of orbital ellipse
261 !! (a):(a/r)**2 [dimensionless]
263  module procedure annual_mean_solar_2d
264  module procedure annual_mean_solar_1d
265  module procedure annual_mean_solar_2level
266 end interface
267 
268 !> \page get_period get_period Interface
269 !!
270 !! ~~~~~~~~~~{.f90}
271 !! call get_period (period)
272 !! ~~~~~~~~~~
273 !!
274 !! Separate routines exist within this interface for integer
275 !! and time_type output:
276 !!
277 !! ~~~~~~~~~~{.f90}
278 !! integer, intent(out) :: period
279 !! type(time_type), intent(out) :: period
280 !! ~~~~~~~~~~
281 !!
282 !! \param [out] <period_out> Length of year for calendar in use
283 interface get_period
284  module procedure get_period_time_type, get_period_integer
285 end interface
286 
287 !> \page set_period set_period Interface
288 !!
289 !! ~~~~~~~~~~{.f90}
290 !! call set_period (period_in)
291 !! ~~~~~~~~~~
292 !!
293 !! Separate routines exist within this interface for integer
294 !! and time_type output:
295 !!
296 !! ~~~~~~~~~~{.f90}
297 !! integer, intent(out) :: period_in
298 !! type(time_type), intent(out) :: period_in
299 !! ~~~~~~~~~~
300 !!
301 !! \param [in] <period_in> Length of year for calendar in use
302 interface set_period
303  module procedure set_period_time_type, set_period_integer
304 end interface
305 
306 
307 private &
308  orbit, & ! Called from astronomy_init and set_orbital_parameters
309  r_inv_squared, & ! Called from diurnal_solar, daily_mean_solar and orbit
310  angle, declination, half_day ! called from diurnal_solar and daily_mean_solar
311 ! half_day, orbital_time, & ! called from diurnal_solar and daily_mean_solar
312 ! universal_time ! called from diurnal_solar:
313 
314 !> \page half_day half_day Interface
315 !!
316 !! ~~~~~~~~~~{.f90}
317 !! half_day (latitude, dec) result (h)
318 !! ~~~~~~~~~~
319 !!
320 !! Separate routines exist within this interface for scalar,
321 !! or 2D input and output fields:
322 !!
323 !! ~~~~~~~~~~{.f90}
324 !! real, intent(in), dimension(:,:) :: latitude
325 !! real, intent(in) :: latitude
326 !!
327 !! real, dimension(size(latitude,1),size(latitude,2)) :: h
328 !! real :: h
329 !! ~~~~~~~~~~
330 !!
331 !! \param [in] <latitude> Latitudes of model grid points [radians]
332 !! \param [in] <dec> Solar declination [radians]
333 !! \param [out] <h> Half of the length of daylight at the given latitude and orbital position (dec); value
334 !! ranges between 0 (all darkness) and pi (all daylight) [dimensionless]
335 interface half_day
336  module procedure half_day_2d, half_day_0d
337 end interface
338 
339 
340 !---------------------------------------------------------------------
341 !-------- namelist ---------
342 
343 real :: ecc = 0.01671 !< Eccentricity of Earth's orbit [dimensionless]
344 real :: obliq = 23.439 !< Obliquity [degrees]
345 real :: per = 102.932 !< Longitude of perihelion with respect
346  !! to autumnal equinox in NH [degrees]
347 integer :: period = 0 !< Specified length of year [seconds];
348  !! must be specified to override default
349  !! value given by length_of_year in
350  !! time_manager_mod
351 integer :: day_ae = 23 !< Day of specified autumnal equinox
352 integer :: month_ae = 9 !< Month of specified autumnal equinox
353 integer :: year_ae = 1998 !< Year of specified autumnal equinox
354 integer :: hour_ae = 5 !< Hour of specified autumnal equinox
355 integer :: minute_ae = 37 !< Minute of specified autumnal equinox
356 integer :: second_ae = 0 !< Second of specified autumnal equinox
357 integer :: num_angles = 3600 !< Number of intervals into which the year
358  !! is divided to compute orbital positions
359 
360 
361 namelist /astronomy_nml/ ecc, obliq, per, period, &
362  year_ae, month_ae, day_ae, &
364  num_angles
365 
366 !--------------------------------------------------------------------
367 !------ public data ----------
368 
369 
370 !--------------------------------------------------------------------
371 !------ private data ----------
372 
373 type(time_type) :: autumnal_eq_ref !< time_type variable containing
374  !! specified time of reference
375  !! NH autumnal equinox
376 
377 type(time_type) :: period_time_type !< time_type variable containing
378  !! period of one orbit
379 
380 real, dimension(:), allocatable :: orb_angle !< table of orbital positions (0 to
381  !! 2*pi) as a function of time used
382  !! to find actual orbital position
383  !! via interpolation
384 
385 real :: seconds_per_day=86400. !< seconds in a day
386 real :: deg_to_rad !< conversion from degrees to radians
387 real :: twopi !< 2 *PI
388 logical :: module_is_initialized=.false. !< has the module been initialized ?
389 
390 real, dimension(:,:), allocatable :: &
391  cosz_ann, & !< annual mean cos of zenith angle
392  solar_ann, & !< annual mean solar factor
393  fracday_ann !< annual mean daylight fraction
394 real :: rrsun_ann !< annual mean earth-sun distance
395 logical :: annual_mean_calculated=.false. !< have the annual mean values been calculated?
396 integer :: num_pts = 0 !< count of grid_boxes for which
397  !! annual mean astronomy values have
398  !! been calculated
399 integer :: total_pts !< number of grid boxes owned by the processor
400 
401 !--------------------------------------------------------------------
402 
403 
404 
405  contains
406 
407 
408 
409 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
410 !
411 ! PUBLIC SUBROUTINES
412 !
413 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
414 
415 !> \brief astronomy_init is the constructor for astronomy_mod.
416 !!
417 !! \throw FATAL, "astronomy_mod ecc must be between 0 and 0.99"
418 !! \throw FATAL, "astronomy_mod obliquity must be between -90 and 90 degrees"
419 !! \throw FATAL, "astronomy_mod perihelion must be between 0 and 360 degrees"
420 subroutine astronomy_init (latb, lonb)
422 real, dimension(:,:), intent(in), optional :: latb !< 2d array of model latitudes at cell corners [radians]
423 real, dimension(:,:), intent(in), optional :: lonb !< 2d array of model longitudes at cell corners [radians]
424 
425 !-------------------------------------------------------------------
426 ! local variables:
427 !-------------------------------------------------------------------
428 integer :: unit, ierr, io, seconds, days, jd, id
429 
430 !-------------------------------------------------------------------
431 ! if module has already been initialized, exit.
432 !-------------------------------------------------------------------
433  if (module_is_initialized) return
434 
435 !-------------------------------------------------------------------
436 !> Verify that modules used by this module have been initialized.
437 !-------------------------------------------------------------------
438  call fms_init
439  call time_manager_init
440  call constants_init
441 
442 !-----------------------------------------------------------------------
443 !> Read namelist.
444 !-----------------------------------------------------------------------
445 #ifdef INTERNAL_FILE_NML
446  read (input_nml_file, astronomy_nml, iostat=io)
447  ierr = check_nml_error(io,'astronomy_nml')
448 #else
449  if ( file_exist('input.nml')) then
450  unit = open_namelist_file( )
451  ierr=1; do while (ierr /= 0)
452  read (unit, nml=astronomy_nml, iostat=io, end=10)
453  ierr = check_nml_error(io,'astronomy_nml')
454  end do
455 10 call close_file (unit)
456  endif
457 #endif
458 !---------------------------------------------------------------------
459 !> Write version number and namelist to logfile.
460 !---------------------------------------------------------------------
461  call write_version_number("ASTRONOMY_MOD", version)
462  if (mpp_pe() == mpp_root_pe() ) then
463  unit = stdlog()
464  write (unit, nml=astronomy_nml)
465  endif
466 !--------------------------------------------------------------------
467 !> Be sure input values are within valid ranges.
468 ! QUESTION : ARE THESE THE RIGHT LIMITS ???
469 !---------------------------------------------------------------------
470  if (ecc < 0.0 .or. ecc > 0.99) &
471  call error_mesg ('astronomy_mod', &
472  'ecc must be between 0 and 0.99', fatal)
473  if (obliq < -90. .or. obliq > 90.) &
474  call error_mesg ('astronomy_mod', &
475  'obliquity must be between -90 and 90 degrees', fatal)
476  if (per < 0.0 .or. per > 360.0) &
477  call error_mesg ('astronomy_mod', &
478  'perihelion must be between 0 and 360 degrees', fatal)
479 
480 !----------------------------------------------------------------------
481 !> Set up time-type variable defining specified time of autumnal equinox.
482 !----------------------------------------------------------------------
485 
486 !---------------------------------------------------------------------
487 !> Set up time-type variable defining length of year.
488 !----------------------------------------------------------------------
489  if (period == 0) then
491  call get_time (period_time_type, seconds, days)
492  period = seconds_per_day*days + seconds
493  else
495  endif
496 
497 !---------------------------------------------------------------------
498 !> Define useful module variables.
499 !----------------------------------------------------------------------
500  twopi = 2.*pi
501  deg_to_rad = twopi/360.
502 
503 !---------------------------------------------------------------------
504 !> Call orbit to define table of orbital angles as function of
505 !! orbital time.
506 !----------------------------------------------------------------------
507 ! wfc moved here from orbit
508  allocate ( orb_angle(0:num_angles) )
509  call orbit
510 
511 !--------------------------------------------------------------------
512 !> If annual mean radiation is desired, then latb will be present.
513 !! allocate arrays to hold the needed astronomical factors. define
514 !! the total number of points that the processor is responsible for.
515 !--------------------------------------------------------------------
516  if (present(latb)) then
517  jd = size(latb,2) - 1
518  id = size(lonb,1) - 1
519  allocate (cosz_ann(id, jd))
520  allocate (solar_ann(id, jd))
521  allocate (fracday_ann(id, jd))
522  total_pts = jd*id
523  endif
524 
525 !---------------------------------------------------------------------
526 !> Mark the module as initialized.
527 !---------------------------------------------------------------------
528  module_is_initialized=.true.
529 
530 !---------------------------------------------------------------------
531 
532 end subroutine astronomy_init
533 
534 
535 !> \brief get_period_integer returns the length of the year as an
536 !! integer number of seconds.
537 !!
538 !! \throw FATAL, "astronomy_mod module has not been initialized"
539 subroutine get_period_integer (period_out)
541 integer, intent(out) :: period_out !< Length of year [seconds]
542 
543 !--------------------------------------------------------------------
544 ! local variables:
545 
546  integer :: seconds, days
547 
548 !---------------------------------------------------------------------
549 ! exit if module has not been initialized.
550 !---------------------------------------------------------------------
551  if (.not. module_is_initialized) &
552  call error_mesg ( 'astronomy_mod', &
553  ' module has not been initialized', fatal)
554 
555 !--------------------------------------------------------------------
556 ! define length of year in seconds.
557 !--------------------------------------------------------------------
558  call get_time (period_time_type, seconds, days)
559  period_out = seconds_per_day*days + seconds
560 
561 
562 end subroutine get_period_integer
563 
564 !> \brief get_period_time_type returns the length of the year as a time_type
565 !! variable.
566 !!
567 !! \throw FATAL, "astronomy_mod module has not been initialized"
568 subroutine get_period_time_type (period_out)
570 type(time_type), intent(inout) :: period_out !< Length of year as time_type variable
571 
572 !---------------------------------------------------------------------
573 ! exit if module has not been initialized.
574 !---------------------------------------------------------------------
575  if (.not. module_is_initialized) &
576  call error_mesg ( 'astronomy_mod', &
577  ' module has not been initialized', fatal)
578 
579 !--------------------------------------------------------------------
580 ! define length of year as a time_type variable.
581 !--------------------------------------------------------------------
582  period_out = period_time_type
583 
584 end subroutine get_period_time_type
585 
586 !> \brief set_period_integer saves as the input length of the year (an
587 !! integer) in a time_type module variable.
588 !!
589 !! \throw FATAL, "astronomy_mod module has not been initialized"
590 subroutine set_period_integer (period_in)
592 integer, intent(in) :: period_in !< Length of year as a time_type
593 
594 !---------------------------------------------------------------------
595 ! exit if module has not been initialized.
596 !---------------------------------------------------------------------
597  if (.not. module_is_initialized) &
598  call error_mesg ( 'astronomy_mod', &
599  ' module has not been initialized', fatal)
600 
601 !---------------------------------------------------------------------
602 ! define time_type variable defining the length of year from the
603 ! input value (integer seconds).
604 !---------------------------------------------------------------------
605  period_time_type = set_time(period_in, 0)
606 
607 end subroutine set_period_integer
608 
609 
610 !> \brief Set_period_time_type saves the length of the year (input as a
611 !! time_type variable) into a time_type module variable.
612 !!
613 !! \throw FATAL, "astronomy_mod module has not been initialized"
614 subroutine set_period_time_type(period_in)
616 type(time_type), intent(in) :: period_in
617 
618 !---------------------------------------------------------------------
619 ! exit if module has not been initialized.
620 !---------------------------------------------------------------------
621  if (.not. module_is_initialized) &
622  call error_mesg ( 'astronomy_mod', &
623  ' module has not been initialized', fatal)
624 
625 !---------------------------------------------------------------------
626 ! define time_type variable defining the length of year from the
627 ! input value (time_type).
628 !---------------------------------------------------------------------
629  period_time_type = period_in
630 
631 
632 end subroutine set_period_time_type
633 
634 !> \brief set_orbital_parameters saves the input values of eccentricity,
635 !! obliquity and perihelion time as module variables for use by
636 !! astronomy_mod.
637 !!
638 !! \throw FATAL, "astronomy_mod module has not been initialized"
639 !! \throw FATAL, "astronomy_mod ecc must be between 0 and 0.99"
640 !! \throw FATAL, "astronomy_mod obliquity must be between -90. and 90. degrees"
641 !! \throw FATAL, "astronomy_mod perihelion must be between 0.0 and 360. degrees"
642 subroutine set_orbital_parameters (ecc_in, obliq_in, per_in)
644 real, intent(in) :: ecc_in !< Eccentricity of orbital ellipse [dimensionless]
645 real, intent(in) :: obliq_in !< Obliquity [degrees]
646 real, intent(in) :: per_in !< Longitude of perihelion with respect to autumnal
647  !! equinox in northern hemisphere [degrees]
648 
649 !---------------------------------------------------------------------
650 ! exit if module has not been initialized.
651 !---------------------------------------------------------------------
652  if (.not. module_is_initialized) &
653  call error_mesg ( 'astronomy_mod', &
654  ' module has not been initialized', fatal)
655 
656 !--------------------------------------------------------------------
657 ! be sure input values are within valid ranges.
658 ! QUESTION : ARE THESE THE RIGHT LIMITS ???
659 !---------------------------------------------------------------------
660  if (ecc_in < 0.0 .or. ecc_in > 0.99) &
661  call error_mesg ('astronomy_mod', &
662  'ecc must be between 0 and 0.99', fatal)
663  if (obliq_in < -90.0 .or. obliq > 90.0) &
664  call error_mesg ('astronomy_mod', &
665  'obliquity must be between -90. and 90. degrees', fatal)
666  if (per_in < 0.0 .or. per_in > 360.0) &
667  call error_mesg ('astronomy_mod', &
668  'perihelion must be between 0.0 and 360. degrees', fatal)
669 
670 !---------------------------------------------------------------------
671 ! save input values into module variables.
672 !---------------------------------------------------------------------
673  ecc = ecc_in
674  obliq = obliq_in
675  per = per_in
676 
677 !---------------------------------------------------------------------
678 ! call orbit to define table of orbital angles as function of
679 ! orbital time using the input values of parameters just supplied.
680 !----------------------------------------------------------------------
681  call orbit
682 
683 !----------------------------------------------------------------------
684 
685 end subroutine set_orbital_parameters
686 
687 !> \brief get_orbital_parameters retrieves the orbital parameters for use
688 !! by another module.
689 !!
690 !! \throw FATAL, "astronomy_mod module has not been initialized"
691 subroutine get_orbital_parameters (ecc_out, obliq_out, per_out)
693 !-------------------------------------------------------------------
694 ! get_orbital_parameters retrieves the orbital parameters for use
695 ! by another module.
696 !--------------------------------------------------------------------
697 
698 real, intent(out) :: ecc_out !< Eccentricity of orbital ellipse [dimensionless]
699 real, intent(out) :: obliq_out !< Obliquity [degrees]
700 real, intent(out) :: per_out !< Longitude of perihelion with respect to autumnal
701  !! equinox in northern hemisphere [degrees]
702 
703 !---------------------------------------------------------------------
704 ! exit if module has not been initialized.
705 !---------------------------------------------------------------------
706  if (.not. module_is_initialized) &
707  call error_mesg ( 'astronomy_mod', &
708  ' module has not been initialized', fatal)
709 
710 !--------------------------------------------------------------------
711 ! fill the output arguments with the eccentricity, obliquity and
712 ! perihelion angle.
713 !--------------------------------------------------------------------
714  ecc_out = ecc
715  obliq_out = obliq
716  per_out = per
717 
718 
719 end subroutine get_orbital_parameters
720 
721 
722 !> \brief set_ref_date_of_ae provides a means of specifying the reference
723 !! date of the NH autumnal equinox for a particular year.
724 !!
725 !! \details set_ref_date_of_ae provides a means of specifying the reference
726 !! date of the NH autumnal equinox for a particular year. It is only
727 !! used if calls are made to the calandar versions of the routines
728 !! diurnal_solar and daily_mean_solar. If the NOLEAP calendar is
729 !! used, then the date of autumnal equinox will be the same every
730 !! year. If JULIAN is used, then the date of autumnal equinox will
731 !! return to the same value every 4th year.
732 !!
733 !! \param [in] <day_in> Day of reference autumnal equinox
734 !! \param [in] <month_in> Month of reference autumnal equinox
735 !! \param [in] <year_in> Year of reference autumnal equinox
736 !! \param [out] <second_in> OPTIONAL: Second of reference autumnal equinox
737 !! \param [out] <minute_in> OPTIONAL: Minute of reference autumnal equinox
738 !! \param [out] <hour_in> OPTIONAL: Hour of reference autumnal equinox
739 !!
740 !! \throw FATAL, "astronomy_mod module has not been initialized"
741 subroutine set_ref_date_of_ae (day_in,month_in,year_in, &
742  second_in,minute_in,hour_in)
744 integer, intent(in) :: day_in, month_in, year_in
745 integer, intent(in), optional :: second_in, minute_in, hour_in
746 
747 !---------------------------------------------------------------------
748 ! exit if module has not been initialized.
749 !---------------------------------------------------------------------
750  if (.not. module_is_initialized) &
751  call error_mesg ( 'astronomy_mod', &
752  ' module has not been initialized', fatal)
753 
754 !--------------------------------------------------------------------
755 ! save the input time of ae specification into a time_type module
756 ! variable autumnal_eq_ref.
757 !--------------------------------------------------------------------
758  day_ae = day_in
759  month_ae = month_in
760  year_ae = year_in
761 
762  if (present(second_in)) then
763  second_ae = second_in
764  minute_ae = minute_in
765  hour_ae = hour_in
766  else
767  second_ae = 0
768  minute_ae = 0
769  hour_ae = 0
770  endif
771 
774 
775 !---------------------------------------------------------------------
776 
777 
778 end subroutine set_ref_date_of_ae
779 
780 
781 !> \brief get_ref_date_of_ae retrieves the reference date of the autumnal
782 !! equinox as integer variables.
783 !!
784 !! \throw FATAL, "astronomy_mod module has not been initialized"
785 !!
786 !! \param [out] <day_out> Day of reference autumnal equinox
787 !! \param [out] <month_out> Month of reference autumnal equinox
788 !! \param [out] <year_out> Year of reference autumnal equinox
789 !! \param [out] <second_out> Second of reference autumnal equinox
790 !! \param [out] <minute_out> Minute of reference autumnal equinox
791 !! \param [out] <hour_out> Hour of reference autumnal equinox
792 subroutine get_ref_date_of_ae (day_out,month_out,year_out,&
793  second_out,minute_out,hour_out)
795 integer, intent(out) :: day_out, month_out, year_out, &
796  second_out, minute_out, hour_out
797 
798 !---------------------------------------------------------------------
799 ! exit if module has not been initialized.
800 !---------------------------------------------------------------------
801  if (.not. module_is_initialized) &
802  call error_mesg ( 'astronomy_mod', &
803  ' module has not been initialized', fatal)
804 
805 !---------------------------------------------------------------------
806 ! fill the output fields with the proper module data.
807 !---------------------------------------------------------------------
808  day_out = day_ae
809  month_out = month_ae
810  year_out = year_ae
811  second_out = second_ae
812  minute_out = minute_ae
813  hour_out = hour_ae
814 
815 
816 end subroutine get_ref_date_of_ae
817 
818 
819 !> \brief diurnal_solar_2d returns 2d fields of cosine of zenith angle,
820 !! daylight fraction and earth-sun distance at the specified latitudes,
821 !! longitudes and time. These values may be instantaneous or averaged
822 !! over a specified time interval.
823 !!
824 !! \param [in] <lat> Latitudes of model grid points
825 !! \param [in] <lon> Longitudes of model grid points
826 !! \param [in] <gmt> Time of day at longitude 0.0; midnight = 0.0, one day = 2 * pi
827 !! \param [in] <time_since_ae> Time of year; autumnal equinox = 0.0, one year = 2 * pi
828 !! \param [out] <cosz> Cosine of solar zenith angle
829 !! \param [out] <fracday> Daylight fraction of time interval
830 !! \param [out] <rrsun> earth-sun distance (r) relative to semi-major axis of orbital ellipse (a):(a/r)**2
831 !! \param [in] <dt> OPTIONAL: time interval after gmt over which the astronomical variables are to be
832 !! averaged. this produces averaged output rather than instantaneous.
833 !! \param [in] <allow_negative_cosz> Allow negative values for cosz?
834 !! \param [out] <half_day_out> half_day_out
835 !!
836 !! \throw FATAL, "astronomy_mod time_since_ae not between 0 and 2pi"
837 !! \throw FATAL, "astronomy_mod gmt not between 0 and 2pi"
838 subroutine diurnal_solar_2d (lat, lon, gmt, time_since_ae, cosz, &
839  fracday, rrsun, dt, allow_negative_cosz, &
840  half_day_out)
842 real, dimension(:,:), intent(in) :: lat, lon
843 real, intent(in) :: gmt, time_since_ae
844 real, dimension(:,:), intent(out) :: cosz, fracday
845 real, intent(out) :: rrsun
846 real, intent(in), optional :: dt
847 logical, intent(in), optional :: allow_negative_cosz
848 real, dimension(:,:), intent(out), optional :: half_day_out
849 
850 
851 !---------------------------------------------------------------------
852 ! local variables
853 
854  real, dimension(size(lat,1),size(lat,2)) :: t, tt, h, aa, bb, &
855  st, stt, sh
856  real :: ang, dec
857  logical :: Lallow_negative
858 
859 !---------------------------------------------------------------------
860 ! local variables
861 !
862 ! t time of day with respect to local noon (2 pi = 1 day)
863 ! [ radians ]
864 ! tt end of averaging period [ radians ]
865 ! h half of the daily period of daylight, centered at noon
866 ! [ radians, -pi --> pi ]
867 ! aa sin(lat) * sin(declination)
868 ! bb cos(lat) * cos(declination)
869 ! st sine of time of day
870 ! stt sine of time of day at end of averaging period
871 ! sh sine of half-day period
872 ! ang position of earth in its orbit wrt autumnal equinox
873 ! [ radians ]
874 ! dec earth's declination [ radians ]
875 !
876 !--------------------------------------------------------------------
877 
878 !--------------------------------------------------------------------
879 ! be sure the time in the annual cycle is legitimate.
880 !---------------------------------------------------------------------
881  if (time_since_ae < 0.0 .or. time_since_ae > twopi) &
882  call error_mesg('astronomy_mod', &
883  'time_since_ae not between 0 and 2pi', fatal)
884 
885 !--------------------------------------------------------------------
886 ! be sure the time at longitude = 0.0 is legitimate.
887 !---------------------------------------------------------------------
888  if (gmt < 0.0 .or. gmt > twopi) &
889  call error_mesg('astronomy_mod', &
890  'gmt not between 0 and 2pi', fatal)
891 
892 !---------------------------------------------------------------------
893 !> define the orbital angle (location in year), solar declination and
894 !! earth sun distance factor. use functions contained in this module.
895 !---------------------------------------------------------------------
896  ang = angle(time_since_ae)
897  dec = declination(ang)
898  rrsun = r_inv_squared(ang)
899 
900 !---------------------------------------------------------------------
901 !> define terms needed in the cosine zenith angle equation.
902 !--------------------------------------------------------------------
903  aa = sin(lat)*sin(dec)
904  bb = cos(lat)*cos(dec)
905 
906 !---------------------------------------------------------------------
907 !> define local time. force it to be between -pi and pi.
908 !--------------------------------------------------------------------
909  t = gmt + lon - pi
910  where(t >= pi) t = t - twopi
911  where(t < -pi) t = t + twopi
912 
913  lallow_negative = .false.
914  if (present(allow_negative_cosz)) then
915  if (allow_negative_cosz) lallow_negative = .true.
916  end if
917 
918 !---------------------------------------------------------------------
919 !> perform a time integration to obtain cosz and fracday if desired.
920 !! output is valid over the period from t to t + dt.
921 !--------------------------------------------------------------------
922  h = half_day(lat,dec)
923 
924  if ( present(half_day_out) ) then
925  half_day_out = h
926  end if
927 
928  if ( present(dt) ) then ! (perform time averaging)
929  tt = t + dt
930  st = sin(t)
931  stt = sin(tt)
932  sh = sin(h)
933  cosz = 0.0
934 
935  if (.not. lallow_negative) then
936 !-------------------------------------------------------------------
937 ! case 1: entire averaging period is before sunrise.
938 !-------------------------------------------------------------------
939  where (t < -h .and. tt < -h) cosz = 0.0
940 
941 !-------------------------------------------------------------------
942 ! case 2: averaging period begins before sunrise, ends after sunrise
943 ! but before sunset
944 !-------------------------------------------------------------------
945  where ( (tt+h) /= 0.0 .and. t < -h .and. abs(tt) <= h) &
946  cosz = aa + bb*(stt + sh)/ (tt + h)
947 
948 !-------------------------------------------------------------------
949 ! case 3: averaging period begins before sunrise, ends after sunset,
950 ! but before the next sunrise. modify if averaging period extends
951 ! past the next day's sunrise, but if averaging period is less than
952 ! a half- day (pi) that circumstance will never occur.
953 !-------------------------------------------------------------------
954  where (t < -h .and. h /= 0.0 .and. h < tt) &
955  cosz = aa + bb*( sh + sh)/(h+h)
956 
957 !-------------------------------------------------------------------
958 ! case 4: averaging period begins after sunrise, ends before sunset.
959 !-------------------------------------------------------------------
960  where ( abs(t) <= h .and. abs(tt) <= h) &
961  cosz = aa + bb*(stt - st)/ (tt - t)
962 
963 !-------------------------------------------------------------------
964 ! case 5: averaging period begins after sunrise, ends after sunset.
965 ! modify when averaging period extends past the next day's sunrise.
966 !-------------------------------------------------------------------
967  where ((h-t) /= 0.0 .and. abs(t) <= h .and. h < tt) &
968  cosz = aa + bb*(sh - st)/(h-t)
969 
970 !-------------------------------------------------------------------
971 ! case 6: averaging period begins after sunrise , ends after the
972 ! next day's sunrise. note that this includes the case when the
973 ! day length is one day (h = pi).
974 !-------------------------------------------------------------------
975  where (twopi - h < tt .and. (tt+h-twopi) /= 0.0 .and. t <= h ) &
976  cosz = (cosz*(h - t) + (aa*(tt + h - twopi) + &
977  bb*(stt + sh))) / ((h - t) + (tt + h - twopi))
978 
979 !-------------------------------------------------------------------
980 ! case 7: averaging period begins after sunset and ends before the
981 ! next day's sunrise
982 !-------------------------------------------------------------------
983  where( h < t .and. twopi - h >= tt ) cosz = 0.0
984 
985 !-------------------------------------------------------------------
986 ! case 8: averaging period begins after sunset and ends after the
987 ! next day's sunrise but before the next day's sunset. if the
988 ! averaging period is less than a half-day (pi) the latter
989 ! circumstance will never occur.
990 !-----------------------------------------------------------------
991  where( h < t .and. twopi - h < tt )
992  cosz = aa + bb*(stt + sh) / (tt + h - twopi)
993  end where
994 
995  else
996  cosz = aa + bb*(stt - st)/ (tt - t)
997  end if
998 
999 
1000 
1001 !-------------------------------------------------------------------
1002 ! day fraction is the fraction of the averaging period contained
1003 ! within the (-h,h) period.
1004 !-------------------------------------------------------------------
1005  where (t < -h .and. tt < -h) fracday = 0.0
1006  where (t < -h .and. abs(tt) <= h) fracday = (tt + h )/dt
1007  where (t < -h .and. h < tt) fracday = ( h + h )/dt
1008  where (abs(t) <= h .and. abs(tt) <= h) fracday = (tt - t )/dt
1009  where (abs(t) <= h .and. h < tt) fracday = ( h - t )/dt
1010  where ( h < t ) fracday = 0.0
1011  where (twopi - h < tt) fracday = fracday + &
1012  (tt + h - &
1013  twopi)/dt
1014 !----------------------------------------------------------------------
1015 !> if instantaneous values are desired, define cosz at time t.
1016 !----------------------------------------------------------------------
1017  else ! (no time averaging)
1018  if (.not. lallow_negative) then
1019  where (abs(t) < h)
1020  cosz = aa + bb*cos(t)
1021  fracday = 1.0
1022  elsewhere
1023  cosz = 0.0
1024  fracday = 0.0
1025  end where
1026  else
1027  cosz = aa + bb*cos(t)
1028  where (abs(t) < h)
1029  fracday = 1.0
1030  elsewhere
1031  fracday = 0.0
1032  end where
1033  end if
1034  end if
1035 
1036 !----------------------------------------------------------------------
1037 !> Check that cosz is not negative, if desired.
1038 !----------------------------------------------------------------------
1039  if (.not. lallow_negative) then
1040  cosz = max(0.0, cosz)
1041  end if
1042 
1043 end subroutine diurnal_solar_2d
1044 
1045 
1046 !> \brief diurnal_solar_1d takes 1-d input fields, adds a second dimension
1047 !! and calls diurnal_solar_2d. on return, the 2d fields are returned
1048 !! to the original 1d fields.
1049 !!
1050 !! \param [in] <lat> Latitudes of model grid points
1051 !! \param [in] <lon> Longitudes of model grid points
1052 !! \param [in] <gmt> Time of day at longitude 0.0; midnight = 0.0, one day = 2 * pi
1053 !! \param [in] <time_since_ae> Time of year; autumnal equinox = 0.0, one year = 2 * pi
1054 !! \param [out] <cosz> Cosine of solar zenith angle
1055 !! \param [out] <fracday> Daylight fraction of time interval
1056 !! \param [out] <rrsun> earth-sun distance (r) relative to semi-major axis of orbital ellipse (a):(a/r)**2
1057 !! \param [in] <dt> OPTIONAL: time interval after gmt over which the astronomical variables are to be
1058 !! averaged. this produces averaged output rather than instantaneous.
1059 !! \param [in] <allow_negative_cosz> Allow negative values for cosz?
1060 !! \param [out] <half_day_out> half_day_out
1061 subroutine diurnal_solar_1d (lat, lon, gmt, time_since_ae, cosz, &
1062  fracday, rrsun, dt, allow_negative_cosz, &
1063  half_day_out)
1065 !---------------------------------------------------------------------
1066 real, dimension(:), intent(in) :: lat, lon
1067 real, intent(in) :: gmt, time_since_ae
1068 real, dimension(:), intent(out) :: cosz, fracday
1069 real, intent(out) :: rrsun
1070 real, intent(in), optional :: dt
1071 logical, intent(in), optional :: allow_negative_cosz
1072 real, dimension(:), intent(out), optional :: half_day_out
1073 
1074 !---------------------------------------------------------------------
1075 ! local variables
1076 !---------------------------------------------------------------------
1077  real, dimension(size(lat),1) :: lat_2d, lon_2d, cosz_2d, &
1078  fracday_2d,halfday_2d
1079 
1080 !--------------------------------------------------------------------
1081 !> define 2-d versions of input data arrays.
1082 !--------------------------------------------------------------------
1083  lat_2d(:,1) = lat
1084  lon_2d(:,1) = lon
1085 
1086 !--------------------------------------------------------------------
1087 !> call diurnal_solar_2d to calculate astronomy fields.
1088 !--------------------------------------------------------------------
1089 ! if (present(dt)) then
1090  call diurnal_solar_2d (lat_2d, lon_2d, gmt, time_since_ae,&
1091  cosz_2d, fracday_2d, rrsun, dt=dt, &
1092  allow_negative_cosz=allow_negative_cosz, &
1093  half_day_out=halfday_2d)
1094 ! else
1095 ! call diurnal_solar_2d (lat_2d, lon_2d, gmt, time_since_ae, &
1096 ! cosz_2d, fracday_2d, rrsun)
1097 ! endif
1098 
1099 !-------------------------------------------------------------------
1100 !> place output fields into 1-d arguments for return to
1101 !! calling routine.
1102 !-------------------------------------------------------------------
1103  fracday = fracday_2d(:,1)
1104  cosz = cosz_2d(:,1)
1105  if (present(half_day_out)) then
1106  half_day_out = halfday_2d(:,1)
1107  end if
1108 
1109 end subroutine diurnal_solar_1d
1110 
1111 
1112 !> \brief diurnal_solar_0d takes scalar input fields, makes them into 2d
1113 !! arrays dimensioned (1,1), and calls diurnal_solar_2d. on return,
1114 !! the 2d fields are converted back to the desired scalar output.
1115 !!
1116 !! \param [in] <lat> Latitudes of model grid points
1117 !! \param [in] <lon> Longitudes of model grid points
1118 !! \param [in] <gmt> Time of day at longitude 0.0; midnight = 0.0, one day = 2 * pi
1119 !! \param [in] <time_since_ae> Time of year; autumnal equinox = 0.0, one year = 2 * pi
1120 !! \param [out] <cosz> Cosine of solar zenith angle
1121 !! \param [out] <fracday> Daylight fraction of time interval
1122 !! \param [out] <rrsun> earth-sun distance (r) relative to semi-major axis of orbital ellipse (a):(a/r)**2
1123 !! \param [in] <dt> OPTIONAL: time interval after gmt over which the astronomical variables are to be
1124 !! averaged. this produces averaged output rather than instantaneous.
1125 !! \param [in] <allow_negative_cosz> Allow negative values for cosz?
1126 !! \param [out] <half_day_out> half_day_out
1127 subroutine diurnal_solar_0d (lat, lon, gmt, time_since_ae, cosz, &
1128  fracday, rrsun, dt, allow_negative_cosz, &
1129  half_day_out)
1131 real, intent(in) :: lat, lon, gmt, time_since_ae
1132 real, intent(out) :: cosz, fracday, rrsun
1133 real, intent(in), optional :: dt
1134 logical,intent(in),optional :: allow_negative_cosz
1135 real, intent(out), optional :: half_day_out
1136 
1137 !--------------------------------------------------------------------
1138 ! local variables:
1139 !--------------------------------------------------------------------
1140  real, dimension(1,1) :: lat_2d, lon_2d, cosz_2d, fracday_2d, halfday_2d
1141 
1142 !---------------------------------------------------------------------
1143 !> create 2d arrays from the scalar input fields.
1144 !---------------------------------------------------------------------
1145  lat_2d = lat
1146  lon_2d = lon
1147 
1148 !--------------------------------------------------------------------
1149 !> call diurnal_solar_2d to calculate astronomy fields.
1150 !--------------------------------------------------------------------
1151 ! if (present(dt)) then
1152  call diurnal_solar_2d (lat_2d, lon_2d, gmt, time_since_ae, &
1153  cosz_2d, fracday_2d, rrsun, dt=dt, &
1154  allow_negative_cosz=allow_negative_cosz, &
1155  half_day_out=halfday_2d)
1156 ! else
1157 ! call diurnal_solar_2d (lat_2d, lon_2d, gmt, time_since_ae, &
1158 ! cosz_2d, fracday_2d, rrsun)
1159 ! end if
1160 
1161 !-------------------------------------------------------------------
1162 !> place output fields into scalars for return to calling routine.
1163 !-------------------------------------------------------------------
1164  fracday = fracday_2d(1,1)
1165  cosz = cosz_2d(1,1)
1166  if (present(half_day_out)) then
1167  half_day_out = halfday_2d(1,1)
1168  end if
1169 
1170 end subroutine diurnal_solar_0d
1171 
1172 
1173 !> \brief diurnal_solar_cal_2d receives time_type inputs, converts
1174 !! them to real variables and then calls diurnal_solar_2d to
1175 !! compute desired astronomical variables.
1176 !!
1177 !! \param [in] <lat> Latitudes of model grid points
1178 !! \param [in] <lon> Longitudes of model grid points
1179 !! \param [in] <gmt> Time of day at longitude 0.0; midnight = 0.0, one day = 2 * pi
1180 !! \param [in] <time> Time of year (time_type)
1181 !! \param [out] <cosz> Cosine of solar zenith angle
1182 !! \param [out] <fracday> Daylight fraction of time interval
1183 !! \param [out] <rrsun> earth-sun distance (r) relative to semi-major axis of orbital ellipse (a):(a/r)**2
1184 !! \param [in] <dt> OPTIONAL: time interval after gmt over which the astronomical variables are to be
1185 !! averaged. this produces averaged output rather than instantaneous.
1186 !! \param [in] <allow_negative_cosz> Allow negative values for cosz?
1187 !! \param [out] <half_day_out> half_day_out
1188 !!
1189 !! \throw FATAL, "astronomy_mod radiation time step must be no longer than 12 hrs"
1190 !! \throw FATAL, "astronomy_mod radiation time step must not be an integral number of days"
1191 subroutine diurnal_solar_cal_2d (lat, lon, time, cosz, fracday, &
1192  rrsun, dt_time, allow_negative_cosz, &
1193  half_day_out)
1195 !-------------------------------------------------------------------
1196 real, dimension(:,:), intent(in) :: lat, lon
1197 type(time_type), intent(in) :: time
1198 real, dimension(:,:), intent(out) :: cosz, fracday
1199 real, intent(out) :: rrsun
1200 type(time_type), intent(in), optional :: dt_time
1201 logical, intent(in), optional :: allow_negative_cosz
1202 real, dimension(:,:), intent(out), optional :: half_day_out
1203 !---------------------------------------------------------------------
1204 
1205 !---------------------------------------------------------------------
1206 ! local variables
1207 !---------------------------------------------------------------------
1208  real :: dt
1209  real :: gmt, time_since_ae
1210 
1211 !---------------------------------------------------------------------
1212 !> Extract time of day (gmt) from time_type variable time with
1213 !! function universal_time.
1214 !---------------------------------------------------------------------
1215  gmt = universal_time(time)
1216 
1217 !---------------------------------------------------------------------
1218 !> Extract the time of year (time_since_ae) from time_type variable
1219 !! time using the function orbital_time.
1220 !---------------------------------------------------------------------
1221  time_since_ae = orbital_time(time)
1222 
1223 !---------------------------------------------------------------------
1224 !> Convert optional time_type variable dt_time (length of averaging
1225 !! period) to a real variable dt with the function universal_time.
1226 !---------------------------------------------------------------------
1227  if (present(dt_time)) then
1228  dt = universal_time(dt_time)
1229  if (dt > pi) then
1230  call error_mesg ( 'astronomy_mod', &
1231  'radiation time step must be no longer than 12 hrs', &
1232  fatal)
1233  endif
1234  if (dt == 0.0) then
1235  call error_mesg ( 'astronomy_mod', &
1236  'radiation time step must not be an integral &
1237  &number of days', fatal)
1238  endif
1239 
1240 !--------------------------------------------------------------------
1241 !> Call diurnal_solar_2d to calculate astronomy fields, with or
1242 !! without the optional argument dt.
1243 !--------------------------------------------------------------------
1244  call diurnal_solar_2d (lat, lon, gmt, time_since_ae, cosz, &
1245  fracday, rrsun, dt=dt, &
1246  allow_negative_cosz=allow_negative_cosz, &
1247  half_day_out=half_day_out)
1248  else
1249  call diurnal_solar_2d (lat, lon, gmt, time_since_ae, cosz, &
1250  fracday, rrsun, &
1251  allow_negative_cosz=allow_negative_cosz, &
1252  half_day_out=half_day_out)
1253  end if
1254 
1255 end subroutine diurnal_solar_cal_2d
1256 
1257 
1258 !> \brief diurnal_solar_cal_1d receives time_type inputs, converts
1259 !! them to real variables and then calls diurnal_solar_2d to
1260 !! compute desired astronomical variables.
1261 !!
1262 !! \param [in] <lat> Latitudes of model grid points
1263 !! \param [in] <lon> Longitudes of model grid points
1264 !! \param [in] <gmt> Time of day at longitude 0.0; midnight = 0.0, one day = 2 * pi
1265 !! \param [in] <time> Time of year (time_type)
1266 !! \param [out] <cosz> Cosine of solar zenith angle
1267 !! \param [out] <fracday> Daylight fraction of time interval
1268 !! \param [out] <rrsun> earth-sun distance (r) relative to semi-major axis of orbital ellipse (a):(a/r)**2
1269 !! \param [in] <dt> OPTIONAL: time interval after gmt over which the astronomical variables are to be
1270 !! averaged. this produces averaged output rather than instantaneous.
1271 !! \param [in] <allow_negative_cosz> Allow negative values for cosz?
1272 !! \param [out] <half_day_out> half_day_out
1273 subroutine diurnal_solar_cal_1d (lat, lon, time, cosz, fracday, &
1274  rrsun, dt_time, allow_negative_cosz, &
1275  half_day_out)
1277 !--------------------------------------------------------------------
1278 real, dimension(:), intent(in) :: lat, lon
1279 type(time_type), intent(in) :: time
1280 real, dimension(:), intent(out) :: cosz, fracday
1281 real, intent(out) :: rrsun
1282 type(time_type), intent(in), optional :: dt_time
1283 logical, intent(in), optional :: allow_negative_cosz
1284 real, dimension(:), intent(out), optional :: half_day_out
1285 !--------------------------------------------------------------------
1286 
1287 !-------------------------------------------------------------------
1288 ! local variables
1289 !-------------------------------------------------------------------
1290  real, dimension(size(lat),1) :: lat_2d, lon_2d, cosz_2d, &
1291  fracday_2d, halfday_2d
1292 
1293 !--------------------------------------------------------------------
1294 !> Define 2-d versions of input data arrays.
1295 !--------------------------------------------------------------------
1296  lat_2d(:,1) = lat
1297  lon_2d(:,1) = lon
1298 
1299 !--------------------------------------------------------------------
1300 !> Call diurnal_solar_cal_2d to convert the time_types to reals and
1301 !! then calculate the astronomy fields.
1302 !--------------------------------------------------------------------
1303  if (present(dt_time)) then
1304  call diurnal_solar_cal_2d (lat_2d, lon_2d, time, cosz_2d, &
1305  fracday_2d, rrsun, dt_time=dt_time, &
1306  allow_negative_cosz=allow_negative_cosz, &
1307  half_day_out=halfday_2d)
1308  else
1309  call diurnal_solar_cal_2d (lat_2d, lon_2d, time, cosz_2d, &
1310  fracday_2d, rrsun, &
1311  allow_negative_cosz=allow_negative_cosz, &
1312  half_day_out=halfday_2d)
1313  end if
1314 
1315 !-------------------------------------------------------------------
1316 !> Place output fields into 1-d arguments for return to
1317 !! calling routine.
1318 !-------------------------------------------------------------------
1319  fracday = fracday_2d(:,1)
1320  cosz = cosz_2d(:,1)
1321  if (present(half_day_out)) then
1322  half_day_out = halfday_2d(:,1)
1323  end if
1324 
1325 end subroutine diurnal_solar_cal_1d
1326 
1327 
1328 !> \brief diurnal_solar_cal_0d receives time_type inputs, converts them to real variables
1329 !! and then calls diurnal_solar_2d to compute desired astronomical variables.
1330 !!
1331 !! \param [in] <lat> Latitudes of model grid points
1332 !! \param [in] <lon> Longitudes of model grid points
1333 !! \param [in] <time> Time of year (time_type)
1334 !! \param [out] <cosz> Cosine of solar zenith angle
1335 !! \param [out] <fracday> Daylight fraction of time interval
1336 !! \param [out] <rrsun> Earth-Sun distance (r) relative to semi-major axis of orbital ellipse (a) : (a/r)**2
1337 !! \param [out] <dt_time> OPTIONAL: time interval after gmt over which the astronomical variables are
1338 !! to be averaged. this produces averaged output rather than instantaneous.
1339 !! \param [in] <allow_negative_cosz> allow_negative_cosz
1340 !! \param [out] <half_day_out> half_day_out
1341 subroutine diurnal_solar_cal_0d (lat, lon, time, cosz, fracday, &
1342  rrsun, dt_time, allow_negative_cosz, &
1343  half_day_out)
1345 !---------------------------------------------------------------------
1346 real, intent(in) :: lat, lon
1347 type(time_type), intent(in) :: time
1348 real, intent(out) :: cosz, fracday, rrsun
1349 type(time_type), intent(in), optional :: dt_time
1350 logical, intent(in), optional :: allow_negative_cosz
1351 real, intent(out), optional :: half_day_out
1352 !---------------------------------------------------------------------
1353 
1354 !---------------------------------------------------------------------
1355 ! local variables
1356 !---------------------------------------------------------------------
1357  real, dimension(1,1) :: lat_2d, lon_2d, cosz_2d, fracday_2d, halfday_2d
1358 
1359 !--------------------------------------------------------------------
1360 !> Define 2-d versions of input data arrays.
1361 !--------------------------------------------------------------------
1362  lat_2d = lat
1363  lon_2d = lon
1364 
1365 !--------------------------------------------------------------------
1366 !> Call diurnal_solar_cal_2d to convert the time_types to reals and
1367 !! then calculate the astronomy fields.
1368 !--------------------------------------------------------------------
1369  if (present(dt_time)) then
1370  call diurnal_solar_cal_2d (lat_2d, lon_2d, time, cosz_2d, &
1371  fracday_2d, rrsun, dt_time=dt_time, &
1372  allow_negative_cosz=allow_negative_cosz, &
1373  half_day_out=halfday_2d)
1374  else
1375  call diurnal_solar_cal_2d (lat_2d, lon_2d, time, cosz_2d, &
1376  fracday_2d, rrsun, &
1377  allow_negative_cosz=allow_negative_cosz, &
1378  half_day_out=halfday_2d)
1379  end if
1380 
1381 !-------------------------------------------------------------------
1382 !> Place output fields into 1-d arguments for return to
1383 !! calling routine.
1384 !-------------------------------------------------------------------
1385  fracday= fracday_2d(1,1)
1386  cosz = cosz_2d(1,1)
1387  if (present(half_day_out)) then
1388  half_day_out = halfday_2d(1,1)
1389  end if
1390 
1391 end subroutine diurnal_solar_cal_0d
1392 
1393 
1394 !> \brief daily_mean_solar_2d computes the daily mean astronomical parameters for
1395 !! the input points at latitude lat and time of year time_since_ae.
1396 !!
1397 !! \param [in] <lat> Latitudes of model grid points
1398 !! \param [in] <time_since_ae> Time of year; autumnal equinox = 0.0, one year = 2 * pi
1399 !! \param [out] <cosz> Cosine of solar zenith angle
1400 !! \param [out] <h_out> 2-d array of half-day lengths at the latitudes
1401 !! \param [out] <rr_out> the inverse of the square of the earth-sun distance relative
1402 !! to the mean distance at angle ang in the earth's orbit.
1403 !!
1404 !! \throw FATAL, "astronomy_mod time_since_ae not between 0 and 2pi"
1405 subroutine daily_mean_solar_2d (lat, time_since_ae, cosz, h_out, rr_out)
1407 !----------------------------------------------------------------------
1408 real, dimension(:,:), intent(in) :: lat
1409 real, intent(in) :: time_since_ae
1410 real, dimension(:,:), intent(out) :: cosz, h_out
1411 real, intent(out) :: rr_out
1412 !----------------------------------------------------------------------
1413 
1414 !--------------------------------------------------------------------
1415 ! local variables
1416 !--------------------------------------------------------------------
1417  real, dimension(size(lat,1),size(lat,2)) :: h
1418  real :: ang, dec, rr
1419 
1420 !--------------------------------------------------------------------
1421 ! be sure the time in the annual cycle is legitimate.
1422 !---------------------------------------------------------------------
1423  if (time_since_ae < 0.0 .or. time_since_ae > twopi) &
1424  call error_mesg('astronomy_mod', &
1425  'time_since_ae not between 0 and 2pi', fatal)
1426 
1427 !---------------------------------------------------------------------
1428 !> Define the orbital angle (location in year), solar declination,
1429 !! half-day length and earth sun distance factor. Use functions
1430 !! contained in this module.
1431 !---------------------------------------------------------------------
1432  ang = angle(time_since_ae)
1433  dec = declination(ang)
1434  h = half_day(lat, dec)
1435  rr = r_inv_squared(ang)
1436 
1437 !---------------------------------------------------------------------
1438 !> Where the entire day is dark, define cosz to be zero. otherwise
1439 !! use the standard formula. Define the daylight fraction and earth-
1440 !! sun distance.
1441 !---------------------------------------------------------------------
1442  where (h == 0.0)
1443  cosz = 0.0
1444  elsewhere
1445  cosz = sin(lat)*sin(dec) + cos(lat)*cos(dec)*sin(h)/h
1446  end where
1447  h_out = h/pi
1448  rr_out = rr
1449 
1450 end subroutine daily_mean_solar_2d
1451 
1452 
1453 !> \brief daily_mean_solar_1d takes 1-d input fields, adds a second dimension
1454 !! and calls daily_mean_solar_2d. on return, the 2d fields are
1455 !! returned to the original 1d fields.
1456 !!
1457 !! \param [in] <lat> Latitudes of model grid points
1458 !! \param [in] <time_since_ae> Time of year; autumnal equinox = 0.0, one year = 2 * pi
1459 !! \param [out] <cosz> Cosine of solar zenith angle
1460 !! \param [out] <h_out> 2-d array of half-day lengths at the latitudes
1461 !! \param [out] <rr_out> the inverse of the square of the earth-sun distance relative
1462 !! to the mean distance at angle ang in the earth's orbit.
1463 subroutine daily_mean_solar_1d (lat, time_since_ae, cosz, h_out, rr_out)
1465 !----------------------------------------------------------------------
1466 real, intent(in), dimension(:) :: lat
1467 real, intent(in) :: time_since_ae
1468 real, intent(out), dimension(size(lat(:))) :: cosz
1469 real, intent(out), dimension(size(lat(:))) :: h_out
1470 real, intent(out) :: rr_out
1471 !----------------------------------------------------------------------
1472 
1473 !----------------------------------------------------------------------
1474 ! local variables
1475 !----------------------------------------------------------------------
1476  real, dimension(size(lat),1) :: lat_2d, cosz_2d, hout_2d
1477 
1478 !--------------------------------------------------------------------
1479 !> define 2-d versions of input data array.
1480 !--------------------------------------------------------------------
1481  lat_2d(:,1) = lat
1482 
1483 !--------------------------------------------------------------------
1484 !> call daily_mean_solar_2d to calculate astronomy fields.
1485 !--------------------------------------------------------------------
1486  call daily_mean_solar_2d (lat_2d, time_since_ae, cosz_2d, &
1487  hout_2d, rr_out)
1488 
1489 !-------------------------------------------------------------------
1490 !> place output fields into 1-d arguments for return to
1491 !! calling routine.
1492 !-------------------------------------------------------------------
1493  h_out = hout_2d(:,1)
1494  cosz = cosz_2d(:,1)
1495 
1496 end subroutine daily_mean_solar_1d
1497 
1498 
1499 !> \brief daily_mean_solar_2level takes 1-d input fields, adds a second
1500 !! dimension and calls daily_mean_solar_2d. on return, the 2d fields
1501 !! are returned to the original 1d fields.
1502 !!
1503 !! \param [in] <lat> Latitudes of model grid points
1504 !! \param [in] <time_since_ae> Time of year; autumnal equinox = 0.0, one year = 2 * pi
1505 !! \param [out] <cosz> Cosine of solar zenith angle
1506 !! \param [out] <solar> Shortwave flux factor: cosine of zenith angle * daylight fraction / (earth-sun distance squared)
1507 subroutine daily_mean_solar_2level (lat, time_since_ae, cosz, solar)
1509 !----------------------------------------------------------------------
1510 real, intent(in), dimension(:) :: lat
1511 real, intent(in) :: time_since_ae
1512 real, intent(out), dimension(size(lat(:))) :: cosz, solar
1513 !----------------------------------------------------------------------
1514 
1515 !----------------------------------------------------------------------
1516 ! local variables
1517 !----------------------------------------------------------------------
1518  real, dimension(size(lat),1) :: lat_2d, cosz_2d, hout_2d
1519  real :: rr_out
1520 
1521 !--------------------------------------------------------------------
1522 !> define 2-d versions of input data array.
1523 !--------------------------------------------------------------------
1524  lat_2d(:,1) = lat
1525 
1526 !--------------------------------------------------------------------
1527 !> call daily_mean_solar_2d to calculate astronomy fields.
1528 !--------------------------------------------------------------------
1529  call daily_mean_solar_2d (lat_2d, time_since_ae, cosz_2d, &
1530  hout_2d, rr_out)
1531 
1532 !-------------------------------------------------------------------
1533 !> place output fields into 1-d arguments for return to
1534 !! calling routine.
1535 !-------------------------------------------------------------------
1536  solar = cosz_2d(:,1)*hout_2d(:,1)*rr_out
1537  cosz = cosz_2d(:,1)
1538 
1539 end subroutine daily_mean_solar_2level
1540 
1541 
1542 !> \brief daily_mean_solar_1d takes 1-d input fields, adds a second dimension
1543 !! and calls daily_mean_solar_2d. on return, the 2d fields are
1544 !! returned to the original 1d fields.
1545 !!
1546 !! \param [in] <lat> Latitudes of model grid points
1547 !! \param [in] <time_since_ae> Time of year; autumnal equinox = 0.0, one year = 2 * pi
1548 !! \param [out] <cosz> Cosine of solar zenith angle
1549 !! \param [out] <h_out> 2-d array of half-day lengths at the latitudes
1550 !! \param [out] <rr_out> the inverse of the square of the earth-sun distance relative to
1551 !! the mean distance at angle ang in the earth's orbit.
1552 subroutine daily_mean_solar_0d (lat, time_since_ae, cosz, h_out, rr_out)
1554 real, intent(in) :: lat, time_since_ae
1555 real, intent(out) :: cosz, h_out, rr_out
1556 
1557 !--------------------------------------------------------------------
1558 ! local variables
1559 !--------------------------------------------------------------------
1560  real, dimension(1,1) :: lat_2d, cosz_2d, hout_2d
1561 
1562 !--------------------------------------------------------------------
1563 !> define 2-d versions of input data array.
1564 !--------------------------------------------------------------------
1565  lat_2d = lat
1566 
1567 !--------------------------------------------------------------------
1568 !> call daily_mean_solar_2d to calculate astronomy fields.
1569 !--------------------------------------------------------------------
1570  call daily_mean_solar_2d (lat_2d, time_since_ae, cosz_2d, &
1571  hout_2d, rr_out)
1572 
1573 !-------------------------------------------------------------------
1574 !> return output fields to scalars for return to calling routine.
1575 !-------------------------------------------------------------------
1576  h_out = hout_2d(1,1)
1577  cosz = cosz_2d(1,1)
1578 
1579 end subroutine daily_mean_solar_0d
1580 
1581 
1582 !> \brief daily_mean_solar_cal_2d receives time_type inputs, converts
1583 !! them to real variables and then calls daily_mean_solar_2d to
1584 !! compute desired astronomical variables.
1585 !!
1586 !! \param [in] <lat> Latitudes of model grid points
1587 !! \param [in] <time> Time of year (time_type)
1588 !! \param [out] <cosz> Cosine of solar zenith angle
1589 !! \param [out] <fracday> Daylight fraction of time interval
1590 !! \param [out] <rrsun> Earth-Sun distance (r) relative to semi-major axis of orbital ellipse (a):(a/r)**2
1591 !!
1592 !! \throw FATAL, "astronomy_mod time_since_ae not between 0 and 2pi"
1593 subroutine daily_mean_solar_cal_2d (lat, time, cosz, fracday, rrsun)
1595 !-------------------------------------------------------------------
1596 real, dimension(:,:), intent(in) :: lat
1597 type(time_type), intent(in) :: time
1598 real, dimension(:,:), intent(out) :: cosz, fracday
1599 real, intent(out) :: rrsun
1600 !-------------------------------------------------------------------
1601 
1602 !-------------------------------------------------------------------
1603 ! local variables
1604 !-------------------------------------------------------------------
1605  real :: time_since_ae
1606 
1607 !--------------------------------------------------------------------
1608 ! be sure the time in the annual cycle is legitimate.
1609 !---------------------------------------------------------------------
1610  time_since_ae = orbital_time(time)
1611  if (time_since_ae < 0.0 .or. time_since_ae > twopi) &
1612  call error_mesg ('astronomy_mod', &
1613  'time_since_ae not between 0 and 2pi', fatal)
1614 
1615 !--------------------------------------------------------------------
1616 ! call daily_mean_solar_2d to calculate astronomy fields.
1617 !--------------------------------------------------------------------
1618  call daily_mean_solar_2d (lat, time_since_ae, cosz, &
1619  fracday, rrsun)
1620 
1621 end subroutine daily_mean_solar_cal_2d
1622 
1623 
1624 !> \brief daily_mean_solar_cal_1d receives time_type inputs, converts
1625 !! them to real, 2d variables and then calls daily_mean_solar_2d to
1626 !! compute desired astronomical variables.
1627 !!
1628 !! \param [in] <lat> Latitudes of model grid points
1629 !! \param [in] <time> Time of year (time_type)
1630 !! \param [out] <cosz> Cosine of solar zenith angle
1631 !! \param [out] <fracday> Daylight fraction of time interval
1632 !! \param [out] <rrsun> Earth-Sun distance (r) relative to semi-major axis of orbital ellipse (a):(a/r)**2
1633 subroutine daily_mean_solar_cal_1d (lat, time, cosz, fracday, rrsun)
1635 real, dimension(:), intent(in) :: lat
1636 type(time_type), intent(in) :: time
1637 real, dimension(:), intent(out) :: cosz, fracday
1638 real, intent(out) :: rrsun
1639 
1640 !---------------------------------------------------------------------
1641 ! local variables
1642 !---------------------------------------------------------------------
1643  real, dimension(size(lat),1) :: lat_2d, cosz_2d, fracday_2d
1644 
1645 
1646 !--------------------------------------------------------------------
1647 !> define 2-d versions of input data array.
1648 !--------------------------------------------------------------------
1649  lat_2d(:,1) = lat
1650 
1651 !--------------------------------------------------------------------
1652 !> call daily_mean_solar_cal_2d to convert the time_types to reals and
1653 !! then calculate the astronomy fields.
1654 !--------------------------------------------------------------------
1655  call daily_mean_solar_cal_2d (lat_2d, time, cosz_2d, &
1656  fracday_2d, rrsun)
1657 
1658 !-------------------------------------------------------------------
1659 !> place output fields into 1-d arguments for return to
1660 !! calling routine.
1661 !-------------------------------------------------------------------
1662  fracday = fracday_2d(:,1)
1663  cosz = cosz_2d(:,1)
1664 
1665 end subroutine daily_mean_solar_cal_1d
1666 
1667 
1668 !> \brief daily_mean_solar_cal_2level receives 1d arrays and time_type input,
1669 !! converts them to real, 2d variables and then calls
1670 !! daily_mean_solar_2d to compute desired astronomical variables.
1671 !!
1672 !! \param [in] <lat> Latitudes of model grid points
1673 !! \param [in] <time> Time of year (time_type)
1674 !! \param [out] <cosz> Cosine of solar zenith angle
1675 !! \param [out] <solar> Shortwave flux factor: cosine of zenith angle * daylight fraction / (earth-sun distance squared)
1676 subroutine daily_mean_solar_cal_2level (lat, time, cosz, solar)
1678 real, dimension(:), intent(in) :: lat
1679 type(time_type), intent(in) :: time
1680 real, dimension(:), intent(out) :: cosz, solar
1681 
1682 !---------------------------------------------------------------------
1683 ! local variables
1684 !---------------------------------------------------------------------
1685  real, dimension(size(lat),1) :: lat_2d, cosz_2d, fracday_2d
1686  real :: rrsun
1687 
1688 
1689 !--------------------------------------------------------------------
1690 !> define 2-d versions of input data array.
1691 !--------------------------------------------------------------------
1692  lat_2d(:,1) = lat
1693 
1694 !--------------------------------------------------------------------
1695 !> call daily_mean_solar_cal_2d to convert the time_types to reals and
1696 !! then calculate the astronomy fields.
1697 !--------------------------------------------------------------------
1698  call daily_mean_solar_cal_2d (lat_2d, time, cosz_2d, &
1699  fracday_2d, rrsun)
1700 
1701 !-------------------------------------------------------------------
1702 !> place output fields into 1-d arguments for return to
1703 !! calling routine.
1704 !-------------------------------------------------------------------
1705  solar = cosz_2d(:,1)*fracday_2d(:,1)*rrsun
1706  cosz = cosz_2d(:,1)
1707 
1708 end subroutine daily_mean_solar_cal_2level
1709 
1710 
1711 !> \brief daily_mean_solar_cal_0d converts scalar input fields to real, 2d variables and
1712 !! then calls daily_mean_solar_2d to compute desired astronomical variables.
1713 !!
1714 !! \param [in] <lat> Latitudes of model grid points
1715 !! \param [in] <time> Time of year (time_type)
1716 !! \param [out] <cosz> Cosine of solar zenith angle
1717 !! \param [out] <fracday> Daylight fraction of time interval
1718 !! \param [out] <rrsun> Earth-Sun distance (r) relative to semi-major axis of orbital ellipse (a):(a/r)**2
1719 subroutine daily_mean_solar_cal_0d (lat, time, cosz, fracday, rrsun)
1721 !--------------------------------------------------------------------
1722 real, intent(in) :: lat
1723 type(time_type), intent(in) :: time
1724 real, intent(out) :: cosz, fracday, rrsun
1725 !--------------------------------------------------------------------
1726 
1727 !--------------------------------------------------------------------
1728 ! local variables
1729 !--------------------------------------------------------------------
1730  real, dimension(1,1) :: lat_2d, cosz_2d, fracday_2d
1731 
1732 !--------------------------------------------------------------------
1733 !> define 2-d versions of input data array.
1734 !--------------------------------------------------------------------
1735  lat_2d = lat
1736 
1737 !--------------------------------------------------------------------
1738 !> call daily_mean_solar_cal_2d to convert the time_types to reals and
1739 !! then calculate the astronomy fields.
1740 !--------------------------------------------------------------------
1741  call daily_mean_solar_cal_2d (lat_2d, time, cosz_2d, &
1742  fracday_2d, rrsun)
1743 
1744 !-------------------------------------------------------------------
1745 !> place output fields into scalar arguments for return to
1746 !! calling routine.
1747 !-------------------------------------------------------------------
1748  fracday = fracday_2d(1,1)
1749  cosz = cosz_2d(1,1)
1750 
1751 end subroutine daily_mean_solar_cal_0d
1752 
1753 
1754 !> \brief annual_mean_solar_2d returns 2d fields of annual mean values of the cosine of
1755 !! zenith angle, daylight fraction and earth-sun distance at the specified latitude.
1756 !!
1757 !! \param [in] <jst> Starting index of latitude window
1758 !! \param [in] <jnd> Ending index of latitude window
1759 !! \param [in] <lat> Latitudes of model grid points
1760 !! \param [out] <cosz> Cosine of solar zenith angle
1761 !! \param [out] <solar> Shortwave flux factor: cosine of zenith angle * daylight fraction / (earth-sun distance squared)
1762 !! \param [out] <fracday> Daylight fraction of time interval
1763 !! \param [out] <rrsun> Earth-Sun distance (r) relative to semi-major axis of orbital ellipse (a):(a/r)**2
1764 subroutine annual_mean_solar_2d (js, je, lat, cosz, solar, fracday, &
1765  rrsun)
1767 !--------------------------------------------------------------------
1768 integer, intent(in) :: js, je
1769 real, dimension(:,:), intent(in) :: lat
1770 real, dimension(:,:), intent(out) :: solar, cosz, fracday
1771 real, intent(out) :: rrsun
1772 !--------------------------------------------------------------------
1773 
1774 !--------------------------------------------------------------------
1775 ! local variables
1776 !--------------------------------------------------------------------
1777  real, dimension(size(lat,1),size(lat,2)) :: s,z
1778  real :: t
1779  integer :: n, i
1780 
1781 !--------------------------------------------------------------------
1782 ! if the calculation has not yet been done, do it here.
1783 !--------------------------------------------------------------------
1784  if (.not. annual_mean_calculated) then
1785 
1786 !----------------------------------------------------------------------
1787 !> determine annual mean values of solar flux and product of cosz
1788 !! and solar flux by integrating the annual cycle in num_angles
1789 !! orbital increments.
1790 !----------------------------------------------------------------------
1791  solar = 0.0
1792  cosz = 0.0
1793  do n =1, num_angles
1794  t = float(n-1)*twopi/float(num_angles)
1795  call daily_mean_solar(lat,t, z, fracday, rrsun)
1796  s = z*rrsun*fracday
1797  solar = solar + s
1798  cosz = cosz + z*s
1799  end do
1800  solar = solar/float(num_angles)
1801  cosz = cosz/float(num_angles)
1802 
1803 !--------------------------------------------------------------------
1804 !> define the flux-weighted annual mean cosine of the zenith angle.
1805 !--------------------------------------------------------------------
1806  where(solar.eq.0.0)
1807  cosz = 0.0
1808  elsewhere
1809  cosz = cosz/solar
1810  end where
1811 
1812 !-------------------------------------------------------------------
1813 !> define avg fracday such as to make the avg flux (solar) equal to
1814 !! the product of the avg cosz * avg fracday * assumed mean avg
1815 !! radius of 1.0. it is unlikely that these avg fracday and avg rr
1816 !! will ever be used.
1817 !--------------------------------------------------------------------
1818  where(solar .eq.0.0)
1819  fracday = 0.0
1820  elsewhere
1821  fracday = solar/cosz
1822  end where
1823  rrsun = 1.00
1824 
1825 !---------------------------------------------------------------------
1826 !> save the values that have been calculated as module variables, if
1827 !! those variables are present; i.e., not the spectral 2-layer model.
1828 !---------------------------------------------------------------------
1829  if (allocated (cosz_ann)) then
1830  cosz_ann = cosz
1831  solar_ann = solar
1832  fracday_ann = fracday
1833  rrsun_ann = rrsun
1834 
1835 !--------------------------------------------------------------------
1836 !> increment the points computed counter. set flag to end execution
1837 !! once values have been calculated for all points owned by the
1838 !! processor.
1839 !--------------------------------------------------------------------
1840  num_pts = num_pts + size(lat,1)*size(lat,2)
1841  if ( num_pts == total_pts) annual_mean_calculated = .true.
1842  endif
1843 
1844 !--------------------------------------------------------------------
1845 !> if the calculation has been done, return the appropriate module
1846 !! variables.
1847 !--------------------------------------------------------------------
1848  else
1849  if (allocated (cosz_ann)) then
1850  cosz = cosz_ann
1851  solar = solar_ann
1852  fracday = fracday_ann
1853  rrsun = rrsun_ann
1854  endif
1855  endif
1856 
1857 end subroutine annual_mean_solar_2d
1858 
1859 
1860 !> \brief annual_mean_solar_1d creates 2-d input fields from 1-d input fields and then calls
1861 !! annual_mean_solar_2d to obtain 2-d output fields which are then stored into 1-d
1862 !! fields for return to the calling subroutine.
1863 !!
1864 !! \param [in] <jst> Starting index of latitude window
1865 !! \param [in] <jnd> Ending index of latitude window
1866 !! \param [in] <lat> Latitudes of model grid points
1867 !! \param [out] <cosz> Cosine of solar zenith angle
1868 !! \param [out] <solar> Shortwave flux factor: cosine of zenith angle * daylight fraction / (earth-sun distance squared)
1869 !! \param [out] <fracday> Daylight fraction of time interval
1870 !! \param [out] <rrsun_out> Earth-Sun distance (r) relative to semi-major axis of orbital ellipse (a):(a/r)**2
1871 subroutine annual_mean_solar_1d (jst, jnd, lat, cosz, solar, &
1872  fracday, rrsun_out)
1874 !---------------------------------------------------------------------
1875 integer, intent(in) :: jst, jnd
1876 real, dimension(:), intent(in) :: lat(:)
1877 real, dimension(:), intent(out) :: cosz, solar, fracday
1878 real, intent(out) :: rrsun_out
1879 !---------------------------------------------------------------------
1880 
1881 !---------------------------------------------------------------------
1882 ! local variables
1883 
1884  real, dimension(size(lat),1) :: lat_2d, solar_2d, cosz_2d, &
1885  fracday_2d
1886  real :: rrsun
1887 
1888 !--------------------------------------------------------------------
1889 ! if the calculation has not been done, do it here.
1890 !--------------------------------------------------------------------
1891  if ( .not. annual_mean_calculated) then
1892 
1893 !--------------------------------------------------------------------
1894 !> define 2-d versions of input data array.
1895 !--------------------------------------------------------------------
1896  lat_2d(:,1) = lat
1897 
1898 !--------------------------------------------------------------------
1899 !> call annual_mean_solar_2d to calculate the astronomy fields.
1900 !--------------------------------------------------------------------
1901  call annual_mean_solar_2d (jst, jnd, lat_2d, cosz_2d, &
1902  solar_2d, fracday_2d, rrsun)
1903 
1904 !-------------------------------------------------------------------
1905 !> place output fields into 1-D arrays for return to calling routine.
1906 !-------------------------------------------------------------------
1907  fracday = fracday_2d(:,1)
1908  rrsun_out = rrsun
1909  solar = solar_2d(:,1)
1910  cosz = cosz_2d(:,1)
1911 
1912 !--------------------------------------------------------------------
1913 !> if the calculation has been done, simply return the module
1914 !! variables contain the results at the desired latitudes.
1915 !--------------------------------------------------------------------
1916  else
1917  cosz(:) = cosz_ann(1,jst:jnd)
1918  solar(:) = solar_ann(1,jst:jnd)
1919  fracday(:) = fracday_ann(1,jst:jnd)
1920  rrsun = rrsun_ann
1921  endif
1922 
1923 end subroutine annual_mean_solar_1d
1924 
1925 
1926 !> \brief annual_mean_solar_2level creates 2-d input fields from 1-d input fields
1927 !! and then calls annual_mean_solar_2d to obtain 2-d output fields which are
1928 !! then stored into 1-d fields for return to the calling subroutine. This
1929 !! subroutine will be called during model initialization.
1930 !!
1931 !! \throw FATAL, "astronomy_mod annual_mean_solar_2level should be called only once"
1932 subroutine annual_mean_solar_2level (lat, cosz, solar)
1934 !---------------------------------------------------------------------
1935 real, dimension(:), intent(in) :: lat !< Latitudes of model grid points
1936 real, dimension(:), intent(out) :: cosz !< Cosine of solar zenith angle
1937 real, dimension(:), intent(out) :: solar !< shortwave flux factor: cosine of zenith angle *
1938  !! daylight fraction / (earth-sun distance squared)
1939 !---------------------------------------------------------------------
1940 
1941 !---------------------------------------------------------------------
1942 ! local variables
1943 
1944  real, dimension(size(lat),1) :: lat_2d, solar_2d, cosz_2d, &
1945  fracday_2d
1946  integer :: jst, jnd
1947  real :: rrsun
1948 
1949 !--------------------------------------------------------------------
1950 ! if the calculation has not been done, do it here.
1951 !--------------------------------------------------------------------
1952  if ( .not. annual_mean_calculated) then
1953 
1954 !--------------------------------------------------------------------
1955 !> define 2-d versions of input data array.
1956 !--------------------------------------------------------------------
1957  lat_2d(:,1) = lat
1958  jst = 1
1959  jnd = size(lat(:))
1960 
1961 !--------------------------------------------------------------------
1962 !> call annual_mean_solar_2d to calculate the astronomy fields.
1963 !--------------------------------------------------------------------
1964  call annual_mean_solar_2d (jst, jnd, lat_2d, cosz_2d, &
1965  solar_2d, fracday_2d, rrsun)
1966 
1967 !-------------------------------------------------------------------
1968 !> place output fields into 1-D arrays for return to calling routine.
1969 !-------------------------------------------------------------------
1970  solar = solar_2d(:,1)
1971  cosz = cosz_2d(:,1)
1972 
1973 !--------------------------------------------------------------------
1974 !> if the calculation has been done, print an error message since
1975 !! this subroutine should be called only once.
1976 !--------------------------------------------------------------------
1977  else
1978  call error_mesg ('astronomy_mod', &
1979  'annual_mean_solar_2level should be called only once', &
1980  fatal)
1981  endif
1982  annual_mean_calculated = .true.
1983 
1984 end subroutine annual_mean_solar_2level
1985 
1986 
1987 !> \brief astronomy_end is the destructor for astronomy_mod.
1988 subroutine astronomy_end
1990 !----------------------------------------------------------------------
1991 !> check if the module has been initialized.
1992 !----------------------------------------------------------------------
1993  if (.not. module_is_initialized) return
1994 ! call error_mesg ( 'astronomy_mod', &
1995 ! ' module has not been initialized', FATAL)
1996 
1997 !----------------------------------------------------------------------
1998 !> deallocate module variables.
1999 !----------------------------------------------------------------------
2000  deallocate (orb_angle)
2001  if (allocated(cosz_ann) ) then
2002  deallocate (cosz_ann)
2003  deallocate (fracday_ann)
2004  deallocate (solar_ann)
2005  endif
2006 
2007 !----------------------------------------------------------------------
2008 !> mark the module as uninitialized.
2009 !----------------------------------------------------------------------
2010  module_is_initialized = .false.
2011 
2012 end subroutine astronomy_end
2013 
2014 
2015 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
2016 !
2017 ! PRIVATE SUBROUTINES
2018 !
2019 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
2020 
2021 !> \brief Orbit computes and stores a table of value of orbital angles as a
2022 !! function of orbital time (both the angle and time are zero at
2023 !! autumnal equinox in the NH, and range from 0 to 2*pi).
2024 subroutine orbit
2026 !---------------------------------------------------------------------
2027 ! local variables
2028 
2029  integer :: n
2030  real :: d1, d2, d3, d4, d5, dt, norm
2031 
2032 !--------------------------------------------------------------------
2033 !> allocate the orbital angle array, sized by the namelist parameter
2034 !! num_angles, defining the annual cycle resolution of the earth's
2035 !! orbit. define some constants to be used.
2036 !--------------------------------------------------------------------
2037 ! wfc moving to astronomy_init
2038 ! allocate ( orb_angle(0:num_angles) )
2039  orb_angle(0) = 0.0
2040  dt = twopi/float(num_angles)
2041  norm = sqrt(1.0 - ecc**2)
2042  dt = dt*norm
2043 
2044 !---------------------------------------------------------------------
2045 !> define the orbital angle at each of the num_angles locations in
2046 !! the orbit.
2047 !---------------------------------------------------------------------
2048  do n = 1, num_angles
2049  d1 = dt*r_inv_squared(orb_angle(n-1))
2050  d2 = dt*r_inv_squared(orb_angle(n-1)+0.5*d1)
2051  d3 = dt*r_inv_squared(orb_angle(n-1)+0.5*d2)
2052  d4 = dt*r_inv_squared(orb_angle(n-1)+d3)
2053  d5 = d1/6.0 + d2/3.0 +d3/3.0 +d4/6.0
2054  orb_angle(n) = orb_angle(n-1) + d5
2055  end do
2056 
2057 end subroutine orbit
2058 
2059 
2060 !> \brief r_inv_squared returns the inverse of the square of the earth-sun
2061 !! distance relative to the mean distance at angle ang in the Earth's orbit.
2062 function r_inv_squared (ang)
2064 !--------------------------------------------------------------------
2065 real, intent(in) :: ang !< angular position of earth in its orbit, relative to a
2066  !! value of 0.0 at the NH autumnal equinox, value between
2067  !! 0.0 and 2 * pi [radians]
2068 !--------------------------------------------------------------------
2069 
2070 !---------------------------------------------------------------------
2071 ! local variables
2072 
2073 real :: r_inv_squared !< The inverse of the square of the earth-sun distance relative
2074  !! to the mean distance [dimensionless]
2075 real :: r !< Earth-Sun distance relative to mean distance [dimensionless]
2076 real :: rad_per !< Angular position of perihelion [radians]
2077 
2078 !--------------------------------------------------------------------
2079 !> define the earth-sun distance (r) and then return the inverse of
2080 !! its square (r_inv_squared) to the calling routine.
2081 !--------------------------------------------------------------------
2082  rad_per = per*deg_to_rad
2083  r = (1. - ecc**2)/(1. + ecc*cos(ang - rad_per))
2084  r_inv_squared = r**(-2)
2085 
2086 
2087 end function r_inv_squared
2088 
2089 
2090 !> \brief angle determines the position within the earth's orbit at time t
2091 !! in the year (t = 0 at NH autumnal equinox) by interpolating
2092 !! into the orbital position table.
2093 function angle (t)
2095 !--------------------------------------------------------------------
2096 real, intent(in) :: t !< time of year (between 0 and 2*pi; t=0 at NH autumnal equinox
2097 !--------------------------------------------------------------------
2098 
2099 !--------------------------------------------------------------------
2100 ! local variables
2101 !--------------------------------------------------------------------
2102  real :: angle !< Orbital position relative to NH autumnal equinox [radians]
2103  real :: norm_time !< Index into orbital table corresponding to input time [dimensionless]
2104  real :: x !< Fractional distance between the orbital table entries bracketing the input time [dimensionless]
2105  integer :: int !< Table index which is lower than actual position, but closest to it [dimensionless]
2106  integer :: int_1 !< Next table index just larger than actual orbital position [dimensionless]
2107 
2108 !--------------------------------------------------------------------
2109 !> Define orbital tables indices bracketing current orbital time
2110 !! (int and int_1). Define table index distance between the lower
2111 !! table value (int) and the actual orbital time (x). Define orbital
2112 !! position as being x of the way between int and int_1. Renormalize
2113 !! angle to be within the range 0 to 2*pi.
2114 !--------------------------------------------------------------------
2115  norm_time = t*float(num_angles)/twopi
2116  int = floor(norm_time)
2117  int = modulo(int,num_angles)
2118  int_1 = int+1
2119  x = norm_time - floor(norm_time)
2120  angle = (1.0 -x)*orb_angle(int) + x*orb_angle(int_1)
2121  angle = modulo(angle, twopi)
2122 
2123 end function angle
2124 
2125 
2126 !> \brief Declination returns the solar declination angle at orbital
2127 !! position ang in earth's orbit.
2128 function declination (ang)
2130 !--------------------------------------------------------------------
2131 real, intent(in) :: ang !< solar orbital position ang in earth's orbit
2132 !--------------------------------------------------------------------
2133 
2134 !--------------------------------------------------------------------
2135 ! local variables
2136 
2137  real :: declination !< Solar declination angle [radians]
2138  real :: rad_obliq !< Obliquity of the ecliptic [radians]
2139  real :: sin_dec !< Sine of the solar declination [dimensionless]
2140 
2141 !---------------------------------------------------------------------
2142 ! compute the solar declination.
2143 !---------------------------------------------------------------------
2144  rad_obliq = obliq*deg_to_rad
2145  sin_dec = - sin(rad_obliq)*sin(ang)
2146  declination = asin(sin_dec)
2147 
2148 end function declination
2149 
2150 
2151 !> \brief half_day_2d returns a 2-d array of half-day lengths at the
2152 !! latitudes and declination provided.
2153 !!
2154 function half_day_2d (latitude, dec) result(h)
2156 !---------------------------------------------------------------------
2157 real, dimension(:,:), intent(in) :: latitude !< Latitutde of view point
2158 real, intent(in) :: dec !< Solar declination angle at view point
2159 real, dimension(size(latitude,1),size(latitude,2)) :: h
2160 !---------------------------------------------------------------------
2161 
2162 !---------------------------------------------------------------------
2163 ! local variables
2164 !---------------------------------------------------------------------
2165  real, dimension (size(latitude,1),size(latitude,2)):: &
2166  cos_half_day, & !< Cosine of half-day length [dimensionless]
2167  lat !< Model latitude, adjusted so that it is never 0.5*pi or -0.5*pi
2168  real :: tan_dec !< tangent of solar declination [dimensionless]
2169  real :: eps = 1.0e-05 !< small increment
2170 
2171 !--------------------------------------------------------------------
2172 !> define tangent of the declination.
2173 !--------------------------------------------------------------------
2174  tan_dec = tan(dec)
2175 
2176 !--------------------------------------------------------------------
2177 !> adjust latitude so that its tangent will be defined.
2178 !--------------------------------------------------------------------
2179  lat = latitude
2180  where (latitude == 0.5*pi) lat= latitude - eps
2181  where (latitude == -0.5*pi) lat= latitude + eps
2182 
2183 !--------------------------------------------------------------------
2184 !> define the cosine of the half-day length. adjust for cases of
2185 !! all daylight or all night.
2186 !--------------------------------------------------------------------
2187  cos_half_day = -tan(lat)*tan_dec
2188  where (cos_half_day <= -1.0) h = pi
2189  where (cos_half_day >= +1.0) h = 0.0
2190  where(cos_half_day > -1.0 .and. cos_half_day < 1.0) &
2191  h = acos(cos_half_day)
2192 
2193 end function half_day_2d
2194 
2195 
2196 !> \brief half_day_0d takes scalar input fields, makes them into 2-d fields
2197 !! dimensioned (1,1), and calls half_day_2d. On return, the 2-d
2198 !! fields are converted to the desired scalar output.
2199 !!
2200 !! \param [in] <latitude> Latitutde of view point
2201 !! \param [in] <dec> Solar declination angle at view point
2202 function half_day_0d(latitude, dec) result(h)
2204 real, intent(in) :: latitude, dec
2205 real :: h
2206 
2207 !----------------------------------------------------------------------
2208 ! local variables
2209 !----------------------------------------------------------------------
2210  real, dimension(1,1) :: lat_2d, h_2d
2211 
2212 !---------------------------------------------------------------------
2213 ! create 2d array from the input latitude field.
2214 !---------------------------------------------------------------------
2215  lat_2d = latitude
2216 
2217 !---------------------------------------------------------------------
2218 ! call half_day with the 2d arguments to calculate half-day length.
2219 !---------------------------------------------------------------------
2220  h_2d = half_day(lat_2d, dec)
2221 
2222 !---------------------------------------------------------------------
2223 ! create scalar from 2d array.
2224 !---------------------------------------------------------------------
2225  h = h_2d(1,1)
2226 
2227 end function half_day_0d
2228 
2229 
2230 !> \brief Orbital time returns the time (1 year = 2*pi) since autumnal
2231 !! equinox
2232 !!
2233 !! \details Orbital time returns the time (1 year = 2*pi) since autumnal
2234 !! equinox; autumnal_eq_ref is a module variable of time_type and
2235 !! will have been defined by default or by a call to
2236 !! set_ref_date_of_ae; length_of_year is available through the time
2237 !! manager and is set at the value approriate for the calandar being used
2238 function orbital_time(time) result(t)
2240 type(time_type), intent(in) :: time !< time (1 year = 2*pi) since autumnal equinox
2241 real :: t
2242 
2243  t = real( (time - autumnal_eq_ref)//period_time_type)
2244  t = twopi*(t - floor(t))
2245  if (time < autumnal_eq_ref) t = twopi - t
2246 
2247 
2248 end function orbital_time
2249 
2250 
2251 !> \brief universal_time returns the time of day at longitude = 0.0
2252 !! (1 day = 2*pi)
2253 function universal_time(time) result(t)
2255 type(time_type), intent(in) :: time !< Time (1 year = 2*pi) since autumnal equinox
2256 real :: t
2257 
2258 !--------------------------------------------------------------------
2259 ! local variables
2260 !--------------------------------------------------------------------
2261  integer :: seconds, days
2262 
2263  call get_time (time, seconds, days)
2264  t = twopi*real(seconds)/seconds_per_day
2265 
2266 end function universal_time
2267 
2268 
2269  end module astronomy_mod
Definition: fms.F90:20
real deg_to_rad
conversion from degrees to radians
Definition: astronomy.F90:386
type(time_type) function, public length_of_year()
integer year_ae
Year of specified autumnal equinox.
Definition: astronomy.F90:353
subroutine diurnal_solar_cal_2d(lat, lon, time, cosz, fracday, rrsun, dt_time, allow_negative_cosz, half_day_out)
diurnal_solar_cal_2d receives time_type inputs, converts them to real variables and then calls diurna...
Definition: astronomy.F90:1194
type(time_type) function, public set_date_julian(year, month, day, hour, minute, second)
subroutine set_period_integer(period_in)
set_period_integer saves as the input length of the year (an integer) in a time_type module variable...
Definition: astronomy.F90:591
subroutine daily_mean_solar_0d(lat, time_since_ae, cosz, h_out, rr_out)
daily_mean_solar_1d takes 1-d input fields, adds a second dimension and calls daily_mean_solar_2d. on return, the 2d fields are returned to the original 1d fields.
Definition: astronomy.F90:1553
subroutine, public astronomy_end
astronomy_end is the destructor for astronomy_mod.
Definition: astronomy.F90:1989
type(time_type) autumnal_eq_ref
time_type variable containing specified time of reference NH autumnal equinox
Definition: astronomy.F90:373
subroutine daily_mean_solar_cal_0d(lat, time, cosz, fracday, rrsun)
daily_mean_solar_cal_0d converts scalar input fields to real, 2d variables and then calls daily_mean_...
Definition: astronomy.F90:1720
subroutine, public get_ref_date_of_ae(day_out, month_out, year_out, second_out, minute_out, hour_out)
get_ref_date_of_ae retrieves the reference date of the autumnal equinox as integer variables...
Definition: astronomy.F90:794
subroutine daily_mean_solar_2level(lat, time_since_ae, cosz, solar)
daily_mean_solar_2level takes 1-d input fields, adds a second dimension and calls daily_mean_solar_2d...
Definition: astronomy.F90:1508
logical module_is_initialized
has the module been initialized ?
Definition: astronomy.F90:388
subroutine, public set_orbital_parameters(ecc_in, obliq_in, per_in)
set_orbital_parameters saves the input values of eccentricity, obliquity and perihelion time as modul...
Definition: astronomy.F90:643
integer minute_ae
Minute of specified autumnal equinox.
Definition: astronomy.F90:355
real, dimension(:,:), allocatable fracday_ann
annual mean daylight fraction
Definition: astronomy.F90:390
integer day_ae
Day of specified autumnal equinox.
Definition: astronomy.F90:351
real function, dimension(size(latitude, 1), size(latitude, 2)) half_day_2d(latitude, dec)
half_day_2d returns a 2-d array of half-day lengths at the latitudes and declination provided...
Definition: astronomy.F90:2155
integer hour_ae
Hour of specified autumnal equinox.
Definition: astronomy.F90:354
Definition: mpp.F90:39
integer function, public check_nml_error(IOSTAT, NML_NAME)
Definition: fms.F90:658
real function, public orbital_time(time)
Orbital time returns the time (1 year = 2*pi) since autumnal equinox.
Definition: astronomy.F90:2239
subroutine get_period_integer(period_out)
get_period_integer returns the length of the year as an integer number of seconds.
Definition: astronomy.F90:540
real function, public universal_time(time)
universal_time returns the time of day at longitude = 0.0 (1 day = 2*pi)
Definition: astronomy.F90:2254
subroutine daily_mean_solar_cal_2d(lat, time, cosz, fracday, rrsun)
daily_mean_solar_cal_2d receives time_type inputs, converts them to real variables and then calls dai...
Definition: astronomy.F90:1594
integer month_ae
Month of specified autumnal equinox.
Definition: astronomy.F90:352
subroutine, public set_ref_date_of_ae(day_in, month_in, year_in, second_in, minute_in, hour_in)
set_ref_date_of_ae provides a means of specifying the reference date of the NH autumnal equinox for a...
Definition: astronomy.F90:743
subroutine set_period_time_type(period_in)
Set_period_time_type saves the length of the year (input as a time_type variable) into a time_type mo...
Definition: astronomy.F90:615
integer period
Specified length of year [seconds]; must be specified to override default value given by length_of_ye...
Definition: astronomy.F90:347
character(len=input_str_length), dimension(:), allocatable, target, public input_nml_file
Definition: mpp.F90:1378
real obliq
Obliquity [degrees].
Definition: astronomy.F90:344
real ecc
Eccentricity of Earth&#39;s orbit [dimensionless].
Definition: astronomy.F90:343
real twopi
2 *PI
Definition: astronomy.F90:387
real, parameter, public pi
Ratio of circle circumference to diameter [N/A].
Definition: constants.F90:74
integer num_angles
Number of intervals into which the year is divided to compute orbital positions.
Definition: astronomy.F90:357
subroutine, public fms_init(localcomm)
Definition: fms.F90:353
subroutine diurnal_solar_1d(lat, lon, gmt, time_since_ae, cosz, fracday, rrsun, dt, allow_negative_cosz, half_day_out)
diurnal_solar_1d takes 1-d input fields, adds a second dimension and calls diurnal_solar_2d. on return, the 2d fields are returned to the original 1d fields.
Definition: astronomy.F90:1064
real function, private r_inv_squared(ang)
r_inv_squared returns the inverse of the square of the earth-sun distance relative to the mean distan...
Definition: astronomy.F90:2063
subroutine, public get_orbital_parameters(ecc_out, obliq_out, per_out)
get_orbital_parameters retrieves the orbital parameters for use by another module.
Definition: astronomy.F90:692
real function, private declination(ang)
Declination returns the solar declination angle at orbital position ang in earth&#39;s orbit...
Definition: astronomy.F90:2129
real function half_day_0d(latitude, dec)
half_day_0d takes scalar input fields, makes them into 2-d fields dimensioned (1,1), and calls half_day_2d. On return, the 2-d fields are converted to the desired scalar output.
Definition: astronomy.F90:2203
real, dimension(:), allocatable orb_angle
table of orbital positions (0 to 2*pi) as a function of time used to find actual orbital position via...
Definition: astronomy.F90:380
subroutine, public time_manager_init()
astronomy_mod provides astronomical variables for use by other modules within fms. The only currently used interface is for determination of astronomical values needed by the shortwave radiation packages.
Definition: astronomy.F90:56
integer num_pts
count of grid_boxes for which annual mean astronomy values have been calculated
Definition: astronomy.F90:396
type(time_type) period_time_type
time_type variable containing period of one orbit
Definition: astronomy.F90:377
subroutine, public astronomy_init(latb, lonb)
astronomy_init is the constructor for astronomy_mod.
Definition: astronomy.F90:421
real, dimension(:,:), allocatable cosz_ann
annual mean cos of zenith angle
Definition: astronomy.F90:390
logical annual_mean_calculated
have the annual mean values been calculated?
Definition: astronomy.F90:395
real rrsun_ann
annual mean earth-sun distance
Definition: astronomy.F90:394
real seconds_per_day
seconds in a day
Definition: astronomy.F90:385
subroutine daily_mean_solar_cal_2level(lat, time, cosz, solar)
daily_mean_solar_cal_2level receives 1d arrays and time_type input, converts them to real...
Definition: astronomy.F90:1677
subroutine diurnal_solar_cal_0d(lat, lon, time, cosz, fracday, rrsun, dt_time, allow_negative_cosz, half_day_out)
diurnal_solar_cal_0d receives time_type inputs, converts them to real variables and then calls diurna...
Definition: astronomy.F90:1344
subroutine annual_mean_solar_2d(js, je, lat, cosz, solar, fracday, rrsun)
annual_mean_solar_2d returns 2d fields of annual mean values of the cosine of zenith angle...
Definition: astronomy.F90:1766
integer total_pts
number of grid boxes owned by the processor
Definition: astronomy.F90:399
subroutine diurnal_solar_0d(lat, lon, gmt, time_since_ae, cosz, fracday, rrsun, dt, allow_negative_cosz, half_day_out)
diurnal_solar_0d takes scalar input fields, makes them into 2d arrays dimensioned (1...
Definition: astronomy.F90:1130
subroutine daily_mean_solar_2d(lat, time_since_ae, cosz, h_out, rr_out)
daily_mean_solar_2d computes the daily mean astronomical parameters for the input points at latitude ...
Definition: astronomy.F90:1406
integer second_ae
Second of specified autumnal equinox.
Definition: astronomy.F90:356
#define max(a, b)
Definition: mosaic_util.h:33
subroutine get_period_time_type(period_out)
get_period_time_type returns the length of the year as a time_type variable.
Definition: astronomy.F90:569
subroutine daily_mean_solar_cal_1d(lat, time, cosz, fracday, rrsun)
daily_mean_solar_cal_1d receives time_type inputs, converts them to real, 2d variables and then calls...
Definition: astronomy.F90:1634
subroutine diurnal_solar_2d(lat, lon, gmt, time_since_ae, cosz, fracday, rrsun, dt, allow_negative_cosz, half_day_out)
diurnal_solar_2d returns 2d fields of cosine of zenith angle, daylight fraction and earth-sun distanc...
Definition: astronomy.F90:841
subroutine annual_mean_solar_2level(lat, cosz, solar)
annual_mean_solar_2level creates 2-d input fields from 1-d input fields and then calls annual_mean_so...
Definition: astronomy.F90:1933
subroutine diurnal_solar_cal_1d(lat, lon, time, cosz, fracday, rrsun, dt_time, allow_negative_cosz, half_day_out)
diurnal_solar_cal_1d receives time_type inputs, converts them to real variables and then calls diurna...
Definition: astronomy.F90:1276
subroutine daily_mean_solar_1d(lat, time_since_ae, cosz, h_out, rr_out)
daily_mean_solar_1d takes 1-d input fields, adds a second dimension and calls daily_mean_solar_2d. on return, the 2d fields are returned to the original 1d fields.
Definition: astronomy.F90:1464
real function, private angle(t)
angle determines the position within the earth&#39;s orbit at time t in the year (t = 0 at NH autumnal eq...
Definition: astronomy.F90:2094
subroutine, public get_time(Time, seconds, days, ticks, err_msg)
subroutine, public error_mesg(routine, message, level)
Definition: fms.F90:529
real per
Longitude of perihelion with respect to autumnal equinox in NH [degrees].
Definition: astronomy.F90:345
subroutine annual_mean_solar_1d(jst, jnd, lat, cosz, solar, fracday, rrsun_out)
annual_mean_solar_1d creates 2-d input fields from 1-d input fields and then calls annual_mean_solar_...
Definition: astronomy.F90:1873
subroutine, private orbit
Orbit computes and stores a table of value of orbital angles as a function of orbital time (both the ...
Definition: astronomy.F90:2025
subroutine, public constants_init
dummy routine.
Definition: constants.F90:141
subroutine, public get_date_julian(time, year, month, day, hour, minute, second)
real, dimension(:,:), allocatable solar_ann
annual mean solar factor
Definition: astronomy.F90:390