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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 57 - (show annotations)
Mon Jan 30 12:54:02 2012 UTC (12 years, 3 months ago) by guez
File size: 4186 byte(s)
Write used namelists to file "" instead of standard output.

Avoid aliasing in "inidissip" in calls to "divgrad2", "divgrad",
"gradiv2", "gradiv", "nxgraro2" and "nxgrarot". Add a degenerate
dimension to arrays so they have rank 3, like the dummy arguments in
"divgrad2", "divgrad", "gradiv2", "gradiv", "nxgraro2" and "nxgrarot".

Extract the initialization part from "bilan_dyn" and make a separate
procedure, "init_dynzon", from it.

Move variables from modules "iniprint" and "logic" to module
"conf_gcm_m".

Promote internal procedures of "fxy" to private procedures of module
"fxy_m".

Extracted documentation from "inigeom". Removed useless "save"
attributes. Removed useless intermediate variables. Extracted
processing of poles from loop on latitudes. Write coordinates to file
"longitude_latitude.txt" instead of standard output.

Do not use ozone tracer for radiative transfer.

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

  ViewVC Help
Powered by ViewVC 1.1.21