Changeset 1538 for trunk/LMDZ.GENERIC


Ignore:
Timestamp:
Apr 8, 2016, 4:57:00 PM (9 years ago)
Author:
emillour
Message:

Generic GCM:

  • Made nsoilmx be no longer a "parameter" and thus added the possibility to define the number of subsurface layers nsoilmx, along with first layer thickness "lay1_soil" and companion coefficient "alpha_soil", in callphys.def at run time. As before (when these were hard-coded), these are such that the depth of soil mid-layers are: mlayer(k)=lay1_soil*alpha_soil(k-1/2), for k=0,nsoil-1

EM

Location:
trunk/LMDZ.GENERIC
Files:
6 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/README

    r1537 r1538  
    11771177- Initialize runoff in hydrol.F90.
    11781178
    1179 
    1180 
     1179== 08/04/2016 == EM
     1180- Made nsoilmx be no longer a "parameter" and thus added the possibility to define the
     1181  number of subsurface layers nsoilmx, along with first layer thickness "lay1_soil" and
     1182  companion coefficient "alpha_soil", in callphys.def at run time.
     1183  As before (when these were hard-coded), these are such that the depth of
     1184  soil mid-layers are:
     1185  mlayer(k)=lay1_soil*alpha_soil**(k-1/2), for k=0,nsoil-1
     1186
     1187
     1188
  • trunk/LMDZ.GENERIC/libf/dynlonlat_phylonlat/phystd/newstart.F

    r1524 r1538  
    3636      USE ener_mod, ONLY: etot0,ptot0,ztot0,stot0,ang0
    3737      use inifis_mod, only: inifis
     38      use control_mod, only: nday
    3839      implicit none
    3940
     
    9596c------------------
    9697      REAL tsurf(ngridmx)       ! surface temperature
    97       REAL tsoil(ngridmx,nsoilmx) ! soil temperature
     98      REAL,ALLOCATABLE :: tsoil(:,:) ! soil temperature
    9899!      REAL co2ice(ngridmx)     ! CO2 ice layer
    99100      REAL emis(ngridmx)        ! surface emissivity
     
    103104!      REAL rnaturfi(ngridmx)
    104105      real alb(iip1,jjp1),albfi(ngridmx) ! albedos
    105       real ith(iip1,jjp1,nsoilmx),ithfi(ngridmx,nsoilmx) ! thermal inertia (3D)
     106      real,allocatable :: ith(:,:,:),ithfi(:,:) ! thermal inertia (3D)
    106107      real surfith(iip1,jjp1),surfithfi(ngridmx) ! surface thermal inertia (2D)
    107108      REAL latfi(ngridmx),lonfi(ngridmx),airefi(ngridmx)
     
    113114
    114115      REAL rnat(ngridmx)
    115       REAL tslab(ngridmx,nsoilmx) ! slab ocean temperature
     116      REAL,ALLOCATABLE :: tslab(:,:) ! slab ocean temperature
    116117      REAL pctsrf_sic(ngridmx) ! sea ice cover
    117118      REAL tsea_ice(ngridmx) ! temperature sea_ice
     
    274275! ... this has to be done before phyetat0
    275276      call initracer(ngridmx,nqtot,tname)
     277
     278      ! default nsoilmx is set in comsoil_h
     279      call getin_p("nsoilmx",nsoilmx)
     280
     281      ! allocate some local arrays now that we know all dimensions
     282      allocate(tsoil(ngridmx,nsoilmx))
     283      allocate(ith(iip1,jjp1,nsoilmx))
     284      allocate(ithfi(ngridmx,nsoilmx))
     285      allocate(tslab(ngridmx,nsoilmx))
    276286
    277287! Take care of arrays in common modules
     
    475485!    (for instance initracer needs to know about some flags, and/or
    476486!      'datafile' path may be changed by user)
    477       call inifis(ngridmx,llm,nqtot,day_ini,daysec,dtphys,
     487      call inifis(ngridmx,llm,nqtot,day_ini,daysec,nday,dtphys,
    478488     &                latfi,lonfi,airefi,rad,g,r,cpp)
    479489
  • trunk/LMDZ.GENERIC/libf/dynlonlat_phylonlat/phystd/start2archive.F

    r1478 r1538  
    1919c=======================================================================
    2020
     21      use ioipsl_getin_p_mod, only: getin_p
    2122      use infotrac, only: infotrac_init, nqtot, tname
    2223      USE comsoil_h
     
    6970c ------------------------------------
    7071      REAL tsurf(ngridmx)       ! Surface temperature
    71       REAL tsoil(ngridmx,nsoilmx) ! Soil temperature
     72      REAL,ALLOCATABLE :: tsoil(:,:) ! Soil temperature
    7273      REAL co2ice(ngridmx)      ! CO2 ice layer
    7374      REAL q2(ngridmx,llm+1)
     
    9394      REAL T(ip1jmp1,llm),us(ip1jmp1,llm),vs(ip1jmp1,llm)
    9495      REAL tsurfS(ip1jmp1)
    95       REAL tsoilS(ip1jmp1,nsoilmx)
    96       REAL ithS(ip1jmp1,nsoilmx) ! Soil Thermal Inertia
     96      REAL,ALLOCATABLE :: tsoilS(:,:)
     97      REAL,ALLOCATABLE :: ithS(:,:) ! Soil Thermal Inertia
    9798      REAL co2iceS(ip1jmp1)
    9899      REAL q2S(ip1jmp1,llm+1)
     
    147148      CALL init_phys_lmdz(iim,jjp1,llm,1,(/(jjm-1)*iim+2/))
    148149      call initcomgeomphy
     150
     151! get value of nsoilmx (default is set in comsoil_h)
     152      call getin_p("nsoilmx",nsoilmx)
     153     
     154      ! allocate some local arrays now that we know all dimensions
     155      allocate(tsoil(ngridmx,nsoilmx))
     156      allocate(tsoilS(ip1jmp1,nsoilmx))
     157      allocate(ithS(ip1jmp1,nsoilmx))
    149158
    150159      ! ALLOCATE ARRAYS IN comgeomfi_h (usually done in inifis)
  • trunk/LMDZ.GENERIC/libf/phystd/comsoil_h.F90

    r1516 r1538  
    55!integer, parameter :: nsoilmx = 18 ! for z1=0.0002 m, depth = 18 m => mars case
    66!integer, parameter :: nsoilmx = 13 ! for z1=0.03 m, depth = 104.8 m => earth case
    7   integer, parameter :: nsoilmx = 18
    8 ! Full soil layer depths are set as: layer(k)=lay1*alpha**(k-1) , k=1,nsoil
    9 ! Mid soil layer depths are set as: mlayer(k)=lay1*alpha**(k-1/2) , k=0,nsoil-1
    10   real,parameter :: lay1=2.e-4 ! depth of first full soil layer
    11   real,parameter :: alpha=2
     7  integer,save :: nsoilmx = 18 ! default, but may be set in callphys.def
     8! Full soil layer depths are set as: layer(k)=lay1_soil*alpha_soil**(k-1) , k=1,nsoil
     9! Mid soil layer depths are set as: mlayer(k)=lay1_soil*alpha_soil**(k-1/2) , k=0,nsoil-1
     10  real,save :: lay1_soil=2.e-4 ! depth (m) of first full soil layer (may be set in callphys.def)
     11  real,save :: alpha_soil=2 ! coefficient for soil layer thickness (may be set in callphys.def)
     12!$OMP THREADPRIVATE(nsoilmx,lay1_soil,alpha_soil)
    1213
    1314  real,save,allocatable,dimension(:) :: layer      ! soil layer depths
  • trunk/LMDZ.GENERIC/libf/phystd/inifis_mod.F90

    r1529 r1538  
    1414  use comdiurn_h, only: sinlat, coslat, sinlon, coslon
    1515  use comgeomfi_h, only: long, lati, area, totarea, totarea_planet
    16   use comsoil_h, only: ini_comsoil_h
     16  use comsoil_h, only: ini_comsoil_h, nsoilmx, lay1_soil, alpha_soil
    1717  use time_phylmdz_mod, only: ecritphy,day_step,iphysiq, &
    1818                              init_time, daysec, dtphys
     
    322322     call getin_p("graybody",graybody)
    323323     write(*,*)" graybody = ",graybody
     324
     325! Soil model
     326     write(*,*)"Number of sub-surface layers for soil scheme?"
     327     ! default value of nsoilmx set in comsoil_h
     328     call getin_p("nsoilmx",nsoilmx)
     329     write(*,*)" nsoilmx=",nsoilmx
     330     
     331     write(*,*)"Thickness of topmost soil layer (m)?"
     332     ! default value of lay1_soil set in comsoil_h
     333     call getin_p("lay1_soil",lay1_soil)
     334     write(*,*)" lay1_soil = ",lay1_soil
     335     
     336     write(*,*)"Coefficient for soil layer thickness distribution?"
     337     ! default value of alpha_soil set in comsoil_h
     338     call getin_p("alpha_soil",alpha_soil)
     339     write(*,*)" alpha_soil = ",alpha_soil
    324340
    325341! Slab Ocean
  • trunk/LMDZ.GENERIC/libf/phystd/soil_settings.F

    r1516 r1538  
    33!      use netcdf
    44      use comsoil_h, only: layer, mlayer, inertiedat, volcapa,
    5      &                     lay1, alpha
     5     &                     lay1_soil, alpha_soil
    66      use iostart, only: inquire_field_ndims, get_var, get_field,
    77     &                   inquire_field, inquire_dimension_length
     
    7171      real,dimension(:),allocatable :: newval
    7272
    73       real malpha,mlay1 ! coefficients for building layers
     73      real malpha,mlay1_soil ! coefficients for building layers
    7474      real xmin,xmax ! to display min and max of a field
    7575
     
    8484
    8585        if (dimlen.ne.nsoil) then
    86           write(*,*)'soil_settings: Interpolation of soil temperature ',
    87      &              'and thermal inertia will be required!'
    8886        ! if dimlen doesn't match nsoil, then interpolation of
    8987        ! soil temperatures and thermal inertia will be requiered
    9088          interpol=.true.
    91         ! allocate olmlayer 
    92           allocate(oldmlayer(dimlen),stat=ierr)
    93           if (ierr.ne.0) then
    94             write(*,*) 'soil_settings: failed allocation of oldmlayer!'
    95             stop
    96           endif
    9789        endif
     90
     91        ! allocate oldmlayer 
     92        allocate(oldmlayer(dimlen),stat=ierr)
     93        if (ierr.ne.0) then
     94          write(*,*) 'soil_settings: failed allocation of oldmlayer!'
     95          stop
     96        endif
     97
     98        ! check if olmlayer distribution matches current one
     99        call get_var("soildepth",oldmlayer,found)
     100        if (found) then
     101          malpha=oldmlayer(2)/oldmlayer(1)
     102          if ((abs(malpha-alpha_soil)/alpha_soil).gt.1.e-6) then
     103            ! alpha values are too different, intepolation needed
     104            interpol=.true.
     105          endif
     106          ! check if loaded mid-layer depth value differs from current one
     107          if (abs((oldmlayer(1)-lay1_soil*alpha_soil**(-1./2.))/
     108     &             (lay1_soil*alpha_soil**(-1./2.))).gt.1.e-6) then
     109            interpol=.true.
     110          endif
     111        endif
     112
     113        if (interpol) then
     114          write(*,*)'soil_settings: Interpolation of soil temperature ',
     115     &              'and thermal inertia will be required!'
     116        endif
     117       
    98118! 1.2 Find out the # of dimensions <inertiedat> was defined as using
    99119      ndims=inquire_field_ndims("inertiedat")
     
    133153      if (interpol) then
    134154      ! default mlayer distribution, following a power law:
    135       !  mlayer(k)=lay1*alpha**(k-1/2)
     155      !  mlayer(k)=lay1_soil*alpha_soil**(k-1/2)
    136156        do iloop=0,nsoil-1
    137           mlayer(iloop)=lay1*(alpha**(iloop-0.5))
     157          mlayer(iloop)=lay1_soil*(alpha_soil**(iloop-0.5))
    138158        enddo
    139159      endif
    140160! 1.5 Build layer(); following the same law as mlayer()
    141161      ! Assuming layer distribution follows mid-layer law:
    142       ! layer(k)=lay1*alpha**(k-1)
    143       mlay1=sqrt(mlayer(0)*mlayer(1))
     162      ! layer(k)=lay1_soil*alpha_soil**(k-1)
     163      mlay1_soil=sqrt(mlayer(0)*mlayer(1))
    144164      malpha=mlayer(1)/mlayer(0)
    145165      ! Check that these values are the same as those prescibed for mlayers
    146       if ((abs(mlay1-lay1)/lay1).gt.1.e-6) then
    147          write(*,*) "soil_settings error: mlay1=",mlay1
    148          write(*,*) " does not match comsoil_h lay1=",lay1
     166      if ((abs(mlay1_soil-lay1_soil)/lay1_soil).gt.1.e-6) then
     167         write(*,*) "soil_settings error: mlay1_soil=",mlay1_soil
     168         write(*,*) " does not match comsoil_h lay1_soil=",lay1_soil
    149169         stop
    150170      endif
    151       if ((abs(malpha-alpha)/alpha).gt.1.e-6) then
     171      if ((abs(malpha-alpha_soil)/alpha_soil).gt.1.e-6) then
    152172         write(*,*) "soil_settings error: malpha=",malpha
    153          write(*,*) " does not match comsoil_h alpha=",alpha
     173         write(*,*) " does not match comsoil_h alpha_soil=",alpha_soil
    154174         stop
    155175      endif
    156176      do iloop=1,nsoil
    157         layer(iloop)=mlay1*(malpha**(iloop-1))
     177        layer(iloop)=mlay1_soil*(malpha**(iloop-1))
    158178      enddo
    159179
Note: See TracChangeset for help on using the changeset viewer.