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

Annotation of /trunk/dyn3d/disvert.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 90 - (hide 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 guez 66 module disvert_m
2 guez 3
3     use dimens_m, only: llm
4    
5     implicit none
6    
7     private llm
8    
9 guez 83 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 guez 3 real, parameter:: preff = 101325. ! in Pa
16    
17     contains
18    
19     SUBROUTINE disvert
20    
21 guez 90 ! From dyn3d/disvert.F, version 1.1.1.1, 2004/05/19 12:53:05
22 guez 66 ! Author: P. Le Van
23 guez 3
24 guez 66 ! This procedure sets the vertical grid. It defines the host
25 guez 83 ! variables "ap", "bp", "presnivs". "pa" should be defined before
26     ! this procedure is called.
27 guez 3
28 guez 66 use jumble, only: new_unit
29 guez 78 use nr_util, only: pi, assert
30 guez 57 use unit_nml_m, only: unit_nml
31 guez 3
32 guez 66 ! Local:
33    
34 guez 3 REAL s(llm+1)
35 guez 22 ! "s(l)" is the atmospheric hybrid sigma-pressure coordinate at
36     ! the interface between layers "l" and "l-1"
37 guez 3
38     real ds(llm)
39 guez 90 ! "ds(l)" : \'epaisseur de la couche "l" dans la coordonn\'ee "s"
40 guez 3
41 guez 42 INTEGER l, unit
42 guez 53 REAL alpha, x(llm)
43 guez 3
44 guez 78 character(len=7):: vert_sampling = "tropo"
45     ! other allowed values are "param", "strato" and "read"
46 guez 3
47     real:: h = 7. ! scale height, in km
48 guez 78 ! used only if vert_sampling == "param"
49 guez 3
50 guez 78 ! These variables are used only in the case vert_sampling == "param":
51 guez 90 real:: deltaz = 0.04 ! \'epaisseur de la premi\`ere couche
52 guez 3 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 guez 78 namelist /disvert_nml/h, deltaz, beta, k0, k1, vert_sampling
57 guez 3
58     !-----------------------------------------------------------------------
59    
60     print *, "Call sequence information: disvert"
61    
62     print *, "Enter namelist 'disvert_nml'."
63     read(unit=*, nml=disvert_nml)
64 guez 57 write(unit_nml, nml=disvert_nml)
65 guez 3
66 guez 78 select case (vert_sampling)
67 guez 3 case ("param")
68 guez 42 s(1) = 1.
69     s(llm+1) = 0.
70 guez 3 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 guez 78 call compute_ab
76 guez 65 case ("tropo")
77 guez 42 s(1) = 1.
78     s(llm+1) = 0.
79 guez 3 forall (l = 1: llm) ds(l) &
80 guez 66 = 1. + 7. * SIN(pi * (REAL(l) - 0.5) / real(llm + 1))**2
81 guez 3 ds = ds / sum(ds)
82    
83     DO l = llm, 2, -1
84     s(l) = s(l+1) + ds(l)
85     ENDDO
86    
87 guez 78 call compute_ab
88 guez 65 case ("strato")
89 guez 3 ! Recommended by F. Lott for a domain including the stratosphere
90 guez 42 s(1) = 1.
91     s(llm+1) = 0.
92 guez 3 forall (l = 1: llm) x(l) = pi * (l - 0.5) / (llm + 1)
93    
94 guez 66 ds = (0.3 + 7. * SIN(x)**2) * (1. - tanh(2 * x / pi - 1.))**2 / 4.
95 guez 3 ds = ds / sum(ds)
96    
97     DO l = llm, 2, -1
98     s(l) = s(l+1) + ds(l)
99     ENDDO
100 guez 78
101     call compute_ab
102 guez 42 case("read")
103 guez 78 ! 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 guez 42 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 guez 78 read(unit, fmt=*) ap(l), bp(l)
112 guez 42 end do
113     close(unit)
114     ! Quick check:
115 guez 78 call assert(ap(1) == 0., ap(llm + 1) == 0., bp(1) == 1., &
116     bp(llm + 1) == 0., "disvert: bad ap or bp values")
117 guez 3 case default
118 guez 78 print *, 'Wrong value for "vert_sampling"'
119 guez 3 stop 1
120     END select
121    
122 guez 32 forall (l = 1: llm) presnivs(l) = 0.5 &
123     * (ap(l) + bp(l) * preff + ap(l+1) + bp(l+1) * preff)
124 guez 3
125 guez 78 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 guez 3 END SUBROUTINE disvert
137    
138 guez 66 end module disvert_m

  ViewVC Help
Powered by ViewVC 1.1.21