Ignore:
Timestamp:
Jan 22, 2019, 5:02:23 PM (6 years ago)
Author:
mvals
Message:

Mars GCM:

  • Update of rocketduststorm_mod.F90 :

We want to separate both parametrizations related to the formation of the detached dust layers. Therefore, rocketduststorm_mod.F90 now only comprises the rocket dust storm scheme, whereas it contained
also before the calculation of the vertical velocity induced by the presence of the sub-grid scale topography. This latter part is under development and will be integrated as a fully independant
parametrization: the aim is to simulate the entrainment of dust by slope winds, from the boundary layer up to the top of the sub-grid scale topography.

  • Addition of initial surface parameters "summit" and "base" to prepare the previously described slope wind parametrization

MV

Location:
trunk/LMDZ.MARS/libf/dynphy_lonlat/phymars
Files:
2 edited

Legend:

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

    r1974 r2079  
    22      SUBROUTINE datareadnc(relief,phisinit,alb,ith,z0,
    33     &                    zmea,zstd,zsig,zgam,zthe,
    4      &                    hmons,summit,zavg)
     4     &                    hmons,summit,base,zavg)
    55c=======================================================================
    66c
     
    7979      REAL,intent(out) :: zthe(imdp1*jmdp1)
    8080      REAL,intent(out) :: hmons(imdp1*jmdp1) !CW17,361*180 hmons
    81       REAL,intent(out) :: summit(imdp1*jmdp1) !CW17,361*180 hmons
    82       REAL,intent(out) :: zavg(imdp1*jmdp1) !CW17,361*180 hmons
     81      REAL,intent(out) :: summit(imdp1*jmdp1)
     82      REAL,intent(out) :: base(imdp1*jmdp1)
     83      REAL,intent(out) :: zavg(imdp1*jmdp1)
    8384     
    8485! Local variables:
     
    327328c     we don't need to put them into restartfi.nc
    328329      string(6) = 'summit' !subgrid summit
    329 c     string(7) = 'base'   !base of summit (not subgrid)
     330      string(7) = 'base'   !base of summit (not subgrid)
    330331c     string(8) = 'bottom' !subgrid base
    331       do k=5,6
     332      do k=5,7
    332333          write(*,*) 'string',k,string(k) 
    333334c-----------------------------------------------------------------------
     
    373374          enddo
    374375
    375           if (k .eq. 5 .or. k .eq. 6) then
     376          if (k .eq. 5 .or. k .eq. 6 .or. k .eq. 7) then
    376377              ! hmons, summit, keep the maximum of subgrids
    377378              call mvc_horiz(imd,jmdp1,iim,jjm,longitude,
     
    404405          if (k.eq.6) then  ! summit
    405406             do i=1,iimp1*jjp1
     407                  if (pfield(i) .ne. -999999.) then             
    406408                    summit(i) = pfield(i)
     409                  else
     410                    summit(i)=0.
     411                  endif
     412             enddo
     413          endif
     414
     415          if (k.eq.7) then  ! base
     416             do i=1,iimp1*jjp1
     417                  if (pfield(i) .ne. -999999.) then             
     418                    base(i) = pfield(i)
     419                  else
     420                    base(i)=0.
     421                  endif
    407422             enddo
    408423          endif
  • trunk/LMDZ.MARS/libf/dynphy_lonlat/phymars/newstart.F

    r2001 r2079  
    2727      use surfdat_h, only: phisfi, z0, zmea, zstd, zsig, zgam, zthe,
    2828     &                     albedodat, z0_default, qsurf, tsurf,
    29      &                     co2ice, emis, hmons
     29     &                     co2ice, emis, hmons, summit, base
    3030      use comsoil_h, only: inertiedat, layer, mlayer, nsoilmx, tsoil
    3131      use control_mod, only: day_step, iphysiq, anneeref, planet_type
     
    100100      real hmonsS(iip1,jjp1)
    101101      real summitS(iip1,jjp1)
     102      real baseS(iip1,jjp1)
    102103      real zavgS(iip1,jjp1)
    103104      real z0S(iip1,jjp1)
     
    385386      CALL datareadnc(relief,phis,alb,surfith,z0S,
    386387     &          zmeaS,zstdS,zsigS,zgamS,ztheS,
    387      &          hmonsS,summitS,zavgS)
     388     &          hmonsS,summitS,baseS,zavgS)
    388389
    389390      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,phis,phisfi)
     
    398399      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,ztheS,zthe)
    399400      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,hmonsS,hmons)
    400 !      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,summitS,summit)
     401      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,summitS,summit)
     402      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,baseS,base)
    401403!      CALL gr_dyn_fi(1,iip1,jjp1,ngridmx,zavgS,zavg)
    402404
     
    16511653     &              nqtot,dtphys,real(day_ini),0.0,cell_area,
    16521654     &              albfi,ithfi,zmea,zstd,zsig,zgam,zthe,
    1653      &              hmons)
     1655     &              hmons,summit,base)
    16541656      call physdem1("restartfi.nc",nsoilmx,ngridmx,llm,nqtot,
    16551657     &              dtphys,hour_ini,
Note: See TracChangeset for help on using the changeset viewer.