/[lmdze]/trunk/phylmd/conf_phys.f90
ViewVC logotype

Contents of /trunk/phylmd/conf_phys.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 339 - (show annotations)
Thu Sep 26 17:08:42 2019 UTC (4 years, 7 months ago) by guez
File size: 1781 byte(s)
Simplify newmicro and rename dummy arguments

Rename dummy argument ucov of procedure advect to ang_3d. Rename dummy
argument radliq of procedure fisrtilp to cldliq. Rename dummy argument
qlwp of procedure newmicro to cldliq. Motivation: same name across
procedures.

Remove useless intermediary local variable rel in procedure
newmicro. The value of rad_chaud can be used instead of 1 in the
computation of cldtau if zflwp is 0: it does not matter.

1 module conf_phys_m
2
3 implicit none
4
5 integer, protected:: iflag_pbl = 1 ! for the planetary boundary layer
6 ! 6 : Mellor and Yamada 2.0
7 ! 8 : Mellor and Yamada 2.5
8
9 REAL, protected:: rad_chau1 = 13., rad_chau2 = 9.
10 real:: epmax = 0.993 ! \'efficacit\'e de pr\'ecipitation
11 integer:: iflag_clw = 0
12
13 contains
14
15 subroutine conf_phys
16
17 ! From phylmd/conf_phys.F90, version 1.7 2005/07/05 07:21:23
18
19 ! Configuration de la "physique" de LMDZ.
20
21 USE clesphys, ONLY: read_clesphys
22 use clesphys2, only: read_clesphys2
23 USE comfisrtilp, ONLY: cld_lc_con, cld_lc_lsc, cld_tau_con, &
24 cld_tau_lsc, coef_eva, ffallv_con, ffallv_lsc, iflag_pdf, reevap_ice
25 use nr_util, only: assert
26 use unit_nml_m, only: unit_nml
27 USE yomcst, ONLY: read_YOMCST
28
29 namelist /conf_phys_nml/ epmax, iflag_clw, cld_lc_lsc, cld_lc_con, &
30 cld_tau_lsc, cld_tau_con, ffallv_lsc, ffallv_con, coef_eva, &
31 reevap_ice, iflag_pdf, iflag_pbl
32
33 namelist /nuagecom/ rad_chau1, rad_chau2
34
35 !-----------------------------------------------------------
36
37 print *, "Call sequence information: conf_phys"
38 call read_clesphys2
39 call read_YOMCST
40
41 cld_lc_lsc = 2.6e-4
42 cld_lc_con = 2.6e-4
43 cld_tau_lsc = 3600.
44 cld_tau_con = 3600.
45 ffallv_lsc = 1.
46 ffallv_con = 1.
47 coef_eva = 2.e-5
48 reevap_ice = .false.
49 iflag_pdf = 0
50
51 print *, "Enter namelist 'conf_phys_nml'."
52 read(unit=*, nml=conf_phys_nml)
53 write(unit_nml, nml=conf_phys_nml)
54
55 call assert(any(iflag_pbl == [0, 1, 6, 8, 9]), &
56 "conf_phys: bad value for iflag_pbl")
57 call read_clesphys
58
59 print *, "Enter namelist 'nuagecom'."
60 read(unit=*, nml=nuagecom)
61 write(unit_nml, nml=nuagecom)
62
63 end subroutine conf_phys
64
65 end module conf_phys_m

  ViewVC Help
Powered by ViewVC 1.1.21