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

Contents of /trunk/dyn3d/disvert.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 78 - (show annotations)
Wed Feb 5 17:51:07 2014 UTC (10 years, 3 months ago) by guez
Original Path: trunk/dyn3d/disvert.f90
File size: 4008 byte(s)
Moved procedure inigeom into module comgeom.

In disvert, renamed s_sampling to vert_sampling, following
LMDZ. Removed choice strato1. In case read, read ap and bp instead of
s (following LMDZ).

Added argument phis to start_init_orog and start_init_dyn, and removed
variable phis of module start_init_orog_m. In etat0 and
start_init_orog, renamed relief to zmea_2d. In start_init_dyn, renamed
psol to ps.

In start_init_orog, renamed relief_hi to relief. No need to set
phis(iim + 1, :) = phis(1, :), already done in grid_noro.

Documentation for massbar out of SVN, in massbar.txt. Documentation
was duplicated in massdair, but not relevant in massdair.

In conflx, no need to initialize pen_[ud] and pde_[ud]. In flxasc,
used intermediary variable fact (following LMDZ).

In grid_noro, added local variable zmea0 for zmea not smoothed and
computed zphi from zmea instead of zmea0 (following LMDZ). This
changes the results of ce0l.

Removed arguments pen_u and pde_d of phytrac and nflxtr, which were
not used.

1 module disvert_m
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 ! Author: P. Le Van
23
24 ! This procedure sets the vertical grid. It defines the host
25 ! variables "ap", "bp", "presnivs", "nivsigs" and "nivsig". "pa"
26 ! should be defined before 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)" : épaisseur de la couche "l" dans la coordonnée "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 ! épaisseur de la première 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 forall (l = 1: llm) nivsigs(l) = REAL(l)
63 forall (l = 1: llm + 1) nivsig(l) = REAL(l)
64
65 print *, "Enter namelist 'disvert_nml'."
66 read(unit=*, nml=disvert_nml)
67 write(unit_nml, nml=disvert_nml)
68
69 select case (vert_sampling)
70 case ("param")
71 s(1) = 1.
72 s(llm+1) = 0.
73 alpha = deltaz / tanh(1./k0) * 2.
74 forall (l = 2: llm) s(l) &
75 = cosh((l - 1) / k0) **(- alpha * k0 / h) &
76 * exp(- alpha / h * tanh((llm - k1) / k0) &
77 * beta **(l - 1 - (llm - k1)) / log(beta))
78 call compute_ab
79 case ("tropo")
80 s(1) = 1.
81 s(llm+1) = 0.
82 forall (l = 1: llm) ds(l) &
83 = 1. + 7. * SIN(pi * (REAL(l) - 0.5) / real(llm + 1))**2
84 ds = ds / sum(ds)
85
86 DO l = llm, 2, -1
87 s(l) = s(l+1) + ds(l)
88 ENDDO
89
90 call compute_ab
91 case ("strato")
92 ! Recommended by F. Lott for a domain including the stratosphere
93 s(1) = 1.
94 s(llm+1) = 0.
95 forall (l = 1: llm) x(l) = pi * (l - 0.5) / (llm + 1)
96
97 ds = (0.3 + 7. * SIN(x)**2) * (1. - tanh(2 * x / pi - 1.))**2 / 4.
98 ds = ds / sum(ds)
99
100 DO l = llm, 2, -1
101 s(l) = s(l+1) + ds(l)
102 ENDDO
103
104 call compute_ab
105 case("read")
106 ! Read "ap" and "bp". First line is skipped (title line). "ap"
107 ! should be in Pa. First couple of values should correspond to
108 ! the surface, that is : "bp" should be in descending order.
109 call new_unit(unit)
110 open(unit, file="hybrid.csv", status="old", action="read", &
111 position="rewind")
112 read(unit, fmt=*) ! skip title line
113 do l = 1, llm + 1
114 read(unit, fmt=*) ap(l), bp(l)
115 end do
116 close(unit)
117 ! Quick check:
118 call assert(ap(1) == 0., ap(llm + 1) == 0., bp(1) == 1., &
119 bp(llm + 1) == 0., "disvert: bad ap or bp values")
120 case default
121 print *, 'Wrong value for "vert_sampling"'
122 stop 1
123 END select
124
125 forall (l = 1: llm) presnivs(l) = 0.5 &
126 * (ap(l) + bp(l) * preff + ap(l+1) + bp(l+1) * preff)
127
128 contains
129
130 subroutine compute_ab
131
132 ! Calcul de "ap" et "bp" :
133 bp(:llm) = EXP(1. - 1. / s(:llm)**2)
134 bp(llm + 1) = 0.
135 ap = pa * (s - bp)
136
137 end subroutine compute_ab
138
139 END SUBROUTINE disvert
140
141 end module disvert_m

  ViewVC Help
Powered by ViewVC 1.1.21