1 | module lwi_mod |
---|
2 | |
---|
3 | implicit none |
---|
4 | |
---|
5 | contains |
---|
6 | |
---|
7 | subroutine lwi (ig0,kdlon,kflev |
---|
8 | . ,psi,zdblay,pdp |
---|
9 | . ,newpcolc ) |
---|
10 | |
---|
11 | |
---|
12 | use dimradmars_mod, only: ndlo2, ndlon, nflev, nir |
---|
13 | use yomlw_h, only: gcp, nlaylte, xi |
---|
14 | USE comcstfi_h, ONLY: g, cpp |
---|
15 | USE time_phylmdz_mod, ONLY: dtphys |
---|
16 | implicit none |
---|
17 | |
---|
18 | include "callkeys.h" |
---|
19 | |
---|
20 | CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC |
---|
21 | C |
---|
22 | C - lwi - |
---|
23 | C |
---|
24 | C PURPOSE: Shema semi - implicite |
---|
25 | C |
---|
26 | C |
---|
27 | CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC |
---|
28 | |
---|
29 | |
---|
30 | c************************************************************************ |
---|
31 | c |
---|
32 | c 0. Declarations |
---|
33 | c ------------ |
---|
34 | c |
---|
35 | c------------------------------------------------------------------------- |
---|
36 | c 0.1 Arguments |
---|
37 | c --------- |
---|
38 | c |
---|
39 | |
---|
40 | integer,intent(in) :: ig0 |
---|
41 | integer,intent(in) :: kdlon |
---|
42 | integer,intent(in) :: kflev |
---|
43 | |
---|
44 | real,intent(in) :: psi(ndlo2,kflev) |
---|
45 | real,intent(in) :: zdblay(ndlo2,nir,kflev) |
---|
46 | real,intent(in) :: pdp(ndlo2,kflev) |
---|
47 | |
---|
48 | real,intent(out) :: newpcolc(ndlo2,kflev) |
---|
49 | |
---|
50 | c------------------------------------------------------------------------- |
---|
51 | c 0.2 local arrays |
---|
52 | c ------------ |
---|
53 | c |
---|
54 | real di(ndlon,nflev) |
---|
55 | . , hi(ndlon,nflev) |
---|
56 | . , bi(ndlon,nflev) |
---|
57 | |
---|
58 | real ci(ndlon,nflev) |
---|
59 | . , ai(ndlon,nflev) |
---|
60 | real deltat |
---|
61 | |
---|
62 | real semit, denom |
---|
63 | |
---|
64 | integer i, jl |
---|
65 | |
---|
66 | c************************************************************************ |
---|
67 | c |
---|
68 | c 1. Initialisations |
---|
69 | c --------------- |
---|
70 | c |
---|
71 | c----------------------------------------------------------------------- |
---|
72 | |
---|
73 | deltat = dtphys * iradia |
---|
74 | c print*,'SEMI = ',semi, '(expl:0 semi-implicite:0.5 impl:1)' |
---|
75 | semit = semi * deltat |
---|
76 | c semi = 0. |
---|
77 | |
---|
78 | c print*,'dtphys,iradia,deltat,semit:',dtphys,iradia,deltat,semit |
---|
79 | c print*,'g,cpp',g,cpp |
---|
80 | |
---|
81 | |
---|
82 | c************************************************************************ |
---|
83 | c |
---|
84 | c 2. |
---|
85 | c --------------- |
---|
86 | c |
---|
87 | c------------------------------------------------------------------------- |
---|
88 | c 2.1 Calcul des di |
---|
89 | c ------------- |
---|
90 | c |
---|
91 | |
---|
92 | |
---|
93 | do i = 1 , nlaylte-1 |
---|
94 | do jl = 1 , kdlon |
---|
95 | c ------------------- |
---|
96 | di(jl,i) = 1 + semit * (g / pdp(jl,i) / cpp) * ( |
---|
97 | . ( xi(ig0+jl,1,i,nlaylte+1) |
---|
98 | . + xi(ig0+jl,1,i,i+1) |
---|
99 | . + xi(ig0+jl,1,i,i-1) ) |
---|
100 | . * zdblay(jl,1,i) |
---|
101 | . + ( xi(ig0+jl,2,i,nlaylte+1) |
---|
102 | . + xi(ig0+jl,2,i,i+1) |
---|
103 | . + xi(ig0+jl,2,i,i-1) ) |
---|
104 | . * zdblay(jl,2,i) |
---|
105 | . ) |
---|
106 | c ------------------- |
---|
107 | enddo |
---|
108 | enddo |
---|
109 | |
---|
110 | c couche nlaylte |
---|
111 | c ------------ |
---|
112 | c , on enleve i,i+1 sinon on a 2 fois le cooling2space |
---|
113 | |
---|
114 | do jl = 1 , kdlon |
---|
115 | c ------------------- |
---|
116 | di(jl,nlaylte) = 1 + semit * (g / pdp(jl,nlaylte) / cpp) * ( |
---|
117 | . ( xi(ig0+jl,1,nlaylte,nlaylte+1) |
---|
118 | . + xi(ig0+jl,1,nlaylte,nlaylte-1) ) |
---|
119 | . * zdblay(jl,1,nlaylte) |
---|
120 | . + ( xi(ig0+jl,2,nlaylte,nlaylte+1) |
---|
121 | . + xi(ig0+jl,2,nlaylte,nlaylte-1) ) |
---|
122 | . * zdblay(jl,2,nlaylte) |
---|
123 | . ) |
---|
124 | c ------------------- |
---|
125 | enddo |
---|
126 | |
---|
127 | c------------------------------------------------------------------------- |
---|
128 | c 2.2 Calcul des hi |
---|
129 | c ------------- |
---|
130 | c |
---|
131 | |
---|
132 | do i = 1 , nlaylte-1 |
---|
133 | do jl = 1 , kdlon |
---|
134 | c ------------------- |
---|
135 | hi(jl,i) = - semit * (g / pdp(jl,i) / cpp) * |
---|
136 | . ( xi(ig0+jl,1,i,i+1) * zdblay(jl,1,i+1) |
---|
137 | . + xi(ig0+jl,2,i,i+1) * zdblay(jl,2,i+1) ) |
---|
138 | c ------------------- |
---|
139 | enddo |
---|
140 | enddo |
---|
141 | |
---|
142 | c------------------------------------------------------------------------- |
---|
143 | c 2.3 Calcul des bi |
---|
144 | c ------------- |
---|
145 | c |
---|
146 | |
---|
147 | |
---|
148 | do i = 2 , nlaylte |
---|
149 | do jl = 1 , kdlon |
---|
150 | c ------------------- |
---|
151 | bi(jl,i) = - semit * (g / pdp(jl,i) / cpp) * |
---|
152 | . ( xi(ig0+jl,1,i,i-1) * zdblay(jl,1,i-1) |
---|
153 | . + xi(ig0+jl,2,i,i-1) * zdblay(jl,2,i-1) ) |
---|
154 | c ------------------- |
---|
155 | enddo |
---|
156 | enddo |
---|
157 | |
---|
158 | |
---|
159 | c couche 1 |
---|
160 | c -------- |
---|
161 | c tant qu'on a pas un calcul propre de zdblay(0) qui tienne compte de |
---|
162 | c la discontinuite de temperature au sol , on met b1 = 0 |
---|
163 | |
---|
164 | |
---|
165 | do jl = 1 , kdlon |
---|
166 | bi(jl,1) = 0 |
---|
167 | enddo |
---|
168 | |
---|
169 | c------------------------------------------------------------------------- |
---|
170 | c 2.4 |
---|
171 | c ------------- |
---|
172 | c |
---|
173 | |
---|
174 | c couche nlaylte |
---|
175 | c ------------ |
---|
176 | |
---|
177 | do jl = 1 , kdlon |
---|
178 | c ------------------- |
---|
179 | ci(jl,nlaylte) = (gcp * psi(jl,nlaylte) / pdp(jl,nlaylte)) |
---|
180 | . / di(jl,nlaylte) |
---|
181 | |
---|
182 | ai(jl,nlaylte) = - bi(jl,nlaylte) / di(jl,nlaylte) |
---|
183 | c ------------------- |
---|
184 | enddo |
---|
185 | |
---|
186 | |
---|
187 | |
---|
188 | do i = nlaylte-1 , 1 , -1 |
---|
189 | do jl = 1 , kdlon |
---|
190 | c ------------------- |
---|
191 | denom = di(jl,i) + hi(jl,i) * ai(jl,i+1) |
---|
192 | |
---|
193 | ci(jl,i) = ( gcp * psi(jl,i) / pdp(jl,i) |
---|
194 | . - hi(jl,i) * ci(jl,i+1) ) / denom |
---|
195 | |
---|
196 | ai(jl,i) = -bi(jl,i) / denom |
---|
197 | c ------------------- |
---|
198 | enddo |
---|
199 | enddo |
---|
200 | |
---|
201 | |
---|
202 | c------------------------------------------------------------------------- |
---|
203 | c 2.5 |
---|
204 | c ------------- |
---|
205 | c |
---|
206 | |
---|
207 | c couche 1 |
---|
208 | c ------- |
---|
209 | do jl = 1 , kdlon |
---|
210 | newpcolc(jl,1) = ci(jl,1) |
---|
211 | enddo |
---|
212 | |
---|
213 | |
---|
214 | do i = 2 , nlaylte |
---|
215 | do jl = 1 , kdlon |
---|
216 | newpcolc(jl,i) = ci(jl,i) + ai(jl,i) * newpcolc(jl,i-1) |
---|
217 | enddo |
---|
218 | enddo |
---|
219 | |
---|
220 | |
---|
221 | |
---|
222 | c------------------------------------------------------------------------- |
---|
223 | |
---|
224 | end subroutine lwi |
---|
225 | |
---|
226 | end module lwi_mod |
---|