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

module atmosphere / SOLARINS More...

Functions/Subroutines

subroutine, public atmos_solarins_setup (basepoint_lon, basepoint_lat, iyear)
 setup solar incidence module More...
 
subroutine, public atmos_solarins_orbit (iyear)
 setup solar incidence module More...
 
subroutine, public atmos_solarins_ecliptic_longitude (Re_factor, sinDEC, cosDEC, hourangle, now_date, offset_year, lambda_out)
 calc factor of Earths solar insolation More...
 
subroutine atmos_solarins_insolation_0d (real_lon, real_lat, now_date, offset_year, solins, cosSZA, Re_factor_out)
 calc factor of Earths solar insolation More...
 
subroutine atmos_solarins_insolation_2d (IA, IS, IE, JA, JS, JE, real_lon, real_lat, now_date, offset_year, solins, cosSZA)
 calc factor of Earths solar insolation More...
 

Variables

real(rp), public atmos_solarins_constant = 1360.250117_RP
 
logical, public atmos_solarins_set_ve = .false.
 
logical, public atmos_solarins_set_ideal = .false.
 
real(rp), public atmos_solarins_obliquity = 23.44_RP
 
real(rp), public atmos_solarins_eccentricity = 0.0167_RP
 
real(rp), public atmos_solarins_perihelion_lon = 282.94719_RP
 
integer, dimension(6), public atmos_solarins_ve_date
 
real(rp), public atmos_solarins_diurnal_sec
 
real(rp), public atmos_solarins_annual_sec
 
logical, public atmos_solarins_fixedlatlon = .false.
 
real(rp), public atmos_solarins_lon
 
real(rp), public atmos_solarins_lat
 
logical, public atmos_solarins_fixeddate = .false.
 
integer, dimension(6), public atmos_solarins_date
 

Detailed Description

module atmosphere / 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
NAMELIST
  • PARAM_ATMOS_SOLARINS
    nametypedefault valuecomment
    ATMOS_SOLARINS_CONSTANT real(RP) 1360.250117_RP Solar constant [W/m2]
    ATMOS_SOLARINS_SET_VE logical .false. Set vernal equinox condition?
    ATMOS_SOLARINS_SET_IDEAL logical .false. Set obliquity and eccentricity?
    ATMOS_SOLARINS_OBLIQUITY real(RP) 23.44_RP Obliquity [deg] cf. https://nssdc.gsfc.nasa.gov/planetary/factsheet/earthfact.html
    ATMOS_SOLARINS_ECCENTRICITY real(RP) 0.0167_RP Eccentricity
    ATMOS_SOLARINS_FIXEDLATLON logical .false. Latitude/Longitude is fixed?
    ATMOS_SOLARINS_LON real(RP) Longitude for radiation [rad]
    ATMOS_SOLARINS_LAT real(RP) Latitude for radiation [rad]
    ATMOS_SOLARINS_FIXEDDATE logical .false. Date is fixed?
    ATMOS_SOLARINS_DATE integer, dimension(6) Date for radiation [Y,M,D,H,M,S]
    ATMOS_SOLARINS_PERIHELION_LON real(RP) 282.94719_RP Longitude of perihelion [deg].
    ATMOS_SOLARINS_VE_DATE integer, dimension(6) Date of first vernal equinox
    ATMOS_SOLARINS_ANNUAL_SEC real(RP) Seconds of the annual period [sec]
    ATMOS_SOLARINS_DIURNAL_SEC real(RP) Seconds of the diurnal period [sec]
    DEBUG logical .false.

History Output
No history output

Function/Subroutine Documentation

◆ atmos_solarins_setup()

subroutine, public scale_atmos_solarins::atmos_solarins_setup ( real(rp), intent(in)  basepoint_lon,
real(rp), intent(in)  basepoint_lat,
integer, intent(in)  iyear 
)

setup solar incidence module

Definition at line 570 of file scale_atmos_solarins.F90.

570  use scale_prc, only: &
571  prc_abort
572  use scale_const, only: &
573  const_d2r, &
575  use scale_calendar, only: &
576  calendar_doi, &
577  calendar_hour, &
578  calendar_min, &
579  calendar_sec, &
583  use scale_time, only: &
585  time_nowday, &
587  implicit none
588  integer, intent(in) :: iyear ! year at setup
589  real(RP), intent(in) :: basepoint_lon
590  real(RP), intent(in) :: basepoint_lat
591 
592  namelist / param_atmos_solarins / &
593  atmos_solarins_constant, &
594  atmos_solarins_set_ve, &
595  atmos_solarins_set_ideal, &
596  atmos_solarins_obliquity, &
597  atmos_solarins_eccentricity, &
599  atmos_solarins_lon, &
603  atmos_solarins_perihelion_lon, &
604  atmos_solarins_ve_date, &
605  atmos_solarins_annual_sec, &
606  atmos_solarins_diurnal_sec, &
607  debug
608 
609  real(RP) :: dyear ! delta t [year]
610  integer :: year ! used year
611  integer :: d
612  integer :: date(6)
613  real(RP) :: Re_factor, sinDEC, cosDEC, hourangle
614  real(DP) :: lambda, subsec
615  character(len=27) :: datechar
616 
617  integer :: ierr
618  !---------------------------------------------------------------------------
619 
620  log_newline
621  log_info("ATMOS_SOLARINS_setup",*) 'Setup'
622 
623  atmos_solarins_lon = basepoint_lon
624  atmos_solarins_lat = basepoint_lat
625  atmos_solarins_date(:) = -1
626  atmos_solarins_diurnal_sec = const_undef
627  atmos_solarins_annual_sec = const_undef
628 
629  !--- read namelist
630  rewind(io_fid_conf)
631  read(io_fid_conf,nml=param_atmos_solarins,iostat=ierr)
632  if( ierr < 0 ) then !--- missing
633  log_info("ATMOS_SOLARINS_setup",*) 'Not found namelist. Default used.'
634  elseif( ierr > 0 ) then !--- fatal error
635  log_error("ATMOS_SOLARINS_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_SOLARINS. Check!'
636  call prc_abort
637  endif
638  log_nml(param_atmos_solarins)
639 
640  if ( atmos_solarins_set_ve ) then
641  if ( .NOT. atmos_solarins_set_ideal ) then
642  atmos_solarins_set_ideal = .true.
643  atmos_solarins_obliquity = 0.0_rp
644  atmos_solarins_eccentricity = 0.0_rp
645  atmos_solarins_perihelion_lon = 0.0_rp
646  endif
647  if ( .NOT. atmos_solarins_fixedlatlon ) then
649  atmos_solarins_lon = 0.0_rp
650  atmos_solarins_lat = 0.0_rp
651  endif
652  if ( .NOT. atmos_solarins_fixeddate ) then
653  atmos_solarins_fixeddate = .true.
654  atmos_solarins_date(1) = ve_date(1)
655  atmos_solarins_date(2) = ve_date(2)
656  atmos_solarins_date(3) = ve_date(3)
657  atmos_solarins_date(4) = 12
658  atmos_solarins_date(5) = 0
659  atmos_solarins_date(6) = 0
660  endif
661  endif
662 
663  arcsec2d = 1.0_rp / (60.0_rp*60.0_rp)
664  arcsec2r = arcsec2d * const_d2r
665 
666  year = iyear
667 
668  if ( atmos_solarins_fixeddate ) then
669  if( atmos_solarins_date(1) >= 0 ) year = atmos_solarins_date(1)
670  endif
671 
672  if ( atmos_solarins_set_ideal ) then
673  ve_date = atmos_solarins_ve_date
674  calyearday = calendar_doi
676 
677  if (atmos_solarins_annual_sec == const_undef ) then
678  solyearsec = calyearday*caldaysec
679  atmos_solarins_annual_sec = solyearsec
680  else
681  solyearsec = atmos_solarins_annual_sec
682  year_calender_synced = .false.
683  endif
684 
685  if (atmos_solarins_diurnal_sec == const_undef ) then
686  soldaysec = caldaysec
687  atmos_solarins_diurnal_sec = soldaysec
688  else
689  soldaysec = atmos_solarins_diurnal_sec
690  day_calender_synced = .false.
691  endif
692 
693  else
694  calyearday = 365.2425_dp
696  solyearsec = calyearday*caldaysec
697  soldaysec = caldaysec
698  endif
699 
700  call atmos_solarins_orbit( year )
701 
702  dyear = real( year-year_ref, kind=rp )
703 
704  !----- report data -----
705  log_newline
706  log_info("ATMOS_SOLARINS_setup",'(1x,A)') 'Solar insolation parameters '
707  log_info("ATMOS_SOLARINS_setup",'(1x,A,I7)') 'Reference year : ', year_ref
708  log_info("ATMOS_SOLARINS_setup",'(1x,A,I7)') 'Current year : ', year
709  log_info("ATMOS_SOLARINS_setup",'(1x,A,I7)') 'Difference from ref. : ', int(dyear)
710  log_newline
711  log_info("ATMOS_SOLARINS_setup",'(1x,A,F12.7)') 'Solar constant [W/m2] : ', atmos_solarins_constant
712  log_info("ATMOS_SOLARINS_setup",'(1x,A,F12.7)') 'Obliquity [deg] : ', obliquity / const_d2r
713  log_info("ATMOS_SOLARINS_setup",'(1x,A,F12.7)') 'Eccentricity : ', e
714  log_info("ATMOS_SOLARINS_setup",'(1x,A,F12.7)') 'Longitude of perihelion [deg] : ', omega / const_d2r
715 ! LOG_INFO("ATMOS_SOLARINS_setup",'(1x,A,F12.7)') 'Longitude at the vernal equinox associated with circle orbit [deg] : ', lambda_m0 / CONST_D2R
716  log_newline
717  if ( year_calender_synced ) then
718  log_info("ATMOS_SOLARINS_setup",'(1x,A )') 'Solar annual cycle is synced with calendar'
719  else
720  log_info("ATMOS_SOLARINS_setup",'(1x,A )') 'Solar annual cycle is not synced with calendar'
721  endif
722  log_info("ATMOS_SOLARINS_setup",'(1x,A,F14.4)') 'Solar annual cycle [sec] : ', solyearsec
723  log_info("ATMOS_SOLARINS_setup",'(1x,A,F14.4)') 'Solar annual cycle [day] : ', solyearsec/caldaysec
724  if ( day_calender_synced ) then
725  log_info("ATMOS_SOLARINS_setup",'(1x,A )') 'Solar diurnal cycle is synced with calendar'
726  else
727  log_info("ATMOS_SOLARINS_setup",'(1x,A )') 'Solar diurnal cycle is not synced with calendar'
728  endif
729  log_info("ATMOS_SOLARINS_setup",'(1x,A,F14.4)') 'Solar diurnal cycle [sec] : ', soldaysec
730  log_newline
731  log_info("ATMOS_SOLARINS_setup",*) 'Latitude/Longitude is fixed? : ', atmos_solarins_fixedlatlon
732  if ( atmos_solarins_fixedlatlon ) then
733  log_info("ATMOS_SOLARINS_setup",*) 'Longitude [deg] : ', atmos_solarins_lon
734  log_info("ATMOS_SOLARINS_setup",*) 'Latitude [deg] : ', atmos_solarins_lat
735  atmos_solarins_lon = atmos_solarins_lon * const_d2r
737  endif
738 
739  log_info("ATMOS_SOLARINS_setup",*) 'Date is fixed? : ', atmos_solarins_fixeddate
740  if ( atmos_solarins_fixeddate ) then
741  log_info("ATMOS_SOLARINS_setup",*) 'Date : ', atmos_solarins_date
742  endif
743 
744  log_newline
745  log_info("ATMOS_SOLARINS_setup",*) 'Orbit info for next 1 solar year'
746  if (io_l) write(io_fid_log, '(A)') 'DATE SOLAR_CONST, Scaled dist from star, Longitude'
747  do d = time_nowday, time_nowday + int(solyearsec/caldaysec)
749  call atmos_solarins_ecliptic_longitude( &
750  re_factor, & ![OUT]
751  sindec, & ![OUT]
752  cosdec, & ![OUT]
753  hourangle, & ![OUT]
754  date, & ![IN]
755  time_offset_year, & ![IN]
756  lambda_out = lambda & ![optional OUT]
757  )
758  call calendar_date2char(datechar, date, subsec)
759  if (lambda < 0.0_rp) then
760  lambda = lambda/const_d2r + 360.0_rp
761  else
762  lambda = lambda/const_d2r
763  endif
764  if (io_l) then
765  write(io_fid_log,'(2A,F11.4,F23.6,F11.4)') datechar, " ", re_factor*atmos_solarins_constant, 1.0_dp/sqrt(re_factor), lambda
766  endif
767  enddo
768 
769  return

References atmos_solarins_annual_sec, atmos_solarins_constant, atmos_solarins_date, atmos_solarins_diurnal_sec, atmos_solarins_eccentricity, atmos_solarins_ecliptic_longitude(), atmos_solarins_fixeddate, atmos_solarins_fixedlatlon, atmos_solarins_lat, atmos_solarins_lon, atmos_solarins_obliquity, atmos_solarins_orbit(), atmos_solarins_perihelion_lon, atmos_solarins_set_ideal, atmos_solarins_set_ve, atmos_solarins_ve_date, scale_calendar::calendar_date2char(), scale_calendar::calendar_date2daysec(), scale_calendar::calendar_daysec2date(), scale_calendar::calendar_doi, scale_calendar::calendar_hour, scale_calendar::calendar_min, scale_calendar::calendar_sec, scale_const::const_d2r, scale_const::const_undef, scale_io::io_fid_conf, scale_io::io_fid_log, scale_io::io_l, scale_prc::prc_abort(), scale_precision::rp, scale_time::time_nowday, scale_time::time_nowsec, and scale_time::time_offset_year.

Referenced by mod_atmos_driver::atmos_driver_setup().

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 777 of file scale_atmos_solarins.F90.

777  use scale_const, only: &
778  pi => const_pi, &
779  d2r => const_d2r,&
781  implicit none
782 
783  integer, intent(in) :: iyear ! year at setup
784 
785  real(RP) :: dyear ! delta t [year]
786  real(RP) :: EsinOMG ! e * sin(w) : ecliptic parameter
787  real(RP) :: EcosOMG ! e * cos(w) : ecliptic parameter
788  real(RP) :: Perih_pi ! pi : longitude of fixed perihelion [rad]
789  real(RP) :: Perih_psi ! psi : general precession [degree]
790  real(RP) :: Perih_omega ! omega_bar : longitude of perihelion [degree]
791 
792  real(RP) :: temp, beta, EcosOMG_mod
793  integer :: i
794  !---------------------------------------------------------------------------
795 
796  if ( .not. atmos_solarins_set_ideal ) then ! real Earth orbit
797  ! time from reference year(1950.0 AD)
798  dyear = real( iyear-year_ref, kind=rp )
799 
800  ! obliquity
801  temp = 0.0_rp
802  do i = 1, nobliq
803  temp = temp + obliq_amp(i)*arcsec2d * cos( obliq_rate(i)*arcsec2r*dyear + obliq_phase(i)*d2r )
804  enddo
805  obliquity = ( obliquity_ref + temp ) * d2r
806 
807  ! eccentricity
808  esinomg = 0.0_rp
809  ecosomg = 0.0_rp
810  do i = 1, neclip
811  esinomg = esinomg + eclip_amp(i) * sin( eclip_rate(i)*arcsec2r*dyear + eclip_phase(i)*d2r )
812  ecosomg = ecosomg + eclip_amp(i) * cos( eclip_rate(i)*arcsec2r*dyear + eclip_phase(i)*d2r )
813  enddo
814  e = sqrt( esinomg*esinomg + ecosomg*ecosomg )
815 
816  ! longitude of fixed perihelion
817  if ( ecosomg == 0.0_rp ) then
818  ecosomg_mod = 1.e-30_rp
819  else
820  ecosomg_mod = ecosomg
821  endif
822  perih_pi = atan(esinomg/ecosomg_mod) + pi * ( 0.5_rp - sign(0.5_rp, ecosomg) )
823 
824  ! if( abs(EcosOMG) < 1.E-8_RP ) then
825  ! if ( EsinOMG == 0.0_RP ) then
826  ! Perih_pi = 0.0_RP
827  ! elseif( EsinOMG < 0.0_RP ) then
828  ! Perih_pi = 1.5_RP * PI
829  ! elseif( EsinOMG > 0.0_RP ) then
830  ! Perih_pi = 0.5_RP * PI
831  ! endif
832 
833  ! elseif( EcosOMG < 0.0_RP ) then
834 
835  ! Perih_pi = atan(EsinOMG/EcosOMG) + PI
836 
837  ! elseif( EcosOMG > 0.0_RP ) then
838 
839  ! if ( EsinOMG < 0.0_RP ) then
840  ! Perih_pi = atan(EsinOMG/EcosOMG) + 2.0_RP * PI
841  ! else
842  ! Perih_pi = atan(EsinOMG/EcosOMG)
843  ! endif
844 
845  ! endif
846  perih_pi = perih_pi / d2r ! [rad]->[degree]
847 
848  ! general precession
849  temp = 0.0_rp
850  do i = 1, nprece
851  temp = temp + prece_amp(i)*arcsec2d * sin( prece_rate(i)*arcsec2r*dyear + prece_phase(i)*d2r )
852  enddo
853  perih_psi = psi_bar * arcsec2d * dyear + zeta + temp
854 
855  ! longitude of perihelion
856  perih_omega = perih_pi + perih_psi
857 
858  do i = 1, 1000
859  if ( perih_omega + 180.0_rp < 0.0_rp ) then
860  perih_omega = perih_omega + 360.0_rp
861  elseif( perih_omega + 180.0_rp >= 360.0_rp ) then
862  perih_omega = perih_omega - 360.0_rp
863  else
864  exit
865  endif
866  enddo
867 
868  ! longitude of perihelion (see berger et al.(1993))
869  omega = ( perih_omega + 180.0_rp ) * d2r
870 
871  atmos_solarins_eccentricity = e
872  atmos_solarins_obliquity = obliquity / d2r
873  atmos_solarins_perihelion_lon = omega / d2r
874 
875  else ! ideal orbit
876 
877  obliquity = atmos_solarins_obliquity * d2r
878  e = atmos_solarins_eccentricity
879  omega = atmos_solarins_perihelion_lon * d2r
880 
881  endif
882 
883  beta = sqrt( 1.0_rp - e*e )
884 
885  ! longitude at the vernal equinox associated with circle orbit
886  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) &
887  - ( 1.0_rp/4.0_rp*e*e ) * ( 1.0_rp/2.0_rp + beta ) * sin(2.0_rp*omega) &
888  + ( 1.0_rp/8.0_rp*e*e*e ) * ( 1.0_rp/3.0_rp + beta ) * sin(3.0_rp*omega) )
889 
890  return

References atmos_solarins_eccentricity, atmos_solarins_obliquity, atmos_solarins_perihelion_lon, atmos_solarins_set_ideal, scale_const::const_d2r, scale_const::const_pi, scale_const::const_undef, and scale_precision::rp.

Referenced by atmos_solarins_setup().

Here is the caller graph for this function:

◆ atmos_solarins_ecliptic_longitude()

subroutine, public scale_atmos_solarins::atmos_solarins_ecliptic_longitude ( real(rp), intent(out)  Re_factor,
real(rp), intent(out)  sinDEC,
real(rp), intent(out)  cosDEC,
real(rp), intent(out)  hourangle,
integer, dimension(6), intent(in)  now_date,
integer, intent(in)  offset_year,
real(dp), intent(out), optional  lambda_out 
)

calc factor of Earths solar insolation

Definition at line 904 of file scale_atmos_solarins.F90.

904  use scale_const, only: &
905  pi => const_pi
906  use scale_calendar, only: &
909  i_year, i_month, i_day, &
910  i_hour, i_min, i_sec
911  implicit none
912 
913  real(RP), intent(out) :: Re_factor ! factor of the distance of Earth from the sun (1/rho2)
914  real(RP), intent(out) :: sinDEC ! sin/cos(solar declination)
915  real(RP), intent(out) :: cosDEC ! sin/cos(solar declination)
916  real(RP), intent(out) :: hourangle ! hour angle: relative longitude of subsolar point (noon=0deg)
917  integer, intent(in) :: now_date(6) ! date(yyyy,mm,dd,hh,mm,ss)
918  integer, intent(in) :: offset_year ! year offset
919  real(DP), intent(out), optional :: lambda_out ! actual longitude from vernal equinox (for output)
920 
921 
922  integer :: date(6)
923  integer :: oyear
924  integer :: absday, absday_ve
925  real(DP) :: abssec, abssec_ve
926  real(DP) :: diffday, diffyr, diffsec
927 
928  real(DP) :: lambda_m ! mean longitude from vernal equinox
929  real(DP) :: lambda ! actual longitude from vernal equinox
930  real(DP) :: nu
931 
932  !---------------------------------------------------------------------------
933 
934  date(:) = now_date(:)
935  oyear = offset_year
936 
937  if ( atmos_solarins_fixeddate ) then
938  if( atmos_solarins_date(1) >= 0 ) oyear = 0
939  if( atmos_solarins_date(1) >= 0 ) date(1) = atmos_solarins_date(1)
940  if( atmos_solarins_date(2) >= 1 ) date(2) = atmos_solarins_date(2)
941  if( atmos_solarins_date(3) >= 1 ) date(3) = atmos_solarins_date(3)
942  if( atmos_solarins_date(4) >= 0 ) date(4) = atmos_solarins_date(4)
943  if( atmos_solarins_date(5) >= 0 ) date(5) = atmos_solarins_date(5)
944  if( atmos_solarins_date(6) >= 0 ) date(6) = atmos_solarins_date(6)
945  endif
946 
947  call calendar_ymd2absday( absday, & ! [OUT]
948  date(i_year), & ! [IN]
949  date(i_month), & ! [IN]
950  date(i_day), & ! [IN]
951  oyear ) ! [IN]
952 
953  call calendar_hms2abssec( abssec, & ! [OUT]
954  date(i_hour), & ! [IN]
955  date(i_min), & ! [IN]
956  date(i_sec), & ! [IN]
957  0.0_dp ) ! [IN]
958 
959  call calendar_ymd2absday( absday_ve, & ! [OUT]
960  ve_date(i_year), & ! [IN]
961  ve_date(i_month), & ! [IN]
962  ve_date(i_day), & ! [IN]
963  oyear ) ! [IN]
964 
965  call calendar_hms2abssec( abssec_ve, & ! [OUT]
966  ve_date(i_hour), & ! [IN]
967  ve_date(i_min), & ! [IN]
968  ve_date(i_sec), & ! [IN]
969  0.0_dp ) ! [IN]
970 
971  if ( year_calender_synced ) then
972  diffday = real(absday-absday_ve,kind=dp) + (abssec-abssec_ve) / caldaysec
973  diffyr = diffday / calyearday - real( nint( diffday / calyearday ),kind=dp)
974  else
975  diffsec = real(absday-absday_ve,kind=dp)*caldaysec + (abssec-abssec_ve)
976  diffyr = diffsec / solyearsec - real( nint( diffsec / solyearsec ),kind=dp)
977  endif
978 
979  lambda_m = lambda_m0 + 2.0_rp * pi * diffyr
980 
981  nu = lambda_m - omega
982 
983  ! 1 / (rho*rho)
984  re_factor = 1.0_dp &
985  + 2.0_dp*e * cos( nu ) &
986  + 0.5_dp*e*e * cos( 2.0_dp*nu ) &
987  + 2.5_dp*e*e &
988  + 4.0_dp*e*e*e * cos( 3.0_dp*nu )
989 
990  ! actual longitude from vernal equinox
991  lambda = lambda_m &
992  + ( 2.0_dp*e - 1.0_dp/ 4.0_dp*e*e*e ) * sin( nu ) &
993  + ( 5.0_dp/ 4.0_dp*e*e ) * sin( 2.0_dp*nu ) &
994  + ( 13.0_dp/12.0_dp*e*e*e ) * sin( 3.0_dp*nu )
995 
996  ! solar declination
997  sindec = sin(lambda) * sin(obliquity)
998  cosdec = sqrt( 1.0_rp - sindec*sindec )
999 
1000  ! hour angle
1001  if ( day_calender_synced ) then
1002  hourangle = pi + 2.0_rp * pi * abssec / caldaysec
1003  else
1004  hourangle = pi + 2.0_rp * pi &
1005  * ( (absday*caldaysec + abssec) / atmos_solarins_diurnal_sec &
1006  - real( nint((absday*caldaysec + abssec)/ soldaysec ),kind=dp) )
1007  endif
1008 
1009  if ( debug ) then
1010  log_info("ATMOS_SOLARINS_ecliptic_longitude",*) lambda_m, nu, lambda, re_factor
1011  endif
1012 
1013  if ( present(lambda_out) ) lambda_out = lambda
1014 
1015  return

References atmos_solarins_date, atmos_solarins_diurnal_sec, atmos_solarins_fixeddate, scale_calendar::calendar_hms2abssec(), scale_calendar::calendar_ymd2absday(), scale_const::const_pi, scale_precision::dp, 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.

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

Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_solarins_insolation_0d()

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

calc factor of Earths solar insolation

Definition at line 1025 of file scale_atmos_solarins.F90.

1025  use scale_const, only: &
1026  eps => const_eps
1027  implicit none
1028  real(RP), intent(in) :: real_lon ! longitude
1029  real(RP), intent(in) :: real_lat ! latitude
1030  integer, intent(in) :: now_date(6) ! date(yyyy,mm,dd,hh,mm,ss)
1031  integer, intent(in) :: offset_year ! year offset
1032 
1033  real(RP), intent(out) :: solins ! solar insolation
1034  real(RP), intent(out) :: cosSZA ! cos(Solar Zenith Angle)
1035  real(RP), intent(out), optional :: Re_factor_out ! factor of the distance of Earth from the sun (1/rho2) [for output]
1036 
1037  real(RP) :: Re_factor ! factor of the distance of Earth from the sun (1/rho2)
1038  real(RP) :: sinDEC, cosDEC ! sin/cos(solar declination)
1039  real(RP) :: hourangle ! hour angle: relative longitude of subsolar point (noon=0deg)
1040 
1041  real(RP) :: lon
1042  real(RP) :: lat
1043  !---------------------------------------------------------------------------
1044 
1045  call atmos_solarins_ecliptic_longitude( re_factor, & ! [OUT]
1046  sindec, & ! [OUT]
1047  cosdec, & ! [OUT]
1048  hourangle, & ! [OUT]
1049  now_date(:), & ! [IN]
1050  offset_year ) ! [IN]
1051 
1052  if ( atmos_solarins_fixedlatlon ) then
1053  lon = atmos_solarins_lon
1054  lat = atmos_solarins_lat
1055  else
1056  lon = real_lon
1057  lat = real_lat
1058  endif
1059 
1060  cossza = sin(lat)*sindec + cos(lat)*cosdec*cos(lon+hourangle)
1061  solins = atmos_solarins_constant * re_factor * ( 0.5_rp + sign(0.5_rp,cossza-eps) )
1062 
1063  if ( present( re_factor_out ) ) then
1064  re_factor_out = re_factor
1065  endif
1066 
1067  return

References atmos_solarins_constant, atmos_solarins_ecliptic_longitude(), atmos_solarins_fixedlatlon, atmos_solarins_lat, atmos_solarins_lon, and scale_const::const_eps.

Here is the call graph for this function:

◆ atmos_solarins_insolation_2d()

subroutine scale_atmos_solarins::atmos_solarins_insolation_2d ( integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
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,
real(rp), dimension (ia,ja), intent(out)  solins,
real(rp), dimension (ia,ja), intent(out)  cosSZA 
)

calc factor of Earths solar insolation

Definition at line 1077 of file scale_atmos_solarins.F90.

1077  use scale_const, only: &
1078  eps => const_eps
1079  implicit none
1080  integer, intent(in) :: IA, IS, IE
1081  integer, intent(in) :: JA, JS, JE
1082 
1083  real(RP), intent(in) :: real_lon(IA,JA) ! longitude [rad]
1084  real(RP), intent(in) :: real_lat(IA,JA) ! latitude [rad]
1085  integer, intent(in) :: now_date(6) ! date(yyyy,mm,dd,hh,mm,ss)
1086  integer, intent(in) :: offset_year ! year offset
1087 
1088  real(RP), intent(out) :: solins (IA,JA) ! solar insolation
1089  real(RP), intent(out) :: cosSZA (IA,JA) ! cos(Solar Zenith Angle)
1090 
1091  real(RP) :: Re_factor ! factor of the distance of Earth from the sun (1/rho2)
1092  real(RP) :: sinDEC, cosDEC ! sin/cos(solar declination)
1093  real(RP) :: hourangle ! hour angle: relative longitude of subsolar point (noon=0deg)
1094 
1095  real(RP) :: lon(IA,JA)
1096  real(RP) :: lat(IA,JA)
1097 
1098  integer :: i, j
1099  !---------------------------------------------------------------------------
1100 
1101  call atmos_solarins_ecliptic_longitude( re_factor, & ! [OUT]
1102  sindec, & ! [OUT]
1103  cosdec, & ! [OUT]
1104  hourangle, & ! [OUT]
1105  now_date(:), & ! [IN]
1106  offset_year ) ! [IN]
1107 
1108  !$acc data create(lon, lat)
1109 
1110  if ( atmos_solarins_fixedlatlon ) then
1111  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1112  !$acc kernels
1113  do j = js, je
1114  do i = is, ie
1115  lon(i,j) = atmos_solarins_lon
1116  lat(i,j) = atmos_solarins_lat
1117  enddo
1118  enddo
1119  !$acc end kernels
1120  else
1121  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1122  !$acc kernels
1123  do j = js, je
1124  do i = is, ie
1125  lon(i,j) = real_lon(i,j)
1126  lat(i,j) = real_lat(i,j)
1127  enddo
1128  enddo
1129  !$acc end kernels
1130  endif
1131 
1132  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1133  !$acc kernels
1134  do j = js, je
1135  do i = is, ie
1136  cossza(i,j) = sin(lat(i,j))*sindec + cos(lat(i,j))*cosdec*cos(lon(i,j)+hourangle)
1137  solins(i,j) = atmos_solarins_constant * re_factor * ( 0.5_rp + sign(0.5_rp,cossza(i,j)-eps) )
1138  enddo
1139  enddo
1140  !$acc end kernels
1141 
1142  !$acc end data
1143 
1144  return

References atmos_solarins_constant, atmos_solarins_ecliptic_longitude(), atmos_solarins_fixedlatlon, atmos_solarins_lat, atmos_solarins_lon, and scale_const::const_eps.

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 43 of file scale_atmos_solarins.F90.

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

Referenced by atmos_solarins_insolation_0d(), atmos_solarins_insolation_2d(), atmos_solarins_setup(), mod_rm_driver::rm_driver(), mod_rm_prep::rm_prep(), and scale_atmos_phy_rd_mm5sw::swrad().

◆ atmos_solarins_set_ve

logical, public scale_atmos_solarins::atmos_solarins_set_ve = .false.

Definition at line 45 of file scale_atmos_solarins.F90.

45  logical, public :: ATMOS_SOLARINS_set_ve = .false. ! Set vernal equinox condition?

Referenced by atmos_solarins_setup().

◆ atmos_solarins_set_ideal

logical, public scale_atmos_solarins::atmos_solarins_set_ideal = .false.

Definition at line 47 of file scale_atmos_solarins.F90.

47  logical, public :: ATMOS_SOLARINS_set_ideal = .false. ! Set obliquity and eccentricity?

Referenced by atmos_solarins_orbit(), and atmos_solarins_setup().

◆ atmos_solarins_obliquity

real(rp), public scale_atmos_solarins::atmos_solarins_obliquity = 23.44_RP

Definition at line 48 of file scale_atmos_solarins.F90.

48  real(RP), public :: ATMOS_SOLARINS_obliquity = 23.44_rp ! Obliquity [deg] cf. https://nssdc.gsfc.nasa.gov/planetary/factsheet/earthfact.html

Referenced by atmos_solarins_orbit(), and atmos_solarins_setup().

◆ atmos_solarins_eccentricity

real(rp), public scale_atmos_solarins::atmos_solarins_eccentricity = 0.0167_RP

Definition at line 49 of file scale_atmos_solarins.F90.

49  real(RP), public :: ATMOS_SOLARINS_eccentricity = 0.0167_rp ! Eccentricity

Referenced by atmos_solarins_orbit(), and atmos_solarins_setup().

◆ atmos_solarins_perihelion_lon

real(rp), public scale_atmos_solarins::atmos_solarins_perihelion_lon = 282.94719_RP

Definition at line 50 of file scale_atmos_solarins.F90.

50  real(RP), public :: ATMOS_SOLARINS_perihelion_lon = 282.94719_rp ! Longitude of perihelion [deg].

Referenced by atmos_solarins_orbit(), and atmos_solarins_setup().

◆ atmos_solarins_ve_date

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

Definition at line 51 of file scale_atmos_solarins.F90.

51  integer, public :: ATMOS_SOLARINS_ve_date(6) ! Date of first vernal equinox

Referenced by atmos_solarins_setup().

◆ atmos_solarins_diurnal_sec

real(rp), public scale_atmos_solarins::atmos_solarins_diurnal_sec

Definition at line 53 of file scale_atmos_solarins.F90.

53  real(RP), public :: ATMOS_SOLARINS_diurnal_sec ! Seconds of the diurnal period [sec]

Referenced by atmos_solarins_ecliptic_longitude(), and atmos_solarins_setup().

◆ atmos_solarins_annual_sec

real(rp), public scale_atmos_solarins::atmos_solarins_annual_sec

Definition at line 54 of file scale_atmos_solarins.F90.

54  real(RP), public :: ATMOS_SOLARINS_annual_sec ! Seconds of the annual period [sec]

Referenced by atmos_solarins_setup().

◆ atmos_solarins_fixedlatlon

logical, public scale_atmos_solarins::atmos_solarins_fixedlatlon = .false.

Definition at line 56 of file scale_atmos_solarins.F90.

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

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

◆ atmos_solarins_lon

real(rp), public scale_atmos_solarins::atmos_solarins_lon

Definition at line 57 of file scale_atmos_solarins.F90.

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

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

◆ atmos_solarins_lat

real(rp), public scale_atmos_solarins::atmos_solarins_lat

Definition at line 58 of file scale_atmos_solarins.F90.

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

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

◆ atmos_solarins_fixeddate

logical, public scale_atmos_solarins::atmos_solarins_fixeddate = .false.

Definition at line 60 of file scale_atmos_solarins.F90.

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

Referenced by scale_atmos_phy_rd_profile::atmos_phy_rd_profile_read(), atmos_solarins_ecliptic_longitude(), and atmos_solarins_setup().

◆ atmos_solarins_date

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

Definition at line 61 of file scale_atmos_solarins.F90.

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

Referenced by scale_atmos_phy_rd_profile::atmos_phy_rd_profile_read(), atmos_solarins_ecliptic_longitude(), and atmos_solarins_setup().

scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_calendar::i_min
integer, parameter, public i_min
[index] minute
Definition: scale_calendar.F90:49
scale_calendar::calendar_daysec2date
subroutine, public calendar_daysec2date(ymdhms, subsec, absday, abssec, offset_year)
Convert from gregorian date to absolute day/second.
Definition: scale_calendar.F90:255
scale_atmos_solarins::atmos_solarins_fixedlatlon
logical, public atmos_solarins_fixedlatlon
Definition: scale_atmos_solarins.F90:56
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_calendar::calendar_hms2abssec
subroutine, public calendar_hms2abssec(abssec, hour, minute, second, subsec)
Hour, minute, second, subsecond -> absolute second.
Definition: scale_calendar.F90:384
scale_calendar
module CALENDAR
Definition: scale_calendar.F90:13
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_time::time_nowsec
real(dp), public time_nowsec
subday part of current time [sec]
Definition: scale_time.F90:70
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_calendar::calendar_min
real(dp), public calendar_min
minutes of hour
Definition: scale_calendar.F90:54
scale_calendar::calendar_ymd2absday
subroutine, public calendar_ymd2absday(absday, gyear, gmonth, gday, oyear)
Convert from gregorian date to absolute day, DAY 0 is AD1/1/1.
Definition: scale_calendar.F90:287
scale_calendar::i_hour
integer, parameter, public i_hour
[index] hour
Definition: scale_calendar.F90:48
scale_calendar::calendar_date2daysec
subroutine, public calendar_date2daysec(absday, abssec, ymdhms, subsec, offset_year)
Convert from gregorian date to absolute day/second.
Definition: scale_calendar.F90:192
scale_calendar::calendar_sec
real(dp), public calendar_sec
seconds of minute
Definition: scale_calendar.F90:55
scale_atmos_solarins::atmos_solarins_lat
real(rp), public atmos_solarins_lat
Definition: scale_atmos_solarins.F90:58
scale_time
module TIME
Definition: scale_time.F90:11
scale_const::const_pi
real(rp), parameter, public const_pi
pi
Definition: scale_const.F90:32
scale_calendar::i_month
integer, parameter, public i_month
[index] month
Definition: scale_calendar.F90:46
scale_calendar::calendar_hour
real(dp), public calendar_hour
hours of day
Definition: scale_calendar.F90:53
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:33
scale_calendar::i_year
integer, parameter, public i_year
[index] year
Definition: scale_calendar.F90:45
scale_calendar::i_day
integer, parameter, public i_day
[index] day
Definition: scale_calendar.F90:47
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_solarins::atmos_solarins_date
integer, dimension(6), public atmos_solarins_date
Definition: scale_atmos_solarins.F90:61
scale_calendar::calendar_doi
real(dp), public calendar_doi
days of year
Definition: scale_calendar.F90:52
scale_time::time_offset_year
integer, public time_offset_year
time offset [year]
Definition: scale_time.F90:76
scale_time::time_nowday
integer, public time_nowday
absolute day of current time [day]
Definition: scale_time.F90:69
scale_calendar::calendar_date2char
subroutine, public calendar_date2char(chardate, ymdhms, subsec)
Convert from gregorian date to absolute day/second.
Definition: scale_calendar.F90:665
scale_calendar::i_sec
integer, parameter, public i_sec
[index] second
Definition: scale_calendar.F90:50
scale_atmos_solarins::atmos_solarins_fixeddate
logical, public atmos_solarins_fixeddate
Definition: scale_atmos_solarins.F90:60