Ignore:
Timestamp:
Mar 9, 2023, 9:38:24 AM (22 months ago)
Author:
romain.vande
Message:

Mars PEM:
New Boolean options for following orbital parameters of ob_ex_lsp.asc: var_obl, var_ex, var_lsp.
If using evol_orbit_pem=true, you can specify which parameter to follow.
True by default: Do you want to change the parameter XXX during the PEM run as prescribed in ob_ex_lsp.asc.
If false, it is set to constant (to the value of the tab_cntrl in the start)
RV

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

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/libf/phymars/comslope_mod.F90

    r2900 r2909  
    2323
    2424CONTAINS
    25   SUBROUTINE ini_comslope_h(ngrid)
     25  SUBROUTINE ini_comslope_h(ngrid,nslope_in)
    2626
    2727  IMPLICIT NONE
    2828  INTEGER, INTENT(IN) :: ngrid
     29  INTEGER, INTENT(IN) :: nslope_in
    2930
    30     allocate(def_slope(nslope+1))
    31     allocate(def_slope_mean(nslope))
    32     allocate(sky_slope(nslope))
    33     allocate(subslope_dist(ngrid,nslope))
     31    allocate(def_slope(nslope_in+1))
     32    allocate(def_slope_mean(nslope_in))
     33    allocate(sky_slope(nslope_in))
     34    allocate(subslope_dist(ngrid,nslope_in))
    3435    allocate(major_slope(ngrid))
    3536  END SUBROUTINE ini_comslope_h
     
    7071      qsurf_meshavg(:,:) = 0.
    7172
     73  if(nslope.eq.1) then
     74      albedo_meshavg(:,:) = albedo_slope(:,:,1)
     75      emis_meshavg(:) = emis_slope(:,1)
     76      tsurf_meshavg(:) = tsurf_slope(:,1)
     77      qsurf_meshavg(:,:) = qsurf_slope(:,:,1)
     78  else
     79
    7280  DO ig = 1,ngrid
    7381      DO islope = 1,nslope
     
    8492  ENDDO
    8593
     94  endif
     95
    8696  END SUBROUTINE compute_meshgridavg
    8797
  • trunk/LMDZ.MARS/libf/phymars/comsoil_h.F90

    r2900 r2909  
    7070  end subroutine end_comsoil_h
    7171
     72  subroutine ini_comsoil_h_slope_var(ngrid,nslope)
     73 
     74  implicit none
     75  integer,intent(in) :: ngrid ! number of atmospheric columns
     76  integer,intent(in) :: nslope ! number of sub grid slopes
     77 
     78    allocate(tsoil(ngrid,nsoilmx,nslope)) ! soil temperatures
     79    allocate(mthermdiff(ngrid,0:nsoilmx-1,nslope))
     80    allocate(thermdiff(ngrid,nsoilmx-1,nslope))
     81    allocate(coefd(ngrid,nsoilmx-1,nslope))
     82    allocate(alph(ngrid,nsoilmx-1,nslope))
     83    allocate(beta(ngrid,nsoilmx-1,nslope))
     84 
     85  end subroutine ini_comsoil_h_slope_var
     86
     87
     88  subroutine end_comsoil_h_slope_var
     89
     90  implicit none
     91
     92    if (allocated(tsoil)) deallocate(tsoil)
     93    if (allocated(mthermdiff)) deallocate(mthermdiff)
     94    if (allocated(thermdiff)) deallocate(thermdiff)
     95    if (allocated(coefd)) deallocate(coefd)
     96    if (allocated(alph)) deallocate(alph)
     97    if (allocated(beta)) deallocate(beta)
     98
     99  end subroutine end_comsoil_h_slope_var
     100
    72101end module comsoil_h
  • trunk/LMDZ.MARS/libf/phymars/dimradmars_mod.F90

    r2900 r2909  
    216216
    217217  end subroutine end_dimradmars_mod
     218
     219  subroutine ini_dimradmars_mod_slope_var(ngrid,nslope)
     220 
     221  implicit none
     222 
     223  integer,intent(in) :: ngrid ! number of atmospheric columns
     224  integer,intent(in) :: nslope ! number of subgrid scale slopes
     225
     226   allocate(albedo(ngrid,2,nslope))
     227   allocate(fluxrad_sky(ngrid,nslope))
     228   allocate(fluxrad(ngrid,nslope))
     229
     230  end subroutine ini_dimradmars_mod_slope_var
     231
     232  subroutine end_dimradmars_mod_slope_var
     233
     234  implicit none
     235
     236   if (allocated(albedo))      deallocate(albedo)
     237   if (allocated(fluxrad_sky)) deallocate(fluxrad_sky)
     238   if (allocated(fluxrad))     deallocate(fluxrad)
     239
     240  end subroutine end_dimradmars_mod_slope_var
    218241
    219242 
  • trunk/LMDZ.MARS/libf/phymars/phyetat0_mod.F90

    r2896 r2909  
    134134endif ! of if (startphy_file)
    135135
    136 if(nslope.ne.1) then
    137   call abort_physic(modname, &
    138                 "phyetat0: For now, nslope should be 1 (set in comslope_mod)",1)
    139 endif
    140 
    141 allocate(default_def_slope(nslope+1))
    142 !Sub-grid scale  subslopes
    143 if (nslope.eq.7) then
    144   default_def_slope(1) = -43.
    145   default_def_slope(2) = -19.
    146   default_def_slope(3) = -9.
    147   default_def_slope(4) = -3.
    148   default_def_slope(5) = 3.
    149   default_def_slope(6) = 9.
    150   default_def_slope(7) = 19.
    151   default_def_slope(8) = 43.
    152 elseif (nslope.eq.5) then
    153   default_def_slope(1) = -43.
    154   default_def_slope(2) = -9.
    155   default_def_slope(3) = -3.
    156   default_def_slope(4) = 3.
    157   default_def_slope(5) = 9.
    158   default_def_slope(6) = 43.
    159 elseif (nslope.eq.1) then
    160   default_def_slope(1) = 0.
    161   default_def_slope(2) = 0.
    162 endif
    163 
    164136if (startphy_file) then
    165137 call get_var("def_slope",def_slope,found)
     
    167139     startphy_slope=.false.
    168140     write(*,*)'slope_settings: Problem while reading <def_slope>'
    169      write(*,*)'default def_slope will be used'
    170      do islope=1,nslope+1
    171        def_slope(islope) = default_def_slope(islope)
    172      enddo
    173      write(*,*)'computing corresponding distribution <subslope_dist>'
    174      write(*,*)'For now, woth nslope=1, subslope_dist is straigforward'
    175      write(*,*)'Later this operation will be done by newstart with a specific routine'
    176      subslope_dist(:,:)=1.
    177      !call subslope_mola(ngrid,nslope,def_slope,subslope_dist)
     141     write(*,*)'Checking that nslope=1'
     142     if(nslope.eq.1) then
     143      write(*,*)'We take default def_slope and subslope dist'
     144      allocate(default_def_slope(nslope+1))
     145      default_def_slope(1) = 0.
     146      default_def_slope(2) = 0.
     147      subslope_dist(:,:)=1.
     148     else
     149       write(*,*)'Variable def_slope is not present in the start and'
     150       write(*,*)'you are trying to run with nslope!=1'
     151       write(*,*)'This is not possible, you have to go through newstart'
     152     endif
    178153   else
    179154     startphy_slope=.true.
     
    181156     if(.not.found) then
    182157       write(*,*)'slope_settings: Problem while reading <subslope_dist>'
    183        write(*,*)'computing a new distribution'
    184        write(*,*)'For now, woth nslope=1, subslope_dist is straigforward'
    185        write(*,*)'Later this operation will be done by newstart with a specific routine'
    186        subslope_dist(:,:)=1.
    187      !call subslope_mola(ngrid,nslope,def_slope,subslope_dist)
     158       write(*,*)'We have to abort.'
     159       write(*,*)'Please check that nslope is coherent, you can modify this parameter with newstart'
     160       call abort_physic(modname, &
     161                "phyetat0: Failed loading <subslope_dist>",1)
    188162     endif
    189163    endif
    190 else ! startphy_file 
    191  do islope=1,nslope+1
    192        def_slope(islope) = default_def_slope(islope)
    193      enddo
    194      write(*,*)'computing corresponding distribution <subslope_dist>'
    195      write(*,*)'For now, woth nslope=1, subslope_dist is straigforward'
    196      write(*,*)'Later this operation will be done by newstart with a specific routine'
    197      subslope_dist(:,:)=1.
    198      !call subslope_mola(ngrid,nslope,def_slope,subslope_dist)
     164else ! (no startphy_file)
     165     if(nslope.eq.1) then
     166      allocate(default_def_slope(nslope+1))
     167      default_def_slope(1) = 0.
     168      default_def_slope(2) = 0.
     169      subslope_dist(:,:)=1.
     170    else
     171     write(*,*)'Without starfi, lets first run with nslope=1'
     172       call abort_physic(modname, &
     173                "phyetat0: No startfi and nslope!=1",1)
     174   endif
    199175endif   
    200176
  • trunk/LMDZ.MARS/libf/phymars/phys_state_var_init_mod.F90

    r2900 r2909  
    6161      use dust_rad_adjust_mod, only: ini_dust_rad_adjust_mod, &
    6262                                     end_dust_rad_adjust_mod
    63       use comslope_mod, ONLY: nslope
     63      use comslope_mod, ONLY: nslope,end_comslope_h,ini_comslope_h
    6464      use netcdf
    6565      USE mod_phys_lmdz_para, ONLY: is_master,bcast
     
    182182      call ini_dust_rad_adjust_mod(ngrid)
    183183
     184      ! allocate arrays in "comslope_mod"
     185      call end_comslope_h
     186      call ini_comslope_h(ngrid,nslope)
     187
    184188      END SUBROUTINE phys_state_var_init
    185189
  • trunk/LMDZ.MARS/libf/phymars/physiq_mod.F

    r2908 r2909  
    577577c        ~~~~~~~~~~~~
    578578#ifndef MESOSCALE
    579 
    580          call ini_comslope_h(ngrid)
    581579
    582580! GCM. Read netcdf initial physical parameters.
  • trunk/LMDZ.MARS/libf/phymars/surfdat_h.F90

    r2900 r2909  
    116116  end subroutine end_surfdat_h
    117117
     118  subroutine ini_surfdat_h_slope_var(ngrid,nq,nslope)
     119 
     120  implicit none
     121  integer,intent(in) :: ngrid ! number of atmospheric columns
     122  integer,intent(in) :: nq ! number of tracers 
     123  integer,intent(in) :: nslope ! number of sub-grid scale slope 
     124    allocate(qsurf(ngrid,nq,nslope))
     125    allocate(tsurf(ngrid,nslope))
     126    allocate(watercap(ngrid,nslope))
     127    allocate(emis(ngrid,nslope))
     128    allocate(capcal(ngrid,nslope))
     129    allocate(fluxgrd(ngrid,nslope))
     130       
     131  end subroutine ini_surfdat_h_slope_var
     132
     133  subroutine end_surfdat_h_slope_var
     134
     135  implicit none
     136
     137    if (allocated(qsurf))         deallocate(qsurf)
     138    if (allocated(tsurf))         deallocate(tsurf)
     139    if (allocated(watercap))      deallocate(watercap)
     140    if (allocated(emis))          deallocate(emis)
     141    if (allocated(capcal))        deallocate(capcal)
     142    if (allocated(fluxgrd))       deallocate(fluxgrd)
     143   
     144  end subroutine end_surfdat_h_slope_var
     145
    118146end module surfdat_h
Note: See TracChangeset for help on using the changeset viewer.