/[lmdze]/trunk/dyn3d/disvert.f
ViewVC logotype

Contents of /trunk/dyn3d/disvert.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 90 - (show annotations)
Wed Mar 12 21:16:36 2014 UTC (10 years, 2 months ago) by guez
File size: 3890 byte(s)
Removed procedures ini_histday, ini_histhf, write_histday and
write_histhf.

Divided file regr_pr_coefoz.f into regr_pr_av.f and
regr_pr_int.f. (Following LMDZ.) Divided module regr_pr_coefoz into
modules regr_pr_av_m and regr_pr_int_m. Renamed regr_pr_av_coefoz to
regr_pr_av and regr_pr_int_coefoz to regr_pr_int. The idea is that
those procedures are more general than Mobidic.

Removed argument dudyn of calfis and physiq. dudyn is not used either
in LMDZ. Removed computation in calfis of unused variable zpsrf (not
used either in LMDZ). Removed useless computation of dqfi in calfis
(part 62): the results were overwritten. (Same in LMDZ.)

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

  ViewVC Help
Powered by ViewVC 1.1.21