Changeset 3040 for trunk/LMDZ.MARS/libf


Ignore:
Timestamp:
Sep 12, 2023, 11:41:22 AM (15 months ago)
Author:
jbclement
Message:

Mars PCM:
The variable 'timeperi' (defined in "planete_h.F90" and computed in "iniorbit.F") is renamed into 'lsperi' and thus slightly changed to be coherent to the solar longitude of perihelion in radian. It can now be used out of the box by other subroutines/programs like the PEM.
JBC

Location:
trunk/LMDZ.MARS/libf/phymars
Files:
5 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/libf/phymars/call_dayperi.F

    r2963 r3040  
    1515c   Input:
    1616c   ------
    17 c   Lsperi      solar longitude (Ls) of perohelion (rad)
    18 c   e_elips       Excentricity
    19 c   real year_day      ! number of sols per Mars yar
     17c   Lsperi        ! Solar longitude (Ls) of perihelion (rad)
     18c   e_elips       ! Excentricity
     19c   real year_day ! Number of sols per Mars year
    2020c
    2121c   output
    2222c   ------
    23 c   dayperi       Martian date at perihelion (sol)
     23c   dayperi       ! Martian date at perihelion (sol)
    2424
    2525c-----------------------------------------------------------------------
     
    2727c arguments:
    2828c ----------
    29        real Lsperi,e_elips,dayperi
    30        real year_day      ! number of sols per Mars yar
    31 
     29       real Lsperi, e_elips, dayperi
     30       real year_day ! Number of sols per Mars year
    3231
    3332c      Local
     
    4342     &   ( 2*atan(x1*tan(0.5*Lsperi))
    4443     &     -x2*sin(Lsperi)/(1+e_elips*cos(Lsperi)) )
    45        if(dayperi.lt.0) dayperi=dayperi+year_day
     44       if (dayperi < 0) dayperi=dayperi+year_day
    4645       return
     46
    4747       end
  • trunk/LMDZ.MARS/libf/phymars/iniorbit.F

    r3037 r3040  
    22     $     (paphelie,pperiheli,pyear_day,pperi_day,pobliq)
    33      use planete_h, only: aphelie, periheli, year_day, peri_day,
    4      &                     obliquit, unitastr, e_elips, p_elips,
    5      &                     timeperi
     4     &                     obliquit, unitastr, e_elips, p_elips, lsperi
    65      use comcstfi_h, only: pi
    76      IMPLICIT NONE
     
    6766      PRINT*,'iniorbit: zx0   ',zx0
    6867
    69       timeperi=2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.))
     68      lsperi=-2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.))
     69      lsperi = modulo(lsperi,2.*pi)
    7070      PRINT*,'iniorbit: Perihelion solar long. Ls (deg)=',
    71      &       modulo(-timeperi*180./pi,360.)
     71     &       lsperi*180./pi
    7272
    7373      END
  • trunk/LMDZ.MARS/libf/phymars/orbite.F

    r1368 r3040  
    11      SUBROUTINE orbite(pls,pdist_sol,pdecli)
    2       USE planete_h, ONLY: e_elips, p_elips, obliquit, timeperi
     2      USE planete_h, ONLY: e_elips, p_elips, obliquit, lsperi
    33      USE comcstfi_h, ONLY: pi
    44      IMPLICIT NONE
     
    1919c   -------
    2020c   pdist_sol     Distance Sun-Planet in UA
    21 c   pdecli        declinaison ( in radians )
     21c   pdecli        Declinaison ( in radians )
    2222c
    2323c=======================================================================
     
    2626c ----------
    2727
    28       REAL,INTENT(IN) :: pls
    29       REAL,INTENT(OUT) :: pdist_sol,pdecli
     28      REAL,INTENT(IN)  :: pls
     29      REAL,INTENT(OUT) :: pdist_sol, pdecli
    3030
    3131c-----------------------------------------------------------------------
    3232
    3333c Distance Sun-Planet
    34 
    35       pdist_sol=p_elips/(1.+e_elips*cos(pls+timeperi))
     34      pdist_sol = p_elips/(1.+e_elips*cos(pls-lsperi))
    3635
    3736c Solar declination
    38 
    39       pdecli= asin (sin(pls)*sin(obliquit*pi/180.))
     37      pdecli = asin(sin(pls)*sin(obliquit*pi/180.))
    4038
    4139      END
  • trunk/LMDZ.MARS/libf/phymars/planete_h.F90

    r2228 r3040  
    22      IMPLICIT NONE
    33
    4       REAL aphelie  ! Aphelion, in Mkm
    5       REAL periheli ! Perihelion, in Mkm
    6       REAL year_day ! number of days in the year
    7       REAL peri_day ! date of perihelion, in days
    8       REAL obliquit ! obliquity of the planet, in degrees
     4      REAL aphelie   ! Aphelion, in Mkm
     5      REAL periheli  ! Perihelion, in Mkm
     6      REAL year_day  ! Number of days in the year
     7      REAL peri_day  ! Date of perihelion, in days
     8      REAL obliquit  ! Obliquity of the planet, in degrees
    99      REAL lmixmin
    1010      REAL emin_turb
    1111      REAL coefvis
    1212      REAL coefir
    13       REAL timeperi ! angle of the perihelion, in rad
    14       REAL e_elips ! orbit eccentricity
    15       REAL p_elips
    16       REAL unitastr ! Astronomical unit AU, in Mkm
     13      REAL lsperi    ! Solar longitude of the perihelion, angle in rad
     14      REAL e_elips   ! Orbit eccentricity
     15      REAL p_elips   ! Ellipse semi-latus rectum
     16      REAL unitastr  ! Astronomical unit AU, in Mkm
    1717
    1818      END MODULE planete_h
  • trunk/LMDZ.MARS/libf/phymars/solarlong.F

    r1226 r3040  
    8484      zteta=2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.))
    8585
    86       psollong=zteta-timeperi
     86      psollong=zteta+lsperi
    8787
    8888      IF(psollong.LT.0.) psollong=psollong+2.*pi
Note: See TracChangeset for help on using the changeset viewer.