Ignore:
Timestamp:
May 24, 2012, 5:02:41 PM (13 years ago)
Author:
tnavarro
Message:

Ice effective radius for output is weighted by ice total surface instead of ice total mass in physiq.F

  • New distribution for permanent ice reservoirs in surfini.F

icelocationmode = 1 allows for realistic ice distribution, no matter what resolution.
icelocationmode = 2 is predefined 32x24 or 64x48 and should be USED BY DEFAULT. It currently overestimates ice.
icelocationmode = 3 is good for quick logical definitions (paleoclimates,stability studies,etc...)

  • Possibility to change perihelion date in solar longitude as well as in sol in tabfi.F for newstarts
  • Improved latent heat release when ground ice sublimates: now a implicit scheme for ground temperature
File:
1 edited

Legend:

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

    r552 r669  
    6767      INTEGER Lmodif
    6868      INTEGER lmax
    69       REAL p_rad,p_omeg,p_g,p_mugaz,p_daysec,time
     69      REAL p_rad,p_omeg,p_g,p_mugaz,p_daysec,time,peri_ls
    7070
    7171c Variables
     
    289289      write(*,*) '(18)        obliquit : planet obliquity (deg)'
    290290      write(*,*) '(17)      peri_day : perihelion date (sol since Ls=0)'
     291      write(*,*) '(  )      peri_ls : perihelion date (Ls since Ls=0)'
    291292      write(*,*) '(15)      periheli : min. sun-mars dist (Mkm)'
    292293      write(*,*) '(16)      aphelie  : max. sun-mars dist (Mkm)'
     
    429430        else if (trim(modif) .eq. 'peri_day') then
    430431          write(*,*) 'current value:',peri_day
    431           write(*,*) 'peri_day should be 485 on current Mars'
     432          write(*,*) 'peri_day should be 485 sols on current Mars'
    432433          write(*,*) 'enter new value:'
    433434 116      read(*,*,iostat=ierr) peri_day
     
    435436          write(*,*)
    436437          write(*,*) ' peri_day (new value):',peri_day
     438         
     439        else if (trim(modif) .eq. 'peri_ls') then
     440          write(*,*) 'peri_ls value is not stored in start files,'
     441          write(*,*) 'but it should be 251 degrees on current Mars'
     442          write(*,*) '(peri_day should be 485 sols on current Mars)'
     443          write(*,*) 'enter new value:'
     444 1160     read(*,*,iostat=ierr) peri_ls
     445          if(ierr.ne.0) goto 1160
     446          write(*,*)
     447          write(*,*) ' peri_ls asked :',peri_ls
     448          call lsp2solp(peri_ls,peri_day)
     449          write(*,*) ' peri_day (new value):',peri_day
     450
    437451
    438452        else if (trim(modif) .eq. 'periheli') then
     
    529543c-----------------------------------------------------------------------
    530544      end
     545
     546
     547
     548     
     549!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
     550! gives sol at perihelion for ls at perihelion (for precession cycles)
     551      subroutine lsp2solp(lsp,solp)
     552
     553      implicit none
     554!  Arguments:
     555      real lsp     ! Input: ls at perihelion
     556      real solp    ! Output: sol at perihelion
     557     
     558!  Local:
     559      double precision zx0 ! eccentric anomaly at Ls=0
     560      double precision year_day
     561      double precision e_elips
     562      double precision pi,degrad
     563     
     564      parameter (year_day=668.6d0) ! number of sols in a martian year
     565      parameter (e_elips=0.0934d0)  ! eccentricity of orbit
     566      parameter (pi=3.14159265358979d0)
     567      parameter (degrad=57.2957795130823d0)
     568     
     569      zx0 = -2.0*datan(dtan(0.5*lsp/degrad)
     570     .          *dsqrt((1.-e_elips)/(1.+e_elips)))
     571      if (zx0 .le. 0.) zx0 = zx0 + 2.*pi
     572     
     573      solp  = year_day*(1.-(zx0-e_elips*dsin(zx0))/(2.*pi))
     574
     575
     576      end subroutine lsp2solp
     577
     578
     579
Note: See TracChangeset for help on using the changeset viewer.