Ignore:
Timestamp:
Nov 21, 2025, 4:28:06 PM (5 weeks ago)
Author:
jbclement
Message:

Mars PCM:
Simplification of the orbital parameters initialization (resolve ticket #109): removing redundant actions in subroutine 'iniorbit' in regard of what is done in "tabfi.F" and moving the subroutine 'lsp2solp' from "tabfi.F" to "planete_h.F90" + some cleanings.
JBC

File:
1 edited

Legend:

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

    r3902 r3974  
    6161      use time_phylmdz_mod, only: daysec, dtphys
    6262      use planete_h, only: aphelie, emin_turb, lmixmin, obliquit,
    63      &                     peri_day, periheli, year_day
     63     &                     peri_day, periheli, year_day,
     64     &                     lsp2solp, iniorbit
     65
    6466      implicit none
    6567 
     
    460462          write(*,*) 'perihelion =',periheli
    461463          write(*,*) 'and',year_day,'sols/year'
    462           call lsp2solp(peri_ls,peri_day,aphelie,periheli,year_day)
     464          call lsp2solp(peri_ls,peri_day)
    463465          write(*,*) 'peri_day (new value):',peri_day
    464466
     
    543545
    544546c-----------------------------------------------------------------------
     547c       Initialization of orbital parameters
     548c-----------------------------------------------------------------------
     549      call iniorbit()
     550
     551c-----------------------------------------------------------------------
    545552c       Case when using a start file from before March 1996 (without iceradius...
    546553c-----------------------------------------------------------------------
     
    555562       end if
    556563
    557 c-----------------------------------------------------------------------
    558564      END SUBROUTINE tabfi
    559565
    560 
    561 
    562      
    563 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    564 ! gives sol at perihelion for ls at perihelion (for precession cycles)
    565       subroutine lsp2solp(lsp,solp,aphelie,periheli,year_day)
    566 
    567       implicit none
    568 !  Arguments:
    569       real lsp     ! Input: ls at perihelion
    570       real solp    ! Output: sol at perihelion
    571       real aphelie,periheli,year_day ! Input: parameters
    572  
    573 !  Local:
    574       double precision zx0 ! eccentric anomaly at Ls=0
    575       double precision e_elips
    576       double precision pi,degrad
    577      
    578       parameter (pi=4.d0*atan(1.d0))
    579       parameter (degrad=180.d0/pi)
    580 
    581       e_elips=(aphelie-periheli)/(aphelie+periheli)     
    582       zx0 = -2.0*datan(dtan(0.5*lsp/degrad)
    583      .          *dsqrt((1.-e_elips)/(1.+e_elips)))
    584       if (zx0 .le. 0.) zx0 = zx0 + 2.*pi
    585      
    586       solp  = year_day*(1.-(zx0-e_elips*dsin(zx0))/(2.*pi))
    587 
    588 
    589       end subroutine lsp2solp
    590 
    591 
    592566      END MODULE tabfi_mod
    593 
Note: See TracChangeset for help on using the changeset viewer.