Ignore:
Timestamp:
Apr 4, 2012, 12:14:59 PM (13 years ago)
Author:
lguez
Message:

In procedure "disvert", the values of "sig" are now chosen according
to the value of "vert_sampling". "ok_strato" no longer appears in
"disvert". The analytical formula that was used when "ok_strato" was
true in "disvert" is now selected when "vert_sampling" is set to
"strato". Motivation: we can now read values of "sig" instead of
using this analytical formula, even when "ok_strato" is true.

Default value of "vert_sampling" with llm=39 gives the vertical
distribution of CMIP5 simulations.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • LMDZ5/trunk/libf/dyn3dpar/disvert.F90

    r1524 r1620  
    44
    55  ! Auteur : P. Le Van
     6
     7  use new_unit_m, only: new_unit
     8  use ioipsl, only: getin
    69
    710  IMPLICIT NONE
     
    2629  real zk, zkm1, dzk1, dzk2, k0, k1
    2730
    28   INTEGER l
     31  INTEGER l, unit
    2932  REAL dsigmin
    3033  REAL alpha, beta, deltaz
    31   INTEGER iostat
    3234  REAL x
    3335  character(len=*),parameter :: modname="disvert"
    3436
     37  character(len=6):: vert_sampling
     38  ! (allowed values are "param", "tropo", "strato" and "read")
     39
    3540  !-----------------------------------------------------------------------
     41
     42  print *, "Call sequence information: disvert"
    3643
    3744  ! default scaleheight is 8km for earth
    3845  scaleheight=8.
    3946
    40   OPEN(99, file='sigma.def', status='old', form='formatted', iostat=iostat)
     47  vert_sampling = "strato" ! default value
     48  call getin('vert_sampling', vert_sampling)
     49  print *, 'vert_sampling = ' // vert_sampling
    4150
    42   IF (iostat == 0) THEN
    43      ! cas 1 on lit les options dans sigma.def:
     51  select case (vert_sampling)
     52  case ("param")
     53     ! On lit les options dans sigma.def:
     54     OPEN(99, file='sigma.def', status='old', form='formatted')
    4455     READ(99, *) scaleheight ! hauteur d'echelle 8.
    4556     READ(99, *) deltaz ! epaiseur de la premiere couche 0.04
     
    6879
    6980     sig(llm+1)=0.
    70 
    71      DO l = 1, llm
    72         dsig(l) = sig(l)-sig(l+1)
    73      end DO
    74   ELSE
    75      if (ok_strato) then
    76         if (llm==39) then
    77            dsigmin=0.3
    78         else if (llm==50) then
    79            dsigmin=1.
    80         else
    81            write(lunout,*) trim(modname), &
    82            ' ATTENTION discretisation z a ajuster'
    83            dsigmin=1.
    84         endif
    85         write(lunout,*) trim(modname), &
    86         ' Discretisation verticale DSIGMIN=',dsigmin
    87      endif
    88 
     81  case("tropo")
    8982     DO l = 1, llm
    9083        x = 2*asin(1.) * (l - 0.5) / (llm + 1)
    91 
    92         IF (ok_strato) THEN
    93            dsig(l) =(dsigmin + 7. * SIN(x)**2) &
    94                 *(0.5*(1.-tanh(1.*(x-asin(1.))/asin(1.))))**2
    95         ELSE
    96            dsig(l) = 1.0 + 7.0 * SIN(x)**2
    97         ENDIF
     84        dsig(l) = 1.0 + 7.0 * SIN(x)**2
    9885     ENDDO
    9986     dsig = dsig / sum(dsig)
     
    10289        sig(l) = sig(l+1) + dsig(l)
    10390     ENDDO
    104   ENDIF
     91  case("strato")
     92     if (llm==39) then
     93        dsigmin=0.3
     94     else if (llm==50) then
     95        dsigmin=1.
     96     else
     97        write(lunout,*) trim(modname), ' ATTENTION discretisation z a ajuster'
     98        dsigmin=1.
     99     endif
     100     WRITE(LUNOUT,*) trim(modname), 'Discretisation verticale DSIGMIN=',dsigmin
     101
     102     DO l = 1, llm
     103        x = 2*asin(1.) * (l - 0.5) / (llm + 1)
     104        dsig(l) =(dsigmin + 7. * SIN(x)**2) &
     105             *(0.5*(1.-tanh(1.*(x-asin(1.))/asin(1.))))**2
     106     ENDDO
     107     dsig = dsig / sum(dsig)
     108     sig(llm+1) = 0.
     109     DO l = llm, 1, -1
     110        sig(l) = sig(l+1) + dsig(l)
     111     ENDDO
     112  case("read")
     113     call new_unit(unit)
     114     open(unit, file="hybrid.txt", status="old", action="read", &
     115          position="rewind")
     116     read(unit, fmt=*) ! skip title line
     117     do l = 1, llm + 1
     118        read(unit, fmt=*) sig(l)
     119     end do
     120     close(unit)
     121  case default
     122     call abort_gcm("disvert", 'Wrong value for "vert_sampling"', 1)
     123  END select
    105124
    106125  DO l=1, llm
Note: See TracChangeset for help on using the changeset viewer.