Index: /trunk/LMDZ.COMMON/libf/evolution/orbit_param_criterion_mod.F90
===================================================================
--- /trunk/LMDZ.COMMON/libf/evolution/orbit_param_criterion_mod.F90	(revision 3039)
+++ /trunk/LMDZ.COMMON/libf/evolution/orbit_param_criterion_mod.F90	(revision 3040)
@@ -22,7 +22,7 @@
 #endif
 #ifndef CPP_STD
-    use planete_h,   only: e_elips, obliquit, timeperi
+    use planete_h,   only: e_elips, obliquit, lsperi
 #else
-    use planete_mod, only: e_elips, obliquit, timeperi
+    use planete_mod, only: e_elips, obliquit, lsperi
 #endif
   use time_evol_mod,  only: year_bp_ini, var_obl, var_ecc, var_lsp, convert_years
@@ -57,6 +57,6 @@
 ! 0. Initializations
 ! **********************************************************************
-Year = year_bp_ini + i_myear        ! We compute the current year
-Lsp = modulo(-timeperi*180./pi,360.) ! We convert in degree
+Year = year_bp_ini + i_myear ! We compute the current year
+Lsp = lsperi*180./pi         ! We convert in degree
 
 call ini_lask_param_mod ! Allocation of variables
Index: /trunk/LMDZ.COMMON/libf/evolution/recomp_orb_param_mod.F90
===================================================================
--- /trunk/LMDZ.COMMON/libf/evolution/recomp_orb_param_mod.F90	(revision 3039)
+++ /trunk/LMDZ.COMMON/libf/evolution/recomp_orb_param_mod.F90	(revision 3040)
@@ -15,7 +15,7 @@
 #ifndef CPP_STD      
     use comconst_mod, only: pi
-    use planete_h,    only: e_elips, obliquit, timeperi, periheli, aphelie, p_elips, peri_day, year_day
+    use planete_h,    only: e_elips, obliquit, lsperi, periheli, aphelie, p_elips, peri_day, year_day
 #else
-    use planete_mod,  only: e_elips, obliquit, timeperi, periastr, apoastr, p_elips, peri_day, year_day
+    use planete_mod,  only: e_elips, obliquit, lsperi, periastr, apoastr, p_elips, peri_day, year_day
     use comcstfi_mod, only: pi
 #endif
@@ -56,6 +56,6 @@
 #endif
 
-Year = year_bp_ini + i_myear        ! We compute the current year
-Lsp = modulo(-timeperi*180./pi,360.) ! We convert in degree
+Year = year_bp_ini + i_myear ! We compute the current year
+Lsp = lsperi*180./pi         ! We convert in degree
 
 write(*,*) 'recomp_orb_param, Old year =', Year - year_iter
@@ -100,5 +100,5 @@
 enddo
 halfaxe = 227.94
-timeperi = Lsp*pi/180.
+Lsp = Lsp*pi/180.
 periheli = halfaxe*(1 - e_elips)
 aphelie =  halfaxe*(1 + e_elips)
@@ -107,5 +107,5 @@
 
 #ifndef CPP_STD  
-    call call_dayperi(timeperi,e_elips,peri_day,year_day)
+    call call_dayperi(Lsp,e_elips,peri_day,year_day)
 #endif
 
Index: /trunk/LMDZ.MARS/changelog.txt
===================================================================
--- /trunk/LMDZ.MARS/changelog.txt	(revision 3039)
+++ /trunk/LMDZ.MARS/changelog.txt	(revision 3040)
@@ -4186,5 +4186,4 @@
 JBC
 
-
 == 11/09/2023 == JBC
 * New management of time in PEM. While definitions in "run.def"/"launching script" and Laskar's data are in Earth years, the PEM works in Martian years. It follows several changes and the new variable 'convert_years' to make the conversion;
@@ -4193,2 +4192,5 @@
 * Laskar's data interpolation has been cleaned further;
 * Some cleaning and renaming of PEM subroutines in the course of the previous modifications.
+
+== 12/09/2023 == JBC
+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.
Index: /trunk/LMDZ.MARS/libf/phymars/call_dayperi.F
===================================================================
--- /trunk/LMDZ.MARS/libf/phymars/call_dayperi.F	(revision 3039)
+++ /trunk/LMDZ.MARS/libf/phymars/call_dayperi.F	(revision 3040)
@@ -15,11 +15,11 @@
 c   Input:
 c   ------
-c   Lsperi      solar longitude (Ls) of perohelion (rad)
-c   e_elips       Excentricity
-c   real year_day      ! number of sols per Mars yar
+c   Lsperi        ! Solar longitude (Ls) of perihelion (rad)
+c   e_elips       ! Excentricity
+c   real year_day ! Number of sols per Mars year
 c
 c   output
 c   ------
-c   dayperi       Martian date at perihelion (sol)
+c   dayperi       ! Martian date at perihelion (sol)
 
 c-----------------------------------------------------------------------
@@ -27,7 +27,6 @@
 c arguments:
 c ----------
-       real Lsperi,e_elips,dayperi
-       real year_day      ! number of sols per Mars yar
-
+       real Lsperi, e_elips, dayperi
+       real year_day ! Number of sols per Mars year
 
 c      Local
@@ -43,5 +42,6 @@
      &   ( 2*atan(x1*tan(0.5*Lsperi))
      &     -x2*sin(Lsperi)/(1+e_elips*cos(Lsperi)) )
-       if(dayperi.lt.0) dayperi=dayperi+year_day
+       if (dayperi < 0) dayperi=dayperi+year_day
        return
+
        end
Index: /trunk/LMDZ.MARS/libf/phymars/iniorbit.F
===================================================================
--- /trunk/LMDZ.MARS/libf/phymars/iniorbit.F	(revision 3039)
+++ /trunk/LMDZ.MARS/libf/phymars/iniorbit.F	(revision 3040)
@@ -2,6 +2,5 @@
      $     (paphelie,pperiheli,pyear_day,pperi_day,pobliq)
       use planete_h, only: aphelie, periheli, year_day, peri_day,
-     &                     obliquit, unitastr, e_elips, p_elips,
-     &                     timeperi
+     &                     obliquit, unitastr, e_elips, p_elips, lsperi
       use comcstfi_h, only: pi
       IMPLICIT NONE
@@ -67,7 +66,8 @@
       PRINT*,'iniorbit: zx0   ',zx0
 
-      timeperi=2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.))
+      lsperi=-2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.))
+      lsperi = modulo(lsperi,2.*pi)
       PRINT*,'iniorbit: Perihelion solar long. Ls (deg)=',
-     &       modulo(-timeperi*180./pi,360.)
+     &       lsperi*180./pi
 
       END
Index: /trunk/LMDZ.MARS/libf/phymars/orbite.F
===================================================================
--- /trunk/LMDZ.MARS/libf/phymars/orbite.F	(revision 3039)
+++ /trunk/LMDZ.MARS/libf/phymars/orbite.F	(revision 3040)
@@ -1,4 +1,4 @@
       SUBROUTINE orbite(pls,pdist_sol,pdecli)
-      USE planete_h, ONLY: e_elips, p_elips, obliquit, timeperi
+      USE planete_h, ONLY: e_elips, p_elips, obliquit, lsperi
       USE comcstfi_h, ONLY: pi
       IMPLICIT NONE
@@ -19,5 +19,5 @@
 c   -------
 c   pdist_sol     Distance Sun-Planet in UA
-c   pdecli        declinaison ( in radians )
+c   pdecli        Declinaison ( in radians )
 c
 c=======================================================================
@@ -26,16 +26,14 @@
 c ----------
 
-      REAL,INTENT(IN) :: pls
-      REAL,INTENT(OUT) :: pdist_sol,pdecli
+      REAL,INTENT(IN)  :: pls
+      REAL,INTENT(OUT) :: pdist_sol, pdecli
 
 c-----------------------------------------------------------------------
 
 c Distance Sun-Planet
-
-      pdist_sol=p_elips/(1.+e_elips*cos(pls+timeperi))
+      pdist_sol = p_elips/(1.+e_elips*cos(pls-lsperi))
 
 c Solar declination
-
-      pdecli= asin (sin(pls)*sin(obliquit*pi/180.))
+      pdecli = asin(sin(pls)*sin(obliquit*pi/180.))
 
       END
Index: /trunk/LMDZ.MARS/libf/phymars/planete_h.F90
===================================================================
--- /trunk/LMDZ.MARS/libf/phymars/planete_h.F90	(revision 3039)
+++ /trunk/LMDZ.MARS/libf/phymars/planete_h.F90	(revision 3040)
@@ -2,17 +2,17 @@
       IMPLICIT NONE
 
-      REAL aphelie  ! Aphelion, in Mkm
-      REAL periheli ! Perihelion, in Mkm
-      REAL year_day ! number of days in the year
-      REAL peri_day ! date of perihelion, in days
-      REAL obliquit ! obliquity of the planet, in degrees
+      REAL aphelie   ! Aphelion, in Mkm
+      REAL periheli  ! Perihelion, in Mkm
+      REAL year_day  ! Number of days in the year
+      REAL peri_day  ! Date of perihelion, in days
+      REAL obliquit  ! Obliquity of the planet, in degrees
       REAL lmixmin
       REAL emin_turb
       REAL coefvis
       REAL coefir
-      REAL timeperi ! angle of the perihelion, in rad
-      REAL e_elips ! orbit eccentricity
-      REAL p_elips
-      REAL unitastr ! Astronomical unit AU, in Mkm
+      REAL lsperi    ! Solar longitude of the perihelion, angle in rad
+      REAL e_elips   ! Orbit eccentricity
+      REAL p_elips   ! Ellipse semi-latus rectum
+      REAL unitastr  ! Astronomical unit AU, in Mkm
 
       END MODULE planete_h
Index: /trunk/LMDZ.MARS/libf/phymars/solarlong.F
===================================================================
--- /trunk/LMDZ.MARS/libf/phymars/solarlong.F	(revision 3039)
+++ /trunk/LMDZ.MARS/libf/phymars/solarlong.F	(revision 3040)
@@ -84,5 +84,5 @@
       zteta=2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.))
 
-      psollong=zteta-timeperi
+      psollong=zteta+lsperi
 
       IF(psollong.LT.0.) psollong=psollong+2.*pi
