SCALE-RM
Data Types | Functions/Subroutines | Variables
scale_atmos_solarins Module Reference

module SOLARINS More...

Functions/Subroutines

subroutine, public atmos_solarins_setup (iyear)
 setup solar incidence module More...
 
subroutine, public atmos_solarins_orbit (iyear)
 setup solar incidence module More...
 
subroutine atmos_solarins_insolation_0d (solins, cosSZA, Re_factor, real_lon, real_lat, now_date, offset_year)
 calc factor of Earths solar insolation More...
 
subroutine atmos_solarins_insolation_2d (solins, cosSZA, real_lon, real_lat, now_date, offset_year)
 calc factor of Earths solar insolation More...
 

Variables

real(rp), public atmos_solarins_constant = 1360.250117_RP
 
logical, public atmos_solarins_fixedlatlon = .false.
 
logical, public atmos_solarins_fixeddate = .false.
 
real(rp), public atmos_solarins_lon
 
real(rp), public atmos_solarins_lat
 
integer, dimension(6), public atmos_solarins_date
 

Detailed Description

module SOLARINS

Description
calculate solar insolation based on berger(1978), berger et al.(1993) The original algorithm is valid to +/- 1,000,000 years from AD1950
Author
Team SCALE
History
  • 2013-01-29 (H.Yashiro) [new]
NAMELIST
  • PARAM_ATMOS_SOLARINS
    nametypedefault valuecomment
    ATMOS_SOLARINS_CONSTANT real(RP) 1360.250117_RP Solar constant [W/m2]
    ATMOS_SOLARINS_FIXEDLATLON logical .false. Latitude/Longitude is fixed?
    ATMOS_SOLARINS_FIXEDDATE logical .false. Date is fixed?
    ATMOS_SOLARINS_LON real(RP) Longitude for radiation [rad]
    ATMOS_SOLARINS_LAT real(RP) Latitude for radiation [rad]
    ATMOS_SOLARINS_DATE integer, dimension(6) Date for radiation [Y,M,D,H,M,S]

History Output
No history output

Function/Subroutine Documentation

◆ atmos_solarins_setup()

subroutine, public scale_atmos_solarins::atmos_solarins_setup ( integer, intent(in)  iyear)

setup solar incidence module

Definition at line 547 of file scale_atmos_sub_solarins.F90.

References atmos_solarins_constant, atmos_solarins_date, atmos_solarins_fixeddate, atmos_solarins_fixedlatlon, atmos_solarins_lat, atmos_solarins_lon, atmos_solarins_orbit(), scale_const::const_d2r, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, scale_process::prc_mpistop(), scale_grid_real::real_basepoint_lat, and scale_grid_real::real_basepoint_lon.

Referenced by mod_atmos_driver::atmos_driver_setup().

547  use scale_process, only: &
549  use scale_const, only: &
550  const_d2r
551  use scale_grid_real, only: &
554  implicit none
555 
556  integer, intent(in) :: iyear ! year at setup
557 
558  namelist / param_atmos_solarins / &
559  atmos_solarins_constant, &
562  atmos_solarins_lon, &
565 
566  real(RP) :: dyear ! delta t [year]
567  integer :: year ! used year
568 
569  integer :: ierr
570  !---------------------------------------------------------------------------
571 
572  if( io_l ) write(io_fid_log,*)
573  if( io_l ) write(io_fid_log,*) '++++++ Module[SOLARINS] / Categ[ATMOS SHARE] / Origin[SCALElib]'
574 
575  atmos_solarins_lon = real_basepoint_lon
577  atmos_solarins_date(:) = -1
578 
579  !--- read namelist
580  rewind(io_fid_conf)
581  read(io_fid_conf,nml=param_atmos_solarins,iostat=ierr)
582  if( ierr < 0 ) then !--- missing
583  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
584  elseif( ierr > 0 ) then !--- fatal error
585  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_SOLARINS. Check!'
586  call prc_mpistop
587  endif
588  if( io_nml ) write(io_fid_nml,nml=param_atmos_solarins)
589 
590  arcsec2d = 1.0_rp / (60.0_rp*60.0_rp)
591  arcsec2r = arcsec2d * const_d2r
592 
593  year = iyear
594 
595  if ( atmos_solarins_fixeddate ) then
596  if( atmos_solarins_date(1) >= 0 ) year = atmos_solarins_date(1)
597  endif
598 
599  call atmos_solarins_orbit( year )
600 
601  dyear = real( year-year_ref, kind=rp )
602 
603  !----- report data -----
604  if( io_l ) write(io_fid_log,*)
605  if( io_l ) write(io_fid_log,'(1x,A)') '*** Solar insolation parameters ***'
606  if( io_l ) write(io_fid_log,'(1x,A,I7)') '*** Reference year : ', year_ref
607  if( io_l ) write(io_fid_log,'(1x,A,I7)') '*** Current year : ', iyear
608  if( io_l ) write(io_fid_log,'(1x,A,I7)') '*** Difference from ref. : ', int(dyear)
609  if( io_l ) write(io_fid_log,*)
610  if( io_l ) write(io_fid_log,'(1x,A,F12.7)') '*** Solar constant [W/m2] : ', atmos_solarins_constant
611  if( io_l ) write(io_fid_log,'(1x,A,F12.7)') '*** Obliquity [deg] : ', obliquity / const_d2r
612  if( io_l ) write(io_fid_log,'(1x,A,F12.7)') '*** Eccentricity : ', e
613  if( io_l ) write(io_fid_log,'(1x,A,F12.7)') '*** Longitude of perihelion [deg] : ', omega / const_d2r
614  if( io_l ) write(io_fid_log,'(1x,A,F12.7)') '*** Longitude at the vernal equinox [deg] : ', lambda_m0 / const_d2r
615  if( io_l ) write(io_fid_log,*)
616 
617  if( io_l ) write(io_fid_log,*) '*** Latitude/Longitude is fixed? : ', atmos_solarins_fixedlatlon
618  if ( atmos_solarins_fixedlatlon ) then
619  if( io_l ) write(io_fid_log,*) '*** Longitude [deg] : ', atmos_solarins_lon
620  if( io_l ) write(io_fid_log,*) '*** Latitude [deg] : ', atmos_solarins_lat
621  atmos_solarins_lon = atmos_solarins_lon * const_d2r
623  endif
624 
625  if( io_l ) write(io_fid_log,*) '*** Date is fixed? : ', atmos_solarins_fixeddate
626  if ( atmos_solarins_fixeddate ) then
627  if( io_l ) write(io_fid_log,*) '*** Date : ', atmos_solarins_date
628  endif
629 
630  return
subroutine, public prc_mpistop
Abort MPI.
real(rp), public real_basepoint_lat
position of base point in real world [rad,-pi,pi]
logical, public atmos_solarins_fixeddate
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
logical, public atmos_solarins_fixedlatlon
module GRID (real space)
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, dimension(6), public atmos_solarins_date
real(rp), public real_basepoint_lon
position of base point in real world [rad,0-2pi]
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_solarins_orbit()

subroutine, public scale_atmos_solarins::atmos_solarins_orbit ( integer, intent(in)  iyear)

setup solar incidence module

Definition at line 638 of file scale_atmos_sub_solarins.F90.

References scale_const::const_d2r, and scale_const::const_pi.

Referenced by atmos_solarins_setup().

638  use scale_const, only: &
639  pi => const_pi, &
640  const_d2r
641  implicit none
642 
643  integer, intent(in) :: iyear ! year at setup
644 
645  real(RP) :: dyear ! delta t [year]
646 
647  real(RP) :: EsinOMG ! e * sin(w) : ecliptic parameter
648  real(RP) :: EcosOMG ! e * cos(w) : ecliptic parameter
649  real(RP) :: Perih_pi ! pi : longitude of fixed perihelion [rad]
650  real(RP) :: Perih_psi ! psi : general precession [degree]
651  real(RP) :: Perih_omega ! omega_bar : longitude of perihelion [degree]
652 
653  real(RP) :: temp, beta, EcosOMG_mod
654  integer :: i
655  !---------------------------------------------------------------------------
656 
657  ! time from reference year(1950.0 AD)
658  dyear = real( iyear-year_ref, kind=rp )
659 
660  ! obliquity
661  temp = 0.0_rp
662  do i = 1, nobliq
663  temp = temp + obliq_amp(i)*arcsec2d * cos( obliq_rate(i)*arcsec2r*dyear + obliq_phase(i)*const_d2r )
664  enddo
665  obliquity = ( obliquity_ref + temp ) * const_d2r
666 
667  ! eccentricity
668  esinomg = 0.0_rp
669  ecosomg = 0.0_rp
670  do i = 1, neclip
671  esinomg = esinomg + eclip_amp(i) * sin( eclip_rate(i)*arcsec2r*dyear + eclip_phase(i)*const_d2r )
672  ecosomg = ecosomg + eclip_amp(i) * cos( eclip_rate(i)*arcsec2r*dyear + eclip_phase(i)*const_d2r )
673  enddo
674  e = sqrt( esinomg*esinomg + ecosomg*ecosomg )
675 
676  ! longitude of fixed perihelion
677  if ( ecosomg == 0.0_rp ) then
678  ecosomg_mod = 1.e-30_rp
679  else
680  ecosomg_mod = ecosomg
681  endif
682  perih_pi = atan(esinomg/ecosomg_mod) + pi * ( 0.5_rp - sign(0.5_rp, ecosomg) )
683 
684 ! if( abs(EcosOMG) < 1.E-8_RP ) then
685 ! if ( EsinOMG == 0.0_RP ) then
686 ! Perih_pi = 0.0_RP
687 ! elseif( EsinOMG < 0.0_RP ) then
688 ! Perih_pi = 1.5_RP * PI
689 ! elseif( EsinOMG > 0.0_RP ) then
690 ! Perih_pi = 0.5_RP * PI
691 ! endif
692 !
693 ! elseif( EcosOMG < 0.0_RP ) then
694 !
695 ! Perih_pi = atan(EsinOMG/EcosOMG) + PI
696 !
697 ! elseif( EcosOMG > 0.0_RP ) then
698 !
699 ! if ( EsinOMG < 0.0_RP ) then
700 ! Perih_pi = atan(EsinOMG/EcosOMG) + 2.0_RP * PI
701 ! else
702 ! Perih_pi = atan(EsinOMG/EcosOMG)
703 ! endif
704 !
705 ! endif
706  perih_pi = perih_pi / const_d2r ! [rad]->[degree]
707 
708  ! general precession
709  temp = 0.0_rp
710  do i = 1, nprece
711  temp = temp + prece_amp(i)*arcsec2d * sin( prece_rate(i)*arcsec2r*dyear + prece_phase(i)*const_d2r )
712  enddo
713  perih_psi = psi_bar * arcsec2d * dyear + zeta + temp
714 
715  ! longitude of perihelion
716  perih_omega = perih_pi + perih_psi
717 
718  do i = 1, 1000
719  if ( perih_omega + 180.0_rp < 0.0_rp ) then
720  perih_omega = perih_omega + 360.0_rp
721  elseif( perih_omega + 180.0_rp >= 360.0_rp ) then
722  perih_omega = perih_omega - 360.0_rp
723  else
724  exit
725  endif
726  enddo
727 
728  ! longitude of perigee (see berger et al.(1993))
729  omega = ( perih_omega + 180.0_rp ) * const_d2r
730 
731  beta = sqrt( 1.0_rp - e*e )
732 
733  ! longitude at the vernal equinox
734  lambda_m0 = 2.0_rp * ( ( 1.0_rp/2.0_rp*e + 1.0_rp/8.0_rp*e*e*e ) * ( 1.0_rp + beta ) * sin( omega) &
735  - ( 1.0_rp/4.0_rp*e*e ) * ( 1.0_rp/2.0_rp + beta ) * sin(2.0_rp*omega) &
736  + ( 1.0_rp/8.0_rp*e*e*e ) * ( 1.0_rp/3.0_rp + beta ) * sin(3.0_rp*omega) )
737 
738  return
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
module CONSTANT
Definition: scale_const.F90:14
real(rp), public const_pi
pi
Definition: scale_const.F90:34
Here is the caller graph for this function:

◆ atmos_solarins_insolation_0d()

subroutine scale_atmos_solarins::atmos_solarins_insolation_0d ( real(rp), intent(out)  solins,
real(rp), intent(out)  cosSZA,
real(rp), intent(out)  Re_factor,
real(rp), intent(in)  real_lon,
real(rp), intent(in)  real_lat,
integer, dimension(6), intent(in)  now_date,
integer, intent(in)  offset_year 
)

calc factor of Earths solar insolation

Definition at line 751 of file scale_atmos_sub_solarins.F90.

References atmos_solarins_constant, atmos_solarins_date, atmos_solarins_fixeddate, atmos_solarins_fixedlatlon, atmos_solarins_lat, atmos_solarins_lon, scale_calendar::calendar_getdayofyear(), scale_calendar::calendar_hms2abssec(), scale_calendar::calendar_ymd2absday(), scale_const::const_eps, scale_const::const_pi, scale_calendar::i_day, scale_calendar::i_hour, scale_calendar::i_min, scale_calendar::i_month, scale_calendar::i_sec, and scale_calendar::i_year.

751  use scale_const, only: &
752  pi => const_pi, &
753  eps => const_eps
754  use scale_calendar, only: &
758  i_year, i_month, i_day, &
759  i_hour, i_min, i_sec
760  implicit none
761 
762  real(RP), intent(out) :: solins ! solar insolation
763  real(RP), intent(out) :: cosSZA ! cos(Solar Zenith Angle)
764  real(RP), intent(out) :: Re_factor ! factor of the distance of Earth from the sun (1/rho2)
765  real(RP), intent(in) :: real_lon ! longitude
766  real(RP), intent(in) :: real_lat ! latitude
767  integer, intent(in) :: now_date(6) ! date(yyyy,mm,dd,hh,mm,ss)
768  integer, intent(in) :: offset_year ! year offset
769 
770  real(RP) :: lon
771  real(RP) :: lat
772  integer :: date(6)
773  integer :: oyear
774 
775  real(RP) :: lambda_m ! mean longitude from vernal equinox
776  real(RP) :: lambda !
777  real(RP) :: sinDEC, cosDEC ! sin/cos(solar declination)
778  real(RP) :: hourangle ! hour angle: relative longitude of subsolar point
779 
780  integer :: absday, absday_ve
781  real(DP) :: DayOfYear, abssec
782  real(RP) :: nu
783  !---------------------------------------------------------------------------
784 
785  lon = real_lon
786  lat = real_lat
787  if ( atmos_solarins_fixedlatlon ) then
788  lon = atmos_solarins_lon
789  lat = atmos_solarins_lat
790  endif
791 
792  date(:) = now_date(:)
793  oyear = offset_year
794  if ( atmos_solarins_fixeddate ) then
795  if( atmos_solarins_date(1) >= 0 ) oyear = 0
796  if( atmos_solarins_date(1) >= 0 ) date(1) = atmos_solarins_date(1)
797  if( atmos_solarins_date(2) >= 1 ) date(2) = atmos_solarins_date(2)
798  if( atmos_solarins_date(3) >= 1 ) date(3) = atmos_solarins_date(3)
799  if( atmos_solarins_date(4) >= 0 ) date(4) = atmos_solarins_date(4)
800  if( atmos_solarins_date(5) >= 0 ) date(5) = atmos_solarins_date(5)
801  if( atmos_solarins_date(6) >= 0 ) date(6) = atmos_solarins_date(6)
802  endif
803 
804  call calendar_getdayofyear( dayofyear, date(i_year) )
805 
806  call calendar_ymd2absday( absday, & ! [OUT]
807  date(i_year), & ! [IN]
808  date(i_month), & ! [IN]
809  date(i_day), & ! [IN]
810  oyear ) ! [IN]
811 
812  call calendar_ymd2absday( absday_ve, & ! [OUT]
813  date(i_year), & ! [IN]
814  ve_date(i_month), & ! [IN]
815  ve_date(i_day), & ! [IN]
816  oyear ) ! [IN]
817 
818  call calendar_hms2abssec( abssec, & ! [OUT]
819  date(i_hour), & ! [IN]
820  date(i_min), & ! [IN]
821  date(i_sec), & ! [IN]
822  0.0_dp ) ! [IN]
823 
824  lambda_m = lambda_m0 + 2.0_rp * pi * real(absday-absday_ve,kind=RP) / DayOfYear
825 
826  nu = lambda_m - omega
827 
828  ! actual longitude from vernal equinox
829  lambda = lambda_m &
830  + ( 2.0_rp*e - 1.0_rp/ 4.0_rp*e*e*e ) * sin( nu ) &
831  + ( 5.0_rp/ 4.0_rp*e*e ) * sin( 2.0_rp*nu ) &
832  + ( 13.0_rp/12.0_rp*e*e*e ) * sin( 3.0_rp*nu )
833 
834  ! solar declination
835  sindec = sin(lambda) * sin(obliquity)
836  cosdec = sqrt( 1.0_rp - sindec*sindec )
837 
838  ! hour angle
839  hourangle = lon + 2.0_rp * pi * abssec / (24.0_rp*60.0_rp*60.0_rp)
840 
841  ! cos(Solar Zenith Angle)
842  cossza = sin(lat)*sindec - cos(lat)*cosdec*cos(hourangle)
843 
844  ! 1 / (rho*rho)
845  re_factor = 1.0_rp &
846  + 2.0_rp*e * cos( nu ) &
847  + 0.5_rp*e*e * cos( 2.0_rp*nu ) &
848  + 2.5_rp*e*e &
849  + 4.0_rp*e*e*e * cos( 3.0_rp*nu )
850 
851 ! solins = ATMOS_SOLARINS_constant * Re_factor * max(cosSZA,0.0_RP)
852  solins = atmos_solarins_constant * re_factor * ( 0.5_rp + sign(0.5_rp,cossza-eps) )
853 
854  return
integer, parameter, public i_month
[index] month
integer, parameter, public i_year
[index] year
subroutine, public calendar_ymd2absday(absday, gyear, gmonth, gday, oyear)
Convert from gregorian date to absolute day, DAY 0 is AD1/1/1.
logical, public atmos_solarins_fixeddate
logical, public atmos_solarins_fixedlatlon
integer, parameter, public i_min
[index] minute
integer, parameter, public i_sec
[index] second
module CONSTANT
Definition: scale_const.F90:14
integer, parameter, public i_hour
[index] hour
integer, parameter, public i_day
[index] day
real(rp), public const_eps
small number
Definition: scale_const.F90:36
real(rp), public const_pi
pi
Definition: scale_const.F90:34
module CALENDAR
integer, dimension(6), public atmos_solarins_date
subroutine, public calendar_getdayofyear(DayOfYear, iyear)
Get day of year.
subroutine, public calendar_hms2abssec(abssec, hour, minute, second, subsec)
Hour, minute, second, subsecond -> absolute second.
Here is the call graph for this function:

◆ atmos_solarins_insolation_2d()

subroutine scale_atmos_solarins::atmos_solarins_insolation_2d ( real(rp), dimension (ia,ja), intent(out)  solins,
real(rp), dimension (ia,ja), intent(out)  cosSZA,
real(rp), dimension(ia,ja), intent(in)  real_lon,
real(rp), dimension(ia,ja), intent(in)  real_lat,
integer, dimension(6), intent(in)  now_date,
integer, intent(in)  offset_year 
)

calc factor of Earths solar insolation

Definition at line 866 of file scale_atmos_sub_solarins.F90.

References atmos_solarins_constant, atmos_solarins_date, atmos_solarins_fixeddate, atmos_solarins_fixedlatlon, atmos_solarins_lat, atmos_solarins_lon, scale_calendar::calendar_getdayofyear(), scale_calendar::calendar_hms2abssec(), scale_calendar::calendar_ymd2absday(), scale_const::const_eps, scale_const::const_pi, scale_calendar::i_day, scale_calendar::i_hour, scale_calendar::i_min, scale_calendar::i_month, scale_calendar::i_sec, scale_calendar::i_year, scale_grid_index::ie, scale_grid_index::is, scale_grid_index::je, and scale_grid_index::js.

866  use scale_grid_index
867  use scale_const, only: &
868  pi => const_pi, &
869  eps => const_eps
870  use scale_calendar, only: &
874  i_year, i_month, i_day, &
875  i_hour, i_min, i_sec
876  implicit none
877 
878  real(RP), intent(out) :: solins (IA,JA) ! solar insolation
879  real(RP), intent(out) :: cosSZA (IA,JA) ! cos(Solar Zenith Angle)
880  real(RP), intent(in) :: real_lon(IA,JA) ! longitude [rad]
881  real(RP), intent(in) :: real_lat(IA,JA) ! latitude [rad]
882  integer, intent(in) :: now_date(6) ! date(yyyy,mm,dd,hh,mm,ss)
883  integer, intent(in) :: offset_year ! year offset
884 
885  real(RP) :: lon(IA,JA)
886  real(RP) :: lat(IA,JA)
887  integer :: date(6)
888  integer :: oyear
889 
890  real(RP) :: lambda_m ! mean longitude from vernal equinox
891  real(RP) :: lambda !
892  real(RP) :: sinDEC, cosDEC ! sin/cos(solar declination)
893  real(RP) :: hourangle(IA,JA) ! hour angle: relative longitude of subsolar point
894 
895  integer :: absday, absday_ve
896  real(DP) :: DayOfYear, abssec
897  real(RP) :: nu
898  real(RP) :: Re_factor ! factor of the distance of Earth from the sun (1/rho2)
899 
900  integer :: i, j
901  !---------------------------------------------------------------------------
902 
903  lon(:,:) = real_lon(:,:)
904  lat(:,:) = real_lat(:,:)
905  if ( atmos_solarins_fixedlatlon ) then
906  lon(:,:) = atmos_solarins_lon
907  lat(:,:) = atmos_solarins_lat
908  endif
909 
910  date(:) = now_date(:)
911  oyear = offset_year
912  if ( atmos_solarins_fixeddate ) then
913  if( atmos_solarins_date(1) >= 0 ) oyear = 0
914  if( atmos_solarins_date(1) >= 0 ) date(1) = atmos_solarins_date(1)
915  if( atmos_solarins_date(2) >= 1 ) date(2) = atmos_solarins_date(2)
916  if( atmos_solarins_date(3) >= 1 ) date(3) = atmos_solarins_date(3)
917  if( atmos_solarins_date(4) >= 0 ) date(4) = atmos_solarins_date(4)
918  if( atmos_solarins_date(5) >= 0 ) date(5) = atmos_solarins_date(5)
919  if( atmos_solarins_date(6) >= 0 ) date(6) = atmos_solarins_date(6)
920  endif
921 
922  call calendar_getdayofyear( dayofyear, date(i_year) )
923 
924  call calendar_ymd2absday( absday, & ! [OUT]
925  date(i_year), & ! [IN]
926  date(i_month), & ! [IN]
927  date(i_day), & ! [IN]
928  oyear ) ! [IN]
929 
930  call calendar_ymd2absday( absday_ve, & ! [OUT]
931  date(i_year), & ! [IN]
932  ve_date(i_month), & ! [IN]
933  ve_date(i_day), & ! [IN]
934  oyear ) ! [IN]
935 
936  call calendar_hms2abssec( abssec, & ! [OUT]
937  date(i_hour), & ! [IN]
938  date(i_min), & ! [IN]
939  date(i_sec), & ! [IN]
940  0.0_dp ) ! [IN]
941 
942  lambda_m = lambda_m0 + 2.0_rp * pi * real(absday-absday_ve,kind=RP) / DayOfYear
943 
944  nu = lambda_m - omega
945 
946  ! 1 / (rho*rho)
947  re_factor = 1.0_rp &
948  + 2.0_rp*e * cos( nu ) &
949  + 0.5_rp*e*e * cos( 2.0_rp*nu ) &
950  + 2.5_rp*e*e &
951  + 4.0_rp*e*e*e * cos( 3.0_rp*nu )
952 
953  ! actual longitude from vernal equinox
954  lambda = lambda_m &
955  + ( 2.0_rp*e - 1.0_rp/ 4.0_rp*e*e*e ) * sin( nu ) &
956  + ( 5.0_rp/ 4.0_rp*e*e ) * sin( 2.0_rp*nu ) &
957  + ( 13.0_rp/12.0_rp*e*e*e ) * sin( 3.0_rp*nu )
958 
959  ! solar declination
960  sindec = sin(lambda) * sin(obliquity)
961  cosdec = sqrt( 1.0_rp - sindec*sindec )
962 
963  do j = js, je
964  do i = is, ie
965  ! hour angle
966  hourangle(i,j) = lon(i,j) + 2.0_rp * pi * abssec / (24.0_rp*60.0_rp*60.0_rp)
967  enddo
968  enddo
969 
970  do j = js, je
971  do i = is, ie
972  ! cos(Solar Zenith Angle)
973  cossza(i,j) = sin(lat(i,j))*sindec - cos(lat(i,j))*cosdec*cos(hourangle(i,j))
974  enddo
975  enddo
976 
977  do j = js, je
978  do i = is, ie
979 ! solins(i,j) = ATMOS_SOLARINS_constant * Re_factor * max(cosSZA(i,j),0.0_RP)
980  solins(i,j) = atmos_solarins_constant * re_factor * ( 0.5_rp + sign(0.5_rp,cossza(i,j)-eps) )
981  enddo
982  enddo
983 
984  return
integer, public is
start point of inner domain: x, local
integer, parameter, public i_month
[index] month
integer, parameter, public i_year
[index] year
integer, public je
end point of inner domain: y, local
subroutine, public calendar_ymd2absday(absday, gyear, gmonth, gday, oyear)
Convert from gregorian date to absolute day, DAY 0 is AD1/1/1.
logical, public atmos_solarins_fixeddate
module grid index
logical, public atmos_solarins_fixedlatlon
integer, parameter, public i_min
[index] minute
integer, public js
start point of inner domain: y, local
integer, parameter, public i_sec
[index] second
module CONSTANT
Definition: scale_const.F90:14
integer, parameter, public i_hour
[index] hour
integer, parameter, public i_day
[index] day
integer, public ie
end point of inner domain: x, local
real(rp), public const_eps
small number
Definition: scale_const.F90:36
real(rp), public const_pi
pi
Definition: scale_const.F90:34
module CALENDAR
integer, dimension(6), public atmos_solarins_date
subroutine, public calendar_getdayofyear(DayOfYear, iyear)
Get day of year.
subroutine, public calendar_hms2abssec(abssec, hour, minute, second, subsec)
Hour, minute, second, subsecond -> absolute second.
Here is the call graph for this function:

Variable Documentation

◆ atmos_solarins_constant

real(rp), public scale_atmos_solarins::atmos_solarins_constant = 1360.250117_RP

Definition at line 44 of file scale_atmos_sub_solarins.F90.

Referenced by atmos_solarins_insolation_0d(), atmos_solarins_insolation_2d(), atmos_solarins_setup(), and scale_atmos_phy_rd_mm5sw::swrad().

44  real(RP), public :: ATMOS_SOLARINS_constant = 1360.250117_rp ! Solar constant [W/m2]

◆ atmos_solarins_fixedlatlon

logical, public scale_atmos_solarins::atmos_solarins_fixedlatlon = .false.

Definition at line 45 of file scale_atmos_sub_solarins.F90.

Referenced by scale_atmos_phy_rd_profile::atmos_phy_rd_profile_read(), atmos_solarins_insolation_0d(), atmos_solarins_insolation_2d(), and atmos_solarins_setup().

45  logical, public :: ATMOS_SOLARINS_fixedlatlon = .false. ! Latitude/Longitude is fixed?

◆ atmos_solarins_fixeddate

logical, public scale_atmos_solarins::atmos_solarins_fixeddate = .false.

Definition at line 46 of file scale_atmos_sub_solarins.F90.

Referenced by scale_atmos_phy_rd_profile::atmos_phy_rd_profile_read(), atmos_solarins_insolation_0d(), atmos_solarins_insolation_2d(), and atmos_solarins_setup().

46  logical, public :: ATMOS_SOLARINS_fixeddate = .false. ! Date is fixed?

◆ atmos_solarins_lon

real(rp), public scale_atmos_solarins::atmos_solarins_lon

Definition at line 47 of file scale_atmos_sub_solarins.F90.

Referenced by atmos_solarins_insolation_0d(), atmos_solarins_insolation_2d(), and atmos_solarins_setup().

47  real(RP), public :: ATMOS_SOLARINS_lon ! Longitude for radiation [rad]

◆ atmos_solarins_lat

real(rp), public scale_atmos_solarins::atmos_solarins_lat

Definition at line 48 of file scale_atmos_sub_solarins.F90.

Referenced by scale_atmos_phy_rd_profile::atmos_phy_rd_profile_read(), atmos_solarins_insolation_0d(), atmos_solarins_insolation_2d(), and atmos_solarins_setup().

48  real(RP), public :: ATMOS_SOLARINS_lat ! Latitude for radiation [rad]

◆ atmos_solarins_date

integer, dimension(6), public scale_atmos_solarins::atmos_solarins_date

Definition at line 49 of file scale_atmos_sub_solarins.F90.

Referenced by scale_atmos_phy_rd_profile::atmos_phy_rd_profile_read(), atmos_solarins_insolation_0d(), atmos_solarins_insolation_2d(), and atmos_solarins_setup().

49  integer, public :: ATMOS_SOLARINS_date(6) ! Date for radiation [Y,M,D,H,M,S]