/[lmdze]/trunk/libf/dyn3d/comvert.f90
ViewVC logotype

Contents of /trunk/libf/dyn3d/comvert.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 47 - (show annotations)
Fri Jul 1 15:00:48 2011 UTC (12 years, 10 months ago) by guez
File size: 4163 byte(s)
Split "thermcell.f" and "cv3_routines.f".
Removed copies of files that are now in "L_util".
Moved "mva9" and "diagetpq" to their own files.
Unified variable names across procedures.

1 module comvert
2
3 use dimens_m, only: llm
4
5 implicit none
6
7 private llm
8
9 real ap(llm+1), pa ! in Pa
10 real bp(llm+1)
11 real presnivs(llm) ! pressions approximatives des milieux de couches, en Pa
12 real, parameter:: preff = 101325. ! in Pa
13 real nivsigs(llm), nivsig(llm+1)
14
15 save
16
17 contains
18
19 SUBROUTINE disvert
20
21 ! From dyn3d/disvert.F, v 1.1.1.1 2004/05/19 12:53:05
22 ! Auteur : P. Le Van
23
24 ! This procedure sets the vertical grid.
25 ! It defines the host variables "ap", "bp", "presnivs", "nivsigs"
26 ! and "nivsig".
27 ! "pa" should be defined before this procedure is called.
28
29 use nr_util, only: pi
30 use l_util, only: new_unit
31
32 REAL s(llm+1)
33 ! "s(l)" is the atmospheric hybrid sigma-pressure coordinate at
34 ! the interface between layers "l" and "l-1"
35
36 real ds(llm)
37 ! "ds(l)" : épaisseur de la couche "l" dans la coordonnée "s"
38
39 INTEGER l, unit
40 REAL alpha, x(llm), trash
41
42 character(len=7):: s_sampling = "LMD5"
43 ! (other allowed values are "param", "strato1", "strato2" and "read")
44
45 real:: h = 7. ! scale height, in km
46 ! (used only if "s_sampling" == "param" or "strato1")
47
48 ! These variables are used only in the case 's_sampling == "param"':
49 real:: deltaz = 0.04 ! épaisseur de la première couche
50 real:: beta = 1.3 ! facteur d'accroissement en haut
51 real:: k0 = 20. ! nombre de couches dans la transition surface
52 real:: k1 = 1.2 ! nombre de couches dans la transition haute
53
54 REAL ZZ(llm + 1), DZ(llm) ! in km
55
56 namelist /disvert_nml/h, deltaz, beta, k0, k1, s_sampling
57
58 !-----------------------------------------------------------------------
59
60 print *, "Call sequence information: disvert"
61
62 forall (l = 1: llm) nivsigs(l) = REAL(l)
63 forall (l = 1: llm + 1) nivsig(l) = REAL(l)
64
65 ! Compute "s":
66
67 print *, "Enter namelist 'disvert_nml'."
68 read(unit=*, nml=disvert_nml)
69 write(unit=*, nml=disvert_nml)
70
71 select case (s_sampling)
72 case ("param")
73 s(1) = 1.
74 s(llm+1) = 0.
75 alpha = deltaz / tanh(1./k0) * 2.
76 forall (l = 2: llm) s(l) &
77 = cosh((l - 1) / k0) **(- alpha * k0 / h) &
78 * exp(- alpha / h * tanh((llm - k1) / k0) &
79 * beta **(l - 1 - (llm - k1)) / log(beta))
80 case ("LMD5")
81 ! Ancienne discrétisation
82 s(1) = 1.
83 s(llm+1) = 0.
84 forall (l = 1: llm) ds(l) &
85 = 1. + 7. * SIN(pi * (REAL(l)-0.5) / real(llm+1))**2
86 ds = ds / sum(ds)
87
88 DO l = llm, 2, -1
89 s(l) = s(l+1) + ds(l)
90 ENDDO
91 case ("strato1")
92 ! F. Lott 70 niveaux et plus
93 s(1) = 1.
94 s(llm+1) = 0.
95 forall (l = 1: llm) dz(l) = 1.56 + TANH(REAL(l - 12) / 5.) &
96 + TANH(REAL(l - llm) / 10.) / 2.
97
98 zz(1) = 0.
99 DO l = 2, llm + 1
100 zz(l) = zz(l - 1) + dz(l - 1)
101 end DO
102
103 s(2:llm) = (exp(- zz(2:llm) / h) - exp(- zz(llm + 1) / h)) &
104 / (1. - exp(- zz(llm + 1) / h))
105 case ("strato2")
106 ! Recommended by F. Lott for a domain including the stratosphere
107 s(1) = 1.
108 s(llm+1) = 0.
109 forall (l = 1: llm) x(l) = pi * (l - 0.5) / (llm + 1)
110
111 ds = (1. + 7. * SIN(x)**2) * (1. - tanh(2 * x / pi - 1.))**2 / 4.
112 ds = ds / sum(ds)
113
114 DO l = llm, 2, -1
115 s(l) = s(l+1) + ds(l)
116 ENDDO
117 case("read")
118 call new_unit(unit)
119 open(unit, file="hybrid.csv", status="old", action="read", &
120 position="rewind")
121 read(unit, fmt=*) ! skip title line
122 do l = 1, llm + 1
123 read(unit, fmt=*) trash, s(l)
124 end do
125 close(unit)
126 ! Quick check:
127 if (s(1) /= 1. .or. s(llm + 1) /= 0. .or. s(1) <= s(2)) then
128 print *, '"s" should be in descending order, from 1 to 0.'
129 stop 1
130 end if
131 case default
132 print *, 'Wrong value for "s_sampling"'
133 stop 1
134 END select
135
136 ! Calcul de "ap" et "bp" :
137 bp(:llm) = EXP(1. - 1. / s(:llm)**2)
138 bp(llm + 1) = 0.
139 ap = pa * (s - bp)
140
141 forall (l = 1: llm) presnivs(l) = 0.5 &
142 * (ap(l) + bp(l) * preff + ap(l+1) + bp(l+1) * preff)
143
144 END SUBROUTINE disvert
145
146 end module comvert

  ViewVC Help
Powered by ViewVC 1.1.21