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 atmos_solarins_ecliptic_longitude (Re_factor, sinDEC, cosDEC, hourangle, now_date, offset_year)
 calc factor of Earths solar insolation More...
 
subroutine atmos_solarins_insolation_0d (real_lon, real_lat, now_date, offset_year, solins, cosSZA)
 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 = 0.0_RP
 
real(rp), public atmos_solarins_eccentricity = 0.0_RP
 
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) 0.0_RP Obliquity [deg]
    ATMOS_SOLARINS_ECCENTRICITY real(RP) 0.0_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]
    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 556 of file scale_atmos_solarins.F90.

556  use scale_prc, only: &
557  prc_abort
558  use scale_const, only: &
559  const_d2r
560  implicit none
561  integer, intent(in) :: iyear ! year at setup
562  real(RP), intent(in) :: basepoint_lon
563  real(RP), intent(in) :: basepoint_lat
564 
565  namelist / param_atmos_solarins / &
566  atmos_solarins_constant, &
567  atmos_solarins_set_ve, &
568  atmos_solarins_set_ideal, &
569  atmos_solarins_obliquity, &
570  atmos_solarins_eccentricity, &
572  atmos_solarins_lon, &
576  debug
577 
578  real(RP) :: dyear ! delta t [year]
579  integer :: year ! used year
580 
581  integer :: ierr
582  !---------------------------------------------------------------------------
583 
584  log_newline
585  log_info("ATMOS_SOLARINS_setup",*) 'Setup'
586 
587  atmos_solarins_lon = basepoint_lon
588  atmos_solarins_lat = basepoint_lat
589  atmos_solarins_date(:) = -1
590 
591  !--- read namelist
592  rewind(io_fid_conf)
593  read(io_fid_conf,nml=param_atmos_solarins,iostat=ierr)
594  if( ierr < 0 ) then !--- missing
595  log_info("ATMOS_SOLARINS_setup",*) 'Not found namelist. Default used.'
596  elseif( ierr > 0 ) then !--- fatal error
597  log_error("ATMOS_SOLARINS_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_SOLARINS. Check!'
598  call prc_abort
599  endif
600  log_nml(param_atmos_solarins)
601 
602  if ( atmos_solarins_set_ve ) then
603  if ( .NOT. atmos_solarins_set_ideal ) then
604  atmos_solarins_set_ideal = .true.
605  atmos_solarins_obliquity = 0.0_rp
606  atmos_solarins_eccentricity = 0.0_rp
607  endif
608  if ( .NOT. atmos_solarins_fixedlatlon ) then
610  atmos_solarins_lon = 0.0_rp
611  atmos_solarins_lat = 0.0_rp
612  endif
613  if ( .NOT. atmos_solarins_fixeddate ) then
614  atmos_solarins_fixeddate = .true.
615  atmos_solarins_date(1) = ve_date(1)
616  atmos_solarins_date(2) = ve_date(2)
617  atmos_solarins_date(3) = ve_date(3)
618  atmos_solarins_date(4) = 12
619  atmos_solarins_date(5) = 0
620  atmos_solarins_date(6) = 0
621  endif
622  endif
623 
624  arcsec2d = 1.0_rp / (60.0_rp*60.0_rp)
625  arcsec2r = arcsec2d * const_d2r
626 
627  year = iyear
628 
629  if ( atmos_solarins_fixeddate ) then
630  if( atmos_solarins_date(1) >= 0 ) year = atmos_solarins_date(1)
631  endif
632 
633  call atmos_solarins_orbit( year )
634 
635  dyear = real( year-year_ref, kind=rp )
636 
637  !----- report data -----
638  log_newline
639  log_info("ATMOS_SOLARINS_setup",'(1x,A)') 'Solar insolation parameters '
640  log_info("ATMOS_SOLARINS_setup",'(1x,A,I7)') 'Reference year : ', year_ref
641  log_info("ATMOS_SOLARINS_setup",'(1x,A,I7)') 'Current year : ', year
642  log_info("ATMOS_SOLARINS_setup",'(1x,A,I7)') 'Difference from ref. : ', int(dyear)
643  log_newline
644  log_info("ATMOS_SOLARINS_setup",'(1x,A,F12.7)') 'Solar constant [W/m2] : ', atmos_solarins_constant
645  log_info("ATMOS_SOLARINS_setup",'(1x,A,F12.7)') 'Obliquity [deg] : ', obliquity / const_d2r
646  log_info("ATMOS_SOLARINS_setup",'(1x,A,F12.7)') 'Eccentricity : ', e
647  log_info("ATMOS_SOLARINS_setup",'(1x,A,F12.7)') 'Longitude of perihelion [deg] : ', omega / const_d2r
648  log_info("ATMOS_SOLARINS_setup",'(1x,A,F12.7)') 'Longitude at the vernal equinox [deg] : ', lambda_m0 / const_d2r
649  log_newline
650 
651  log_info("ATMOS_SOLARINS_setup",*) 'Latitude/Longitude is fixed? : ', atmos_solarins_fixedlatlon
652  if ( atmos_solarins_fixedlatlon ) then
653  log_info("ATMOS_SOLARINS_setup",*) 'Longitude [deg] : ', atmos_solarins_lon
654  log_info("ATMOS_SOLARINS_setup",*) 'Latitude [deg] : ', atmos_solarins_lat
655  atmos_solarins_lon = atmos_solarins_lon * const_d2r
657  endif
658 
659  log_info("ATMOS_SOLARINS_setup",*) 'Date is fixed? : ', atmos_solarins_fixeddate
660  if ( atmos_solarins_fixeddate ) then
661  log_info("ATMOS_SOLARINS_setup",*) 'Date : ', atmos_solarins_date
662  endif
663 
664  return

References atmos_solarins_constant, atmos_solarins_date, atmos_solarins_eccentricity, atmos_solarins_fixeddate, atmos_solarins_fixedlatlon, atmos_solarins_lat, atmos_solarins_lon, atmos_solarins_obliquity, atmos_solarins_orbit(), atmos_solarins_set_ideal, atmos_solarins_set_ve, scale_const::const_d2r, scale_io::io_fid_conf, scale_prc::prc_abort(), and scale_precision::rp.

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

672  use scale_const, only: &
673  pi => const_pi, &
674  d2r => const_d2r
675  implicit none
676 
677  integer, intent(in) :: iyear ! year at setup
678 
679  real(RP) :: dyear ! delta t [year]
680  real(RP) :: EsinOMG ! e * sin(w) : ecliptic parameter
681  real(RP) :: EcosOMG ! e * cos(w) : ecliptic parameter
682  real(RP) :: Perih_pi ! pi : longitude of fixed perihelion [rad]
683  real(RP) :: Perih_psi ! psi : general precession [degree]
684  real(RP) :: Perih_omega ! omega_bar : longitude of perihelion [degree]
685 
686  real(RP) :: temp, beta, EcosOMG_mod
687  integer :: i
688  !---------------------------------------------------------------------------
689 
690  ! time from reference year(1950.0 AD)
691  dyear = real( iyear-year_ref, kind=rp )
692 
693  ! obliquity
694  temp = 0.0_rp
695  do i = 1, nobliq
696  temp = temp + obliq_amp(i)*arcsec2d * cos( obliq_rate(i)*arcsec2r*dyear + obliq_phase(i)*d2r )
697  enddo
698  obliquity = ( obliquity_ref + temp ) * d2r
699 
700  ! eccentricity
701  esinomg = 0.0_rp
702  ecosomg = 0.0_rp
703  do i = 1, neclip
704  esinomg = esinomg + eclip_amp(i) * sin( eclip_rate(i)*arcsec2r*dyear + eclip_phase(i)*d2r )
705  ecosomg = ecosomg + eclip_amp(i) * cos( eclip_rate(i)*arcsec2r*dyear + eclip_phase(i)*d2r )
706  enddo
707  e = sqrt( esinomg*esinomg + ecosomg*ecosomg )
708 
709  if ( atmos_solarins_set_ideal ) then
710  obliquity = atmos_solarins_obliquity * d2r
711  e = atmos_solarins_eccentricity
712  endif
713 
714  ! longitude of fixed perihelion
715  if ( ecosomg == 0.0_rp ) then
716  ecosomg_mod = 1.e-30_rp
717  else
718  ecosomg_mod = ecosomg
719  endif
720  perih_pi = atan(esinomg/ecosomg_mod) + pi * ( 0.5_rp - sign(0.5_rp, ecosomg) )
721 
722 ! if( abs(EcosOMG) < 1.E-8_RP ) then
723 ! if ( EsinOMG == 0.0_RP ) then
724 ! Perih_pi = 0.0_RP
725 ! elseif( EsinOMG < 0.0_RP ) then
726 ! Perih_pi = 1.5_RP * PI
727 ! elseif( EsinOMG > 0.0_RP ) then
728 ! Perih_pi = 0.5_RP * PI
729 ! endif
730 !
731 ! elseif( EcosOMG < 0.0_RP ) then
732 !
733 ! Perih_pi = atan(EsinOMG/EcosOMG) + PI
734 !
735 ! elseif( EcosOMG > 0.0_RP ) then
736 !
737 ! if ( EsinOMG < 0.0_RP ) then
738 ! Perih_pi = atan(EsinOMG/EcosOMG) + 2.0_RP * PI
739 ! else
740 ! Perih_pi = atan(EsinOMG/EcosOMG)
741 ! endif
742 !
743 ! endif
744  perih_pi = perih_pi / d2r ! [rad]->[degree]
745 
746  ! general precession
747  temp = 0.0_rp
748  do i = 1, nprece
749  temp = temp + prece_amp(i)*arcsec2d * sin( prece_rate(i)*arcsec2r*dyear + prece_phase(i)*d2r )
750  enddo
751  perih_psi = psi_bar * arcsec2d * dyear + zeta + temp
752 
753  ! longitude of perihelion
754  perih_omega = perih_pi + perih_psi
755 
756  do i = 1, 1000
757  if ( perih_omega + 180.0_rp < 0.0_rp ) then
758  perih_omega = perih_omega + 360.0_rp
759  elseif( perih_omega + 180.0_rp >= 360.0_rp ) then
760  perih_omega = perih_omega - 360.0_rp
761  else
762  exit
763  endif
764  enddo
765 
766  ! longitude of perigee (see berger et al.(1993))
767  omega = ( perih_omega + 180.0_rp ) * d2r
768 
769  beta = sqrt( 1.0_rp - e*e )
770 
771  ! longitude at the vernal equinox
772  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) &
773  - ( 1.0_rp/4.0_rp*e*e ) * ( 1.0_rp/2.0_rp + beta ) * sin(2.0_rp*omega) &
774  + ( 1.0_rp/8.0_rp*e*e*e ) * ( 1.0_rp/3.0_rp + beta ) * sin(3.0_rp*omega) )
775 
776  return

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

Referenced by atmos_solarins_setup().

Here is the caller graph for this function:

◆ atmos_solarins_ecliptic_longitude()

subroutine 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 
)

calc factor of Earths solar insolation

Definition at line 788 of file scale_atmos_solarins.F90.

788  use scale_const, only: &
789  pi => const_pi
790  use scale_calendar, only: &
794  i_year, i_month, i_day, &
795  i_hour, i_min, i_sec
796  implicit none
797 
798  real(RP), intent(out) :: Re_factor ! factor of the distance of Earth from the sun (1/rho2)
799  real(RP), intent(out) :: sinDEC ! sin/cos(solar declination)
800  real(RP), intent(out) :: cosDEC ! sin/cos(solar declination)
801  real(RP), intent(out) :: hourangle ! hour angle: relative longitude of subsolar point
802  integer, intent(in) :: now_date(6) ! date(yyyy,mm,dd,hh,mm,ss)
803  integer, intent(in) :: offset_year ! year offset
804 
805  integer :: date(6)
806  integer :: oyear
807  integer :: absday, absday_ve
808  real(DP) :: DayOfYear, abssec
809 
810  real(RP) :: lambda_m ! mean longitude from vernal equinox
811  real(RP) :: lambda ! actual longitude from vernal equinox
812  real(RP) :: nu
813  !---------------------------------------------------------------------------
814 
815  date(:) = now_date(:)
816  oyear = offset_year
817 
818  if ( atmos_solarins_fixeddate ) then
819  if( atmos_solarins_date(1) >= 0 ) oyear = 0
820  if( atmos_solarins_date(1) >= 0 ) date(1) = atmos_solarins_date(1)
821  if( atmos_solarins_date(2) >= 1 ) date(2) = atmos_solarins_date(2)
822  if( atmos_solarins_date(3) >= 1 ) date(3) = atmos_solarins_date(3)
823  if( atmos_solarins_date(4) >= 0 ) date(4) = atmos_solarins_date(4)
824  if( atmos_solarins_date(5) >= 0 ) date(5) = atmos_solarins_date(5)
825  if( atmos_solarins_date(6) >= 0 ) date(6) = atmos_solarins_date(6)
826  endif
827 
828  call calendar_getdayofyear( dayofyear, date(i_year) )
829 
830  call calendar_ymd2absday( absday, & ! [OUT]
831  date(i_year), & ! [IN]
832  date(i_month), & ! [IN]
833  date(i_day), & ! [IN]
834  oyear ) ! [IN]
835 
836  call calendar_ymd2absday( absday_ve, & ! [OUT]
837  date(i_year), & ! [IN]
838  ve_date(i_month), & ! [IN]
839  ve_date(i_day), & ! [IN]
840  oyear ) ! [IN]
841 
842  call calendar_hms2abssec( abssec, & ! [OUT]
843  date(i_hour), & ! [IN]
844  date(i_min), & ! [IN]
845  date(i_sec), & ! [IN]
846  0.0_dp ) ! [IN]
847 
848  lambda_m = lambda_m0 + 2.0_rp * pi * real(absday-absday_ve,kind=rp) / dayofyear
849 
850  nu = lambda_m - omega
851 
852  ! 1 / (rho*rho)
853  re_factor = 1.0_rp &
854  + 2.0_rp*e * cos( nu ) &
855  + 0.5_rp*e*e * cos( 2.0_rp*nu ) &
856  + 2.5_rp*e*e &
857  + 4.0_rp*e*e*e * cos( 3.0_rp*nu )
858 
859  ! actual longitude from vernal equinox
860  lambda = lambda_m &
861  + ( 2.0_rp*e - 1.0_rp/ 4.0_rp*e*e*e ) * sin( nu ) &
862  + ( 5.0_rp/ 4.0_rp*e*e ) * sin( 2.0_rp*nu ) &
863  + ( 13.0_rp/12.0_rp*e*e*e ) * sin( 3.0_rp*nu )
864 
865  ! solar declination
866  sindec = sin(lambda) * sin(obliquity)
867  cosdec = sqrt( 1.0_rp - sindec*sindec )
868 
869  ! hour angle
870  hourangle = 2.0_rp * pi * abssec / (24.0_rp*60.0_rp*60.0_rp)
871 
872  if ( debug ) then
873  log_info("ATMOS_SOLARINS_ecliptic_longitude",*) lambda_m, nu, lambda
874  endif
875 
876  return

References atmos_solarins_date, atmos_solarins_fixeddate, scale_calendar::calendar_getdayofyear(), scale_calendar::calendar_hms2abssec(), scale_calendar::calendar_ymd2absday(), 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.

Referenced by atmos_solarins_insolation_0d(), and atmos_solarins_insolation_2d().

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 
)

calc factor of Earths solar insolation

Definition at line 885 of file scale_atmos_solarins.F90.

885  use scale_const, only: &
886  eps => const_eps
887  implicit none
888  real(RP), intent(in) :: real_lon ! longitude
889  real(RP), intent(in) :: real_lat ! latitude
890  integer, intent(in) :: now_date(6) ! date(yyyy,mm,dd,hh,mm,ss)
891  integer, intent(in) :: offset_year ! year offset
892 
893  real(RP), intent(out) :: solins ! solar insolation
894  real(RP), intent(out) :: cosSZA ! cos(Solar Zenith Angle)
895 
896  real(RP) :: Re_factor ! factor of the distance of Earth from the sun (1/rho2)
897  real(RP) :: sinDEC, cosDEC ! sin/cos(solar declination)
898  real(RP) :: hourangle ! hour angle: relative longitude of subsolar point
899 
900  real(RP) :: lon
901  real(RP) :: lat
902  !---------------------------------------------------------------------------
903 
904  call atmos_solarins_ecliptic_longitude( re_factor, & ! [OUT]
905  sindec, & ! [OUT]
906  cosdec, & ! [OUT]
907  hourangle, & ! [OUT]
908  now_date(:), & ! [IN]
909  offset_year ) ! [IN]
910 
911  if ( atmos_solarins_fixedlatlon ) then
912  lon = atmos_solarins_lon
913  lat = atmos_solarins_lat
914  else
915  lon = real_lon
916  lat = real_lat
917  endif
918 
919  cossza = sin(lat)*sindec - cos(lat)*cosdec*cos(lon+hourangle)
920  solins = atmos_solarins_constant * re_factor * ( 0.5_rp + sign(0.5_rp,cossza-eps) )
921 
922  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 932 of file scale_atmos_solarins.F90.

932  use scale_const, only: &
933  eps => const_eps
934  implicit none
935  integer, intent(in) :: IA, IS, IE
936  integer, intent(in) :: JA, JS, JE
937 
938  real(RP), intent(in) :: real_lon(IA,JA) ! longitude [rad]
939  real(RP), intent(in) :: real_lat(IA,JA) ! latitude [rad]
940  integer, intent(in) :: now_date(6) ! date(yyyy,mm,dd,hh,mm,ss)
941  integer, intent(in) :: offset_year ! year offset
942 
943  real(RP), intent(out) :: solins (IA,JA) ! solar insolation
944  real(RP), intent(out) :: cosSZA (IA,JA) ! cos(Solar Zenith Angle)
945 
946  real(RP) :: Re_factor ! factor of the distance of Earth from the sun (1/rho2)
947  real(RP) :: sinDEC, cosDEC ! sin/cos(solar declination)
948  real(RP) :: hourangle ! hour angle: relative longitude of subsolar point
949 
950  real(RP) :: lon(IA,JA)
951  real(RP) :: lat(IA,JA)
952 
953  integer :: i, j
954  !---------------------------------------------------------------------------
955 
956  call atmos_solarins_ecliptic_longitude( re_factor, & ! [OUT]
957  sindec, & ! [OUT]
958  cosdec, & ! [OUT]
959  hourangle, & ! [OUT]
960  now_date(:), & ! [IN]
961  offset_year ) ! [IN]
962 
963  if ( atmos_solarins_fixedlatlon ) then
964  !$omp parallel do OMP_SCHEDULE_ collapse(2)
965  do j = js, je
966  do i = is, ie
967  lon(i,j) = atmos_solarins_lon
968  lat(i,j) = atmos_solarins_lat
969  enddo
970  enddo
971  else
972  !$omp parallel do OMP_SCHEDULE_ collapse(2)
973  do j = js, je
974  do i = is, ie
975  lon(i,j) = real_lon(i,j)
976  lat(i,j) = real_lat(i,j)
977  enddo
978  enddo
979  endif
980 
981  !$omp parallel do OMP_SCHEDULE_ collapse(2)
982  do j = js, je
983  do i = is, ie
984  cossza(i,j) = sin(lat(i,j))*sindec - cos(lat(i,j))*cosdec*cos(lon(i,j)+hourangle)
985  solins(i,j) = atmos_solarins_constant * re_factor * ( 0.5_rp + sign(0.5_rp,cossza(i,j)-eps) )
986  enddo
987  enddo
988 
989  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 42 of file scale_atmos_solarins.F90.

42  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(), and scale_atmos_phy_rd_mm5sw::swrad().

◆ atmos_solarins_set_ve

logical, public scale_atmos_solarins::atmos_solarins_set_ve = .false.

Definition at line 44 of file scale_atmos_solarins.F90.

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

46  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 = 0.0_RP

Definition at line 47 of file scale_atmos_solarins.F90.

47  real(RP), public :: ATMOS_SOLARINS_obliquity = 0.0_rp ! Obliquity [deg]

Referenced by atmos_solarins_orbit(), and atmos_solarins_setup().

◆ atmos_solarins_eccentricity

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

Definition at line 48 of file scale_atmos_solarins.F90.

48  real(RP), public :: ATMOS_SOLARINS_eccentricity = 0.0_rp ! Eccentricity

Referenced by atmos_solarins_orbit(), and atmos_solarins_setup().

◆ atmos_solarins_fixedlatlon

logical, public scale_atmos_solarins::atmos_solarins_fixedlatlon = .false.

Definition at line 50 of file scale_atmos_solarins.F90.

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

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

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

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

55  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:342
scale_calendar::i_min
integer, parameter, public i_min
[index] minute
Definition: scale_calendar.F90:49
scale_atmos_solarins::atmos_solarins_fixedlatlon
logical, public atmos_solarins_fixedlatlon
Definition: scale_atmos_solarins.F90:50
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:33
scale_calendar::calendar_getdayofyear
subroutine, public calendar_getdayofyear(DayOfYear, iyear)
Get day of year.
Definition: scale_calendar.F90:139
scale_calendar::calendar_hms2abssec
subroutine, public calendar_hms2abssec(abssec, hour, minute, second, subsec)
Hour, minute, second, subsecond -> absolute second.
Definition: scale_calendar.F90:322
scale_calendar
module CALENDAR
Definition: scale_calendar.F90:13
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_const::const_pi
real(rp), public const_pi
pi
Definition: scale_const.F90:31
scale_const
module CONSTANT
Definition: scale_const.F90:11
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:223
scale_calendar::i_hour
integer, parameter, public i_hour
[index] hour
Definition: scale_calendar.F90:48
scale_atmos_solarins::atmos_solarins_lat
real(rp), public atmos_solarins_lat
Definition: scale_atmos_solarins.F90:52
scale_calendar::i_month
integer, parameter, public i_month
[index] month
Definition: scale_calendar.F90:46
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
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_atmos_solarins::atmos_solarins_date
integer, dimension(6), public atmos_solarins_date
Definition: scale_atmos_solarins.F90:55
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:54