/[lmdze]/trunk/phylmd/clesphys2.f
ViewVC logotype

Annotation of /trunk/phylmd/clesphys2.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 267 - (hide annotations)
Thu May 3 16:14:08 2018 UTC (6 years ago) by guez
File size: 1224 byte(s)
Rename procedure clmain to pbl_surface (following LMDZ).

Remove choice soil_model = f. This choice made the algorithm unclear
in interfsurf_hq. Also soil_model = f is never used in LMDZ. radsol
was intent inout in clqh because of its possible modification in
interfsurf_hq, but the corresponding actual argument yrads is not used
in pbl_surface. The modification of radsol in interfsurf_hq was a bad
idea. Now we can more clearly make radsol an intent in argument of
clqh and interfsurf_hq.

1 guez 12 module clesphys2
2    
3 guez 51 ! From version 1.3 2005/06/06 13:16:33
4 guez 12
5     implicit none
6    
7 guez 51 LOGICAL:: new_oliq = .TRUE.
8     ! Permet de mettre en route la nouvelle parametrisation de l'eau liquide
9 guez 12
10 guez 69 ! Pour l'orographie:
11     LOGICAL:: ok_orodr = .TRUE.
12     LOGICAL:: ok_orolf = .TRUE.
13 guez 12
14 guez 51 LOGICAL:: ok_limitvrai = .FALSE.
15     ! On peut forcer le modele a lire le fichier SST de la bonne
16 guez 69 ! annee.
17 guez 12
18 guez 51 INTEGER:: nbapp_rad = 12
19 guez 69 ! nombre d'appels des routines de rayonnements par jour
20 guez 12
21 guez 182 logical:: conv_emanuel = .true. ! convection scheme of Emanuel, else Tiedtke
22 guez 12
23 guez 13 contains
24    
25     subroutine read_clesphys2
26    
27 guez 57 use unit_nml_m, only: unit_nml
28 guez 99 use nr_util, only: assert
29 guez 154 use conf_gcm_m, only: day_step, iphysiq
30 guez 57
31 guez 267 namelist /clesphys2_nml/new_oliq, ok_orodr, ok_orolf, ok_limitvrai, &
32     nbapp_rad, conv_emanuel
33 guez 13
34     !------------------------------------
35    
36     print *, "Enter namelist 'clesphys2_nml'."
37     read(unit=*, nml=clesphys2_nml)
38 guez 57 write(unit_nml, nml=clesphys2_nml)
39 guez 154 call assert(mod(day_step / iphysiq, nbapp_rad) == 0, &
40     "read_clesphys2 nbapp_rad")
41 guez 209 call assert(nbapp_rad >= 4, &
42     "read_clesphys2: minimum 4 calls to radiative transfer per day")
43 guez 13
44     end subroutine read_clesphys2
45    
46 guez 12 end module clesphys2

  ViewVC Help
Powered by ViewVC 1.1.21