/[lmdze]/trunk/Sources/filtrez/inifilr.f
ViewVC logotype

Contents of /trunk/Sources/filtrez/inifilr.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 164 - (show annotations)
Tue Jul 28 14:53:31 2015 UTC (8 years, 9 months ago) by guez
File size: 8215 byte(s)
In procedure inifilr, coefilu2 and coefilv2 were not used. coefilu and
coefilv were defined and used only at filtered latitudes so split them
into north and south arrays. Values in eignvl are necessarily
negative. Simplified the computation of eignft.

1 module inifilr_m
2
3 IMPLICIT NONE
4
5 ! North:
6
7 INTEGER jfiltnu, jfiltnv
8 ! index of the last scalar line filtered in northern hemisphere
9
10 real, allocatable:: matriceun(:, :, :), matrinvn(:, :, :)
11 ! (iim, iim, 2:jfiltnu)
12
13 real, allocatable:: matricevn(:, :, :) ! (iim, iim, jfiltnv)
14
15 ! South:
16
17 integer jfiltsu, jfiltsv
18 ! index of the first line filtered in southern hemisphere
19
20 real, allocatable:: matriceus(:, :, :), matrinvs(:, :, :)
21 ! (iim, iim, jfiltsu:jjm)
22
23 real, allocatable:: matricevs(:, :, :) ! (iim, iim, jfiltsv:jjm)
24
25 contains
26
27 SUBROUTINE inifilr
28
29 ! From filtrez/inifilr.F, version 1.1.1.1, 2004/05/19 12:53:09
30 ! H. Upadhyaya, O. Sharma
31
32 ! This procedure computes the filtering coefficients for scalar
33 ! lines and meridional wind v lines. The modes are filtered from
34 ! modfrst to iim. We filter all those latitude lines where coefil
35 ! < 1. No filtering at poles. colat0 is to be used when alpha
36 ! (stretching coefficient) is set equal to zero for the regular
37 ! grid case.
38
39 USE dimens_m, ONLY : iim, jjm
40 USE dynetat0_m, ONLY : rlatu, rlatv, xprimu, grossismx
41 use inifgn_m, only: inifgn
42 use jumble, only: new_unit
43 use nr_util, only: pi
44
45 ! Local:
46
47 REAL dlatu(jjm)
48 REAL rlamda(2: iim)
49 real eignvl(iim) ! eigenvalues sorted in descending order (<= 0)
50 INTEGER i, j, unit
51 REAL colat0 ! > 0
52 REAL eignft(iim, iim)
53
54 real eignfnu(iim, iim), eignfnv(iim, iim)
55 ! eigenvectors of the discrete second derivative with respect to longitude
56
57 ! Filtering coefficients (lamda_max * cos(rlat) / lamda):
58 real, allocatable:: coefilnu(:, :) ! (iim, 2:jfiltnu)
59 real, allocatable:: coefilsu(:, :) ! (iim, jfiltsu:jjm)
60 real, allocatable:: coefilnv(:, :) ! (iim, jfiltnv)
61 real, allocatable:: coefilsv(:, :) ! (iim, jfiltsv:jjm)
62
63 ! Index of the mode from where modes are filtered:
64 integer, allocatable:: modfrstnu(:) ! (2:jfiltnu)
65 integer, allocatable:: modfrstsu(:) ! (jfiltsu:jjm)
66 integer, allocatable:: modfrstnv(:) ! (jfiltnv)
67 integer, allocatable:: modfrstsv(:) ! (jfiltsv:jjm)
68
69 !-----------------------------------------------------------
70
71 print *, "Call sequence information: inifilr"
72
73 CALL inifgn(eignvl, eignfnu, eignfnv)
74
75 ! Calcul de colat0
76 forall (j = 1:jjm) dlatu(j) = rlatu(j) - rlatu(j + 1)
77 colat0 = min(0.5, minval(dlatu) / minval(xprimu(:iim)))
78 PRINT *, 'colat0 = ', colat0
79
80 rlamda = iim / (pi * colat0 / grossismx) / sqrt(- eignvl(2: iim))
81
82 ! Determination de jfiltnu, jfiltsu, jfiltnv, jfiltsv
83
84 jfiltnu = (jjm + 1) / 2
85 do while (cos(rlatu(jfiltnu)) >= colat0 &
86 .or. rlamda(iim) * cos(rlatu(jfiltnu)) >= 1.)
87 jfiltnu = jfiltnu - 1
88 end do
89
90 jfiltsu = jjm / 2 + 2
91 do while (cos(rlatu(jfiltsu)) >= colat0 &
92 .or. rlamda(iim) * cos(rlatu(jfiltsu)) >= 1.)
93 jfiltsu = jfiltsu + 1
94 end do
95
96 jfiltnv = jjm / 2
97 do while ((cos(rlatv(jfiltnv)) >= colat0 &
98 .or. rlamda(iim) * cos(rlatv(jfiltnv)) >= 1.) .and. jfiltnv >= 2)
99 jfiltnv = jfiltnv - 1
100 end do
101
102 if (cos(rlatv(jfiltnv)) >= colat0 &
103 .or. rlamda(iim) * cos(rlatv(jfiltnv)) >= 1.) then
104 ! {jfiltnv == 1}
105 PRINT *, 'Could not find jfiltnv.'
106 STOP 1
107 END IF
108
109 jfiltsv = (jjm + 1)/ 2 + 1
110 do while ((cos(rlatv(jfiltsv)) >= colat0 &
111 .or. rlamda(iim) * cos(rlatv(jfiltsv)) >= 1.) .and. jfiltsv <= jjm - 1)
112 jfiltsv = jfiltsv + 1
113 end do
114
115 IF (cos(rlatv(jfiltsv)) >= colat0 &
116 .or. rlamda(iim) * cos(rlatv(jfiltsv)) >= 1.) THEN
117 ! {jfiltsv == jjm}
118 PRINT *, 'Could not find jfiltsv.'
119 STOP 1
120 END IF
121
122 PRINT *, 'jfiltnu =', jfiltnu
123 PRINT *, 'jfiltsu =', jfiltsu
124 PRINT *, 'jfiltnv =', jfiltnv
125 PRINT *, 'jfiltsv =', jfiltsv
126
127 ! D\'etermination de coefil[ns][uv], modfrst[ns][uv]:
128
129 allocate(modfrstnu(2:jfiltnu), modfrstsu(jfiltsu:jjm))
130 allocate(modfrstnv(jfiltnv), modfrstsv(jfiltsv:jjm))
131 allocate(coefilnu(iim, 2:jfiltnu), coefilsu(iim, jfiltsu:jjm))
132 allocate(coefilnv(iim, jfiltnv), coefilsv(iim, jfiltsv:jjm))
133
134 coefilnu = 0.
135 coefilnv = 0.
136 coefilsu = 0.
137 coefilsv = 0.
138
139 DO j = 2, jfiltnu
140 modfrstnu(j) = 2
141 do while (rlamda(modfrstnu(j)) * cos(rlatu(j)) >= 1. &
142 .and. modfrstnu(j) <= iim - 1)
143 modfrstnu(j) = modfrstnu(j) + 1
144 end do
145
146 if (rlamda(modfrstnu(j)) * cos(rlatu(j)) < 1.) then
147 DO i = modfrstnu(j), iim
148 coefilnu(i, j) = rlamda(i) * cos(rlatu(j)) - 1.
149 end DO
150 end if
151 END DO
152
153 DO j = 1, jfiltnv
154 modfrstnv(j) = 2
155 do while (rlamda(modfrstnv(j)) * cos(rlatv(j)) >= 1. &
156 .and. modfrstnv(j) <= iim - 1)
157 modfrstnv(j) = modfrstnv(j) + 1
158 end do
159
160 if (rlamda(modfrstnv(j)) * cos(rlatv(j)) < 1.) then
161 DO i = modfrstnv(j), iim
162 coefilnv(i, j) = rlamda(i) * cos(rlatv(j)) - 1.
163 end DO
164 end if
165 end DO
166
167 DO j = jfiltsu, jjm
168 modfrstsu(j) = 2
169 do while (rlamda(modfrstsu(j)) * cos(rlatu(j)) >= 1. &
170 .and. modfrstsu(j) <= iim - 1)
171 modfrstsu(j) = modfrstsu(j) + 1
172 end do
173
174 if (rlamda(modfrstsu(j)) * cos(rlatu(j)) < 1.) then
175 DO i = modfrstsu(j), iim
176 coefilsu(i, j) = rlamda(i) * cos(rlatu(j)) - 1.
177 end DO
178 end if
179 end DO
180
181 DO j = jfiltsv, jjm
182 modfrstsv(j) = 2
183 do while (rlamda(modfrstsv(j)) * cos(rlatv(j)) >= 1. &
184 .and. modfrstsv(j) <= iim - 1)
185 modfrstsv(j) = modfrstsv(j) + 1
186 end do
187
188 if (rlamda(modfrstsv(j)) * cos(rlatv(j)) < 1.) then
189 DO i = modfrstsv(j), iim
190 coefilsv(i, j) = rlamda(i) * cos(rlatv(j)) - 1.
191 end DO
192 end if
193 END DO
194
195 call new_unit(unit)
196 open(unit, file = "inifilr_out.txt", status = "replace", action = "write")
197 write(unit, fmt = *) '"EIGNVL"', eignvl
198 write(unit, fmt = *) '"modfrstnu"', modfrstnu
199 write(unit, fmt = *) '"modfrstsu"', modfrstsu
200 write(unit, fmt = *) '"modfrstnv"', modfrstnv
201 write(unit, fmt = *) '"modfrstsv"', modfrstsv
202 close(unit)
203
204 allocate(matriceun(iim, iim, 2:jfiltnu), matrinvn(iim, iim, 2:jfiltnu))
205 allocate(matricevn(iim, iim, jfiltnv))
206 allocate(matricevs(iim, iim, jfiltsv:jjm))
207 allocate(matriceus(iim, iim, jfiltsu:jjm), matrinvs(iim, iim, jfiltsu:jjm))
208
209 ! Calcul de la matrice filtre 'matriceu' pour les champs situes
210 ! sur la grille scalaire
211
212 DO j = 2, jfiltnu
213 eignft(:modfrstnu(j) - 1, :) = 0.
214 forall (i = modfrstnu(j):iim) eignft(i, :) = eignfnv(:, i) &
215 * coefilnu(i, j)
216 matriceun(:, :, j) = matmul(eignfnv, eignft)
217 END DO
218
219 DO j = jfiltsu, jjm
220 eignft(:modfrstsu(j) - 1, :) = 0.
221 forall (i = modfrstsu(j):iim) eignft(i, :) = eignfnv(:, i) &
222 * coefilsu(i, j)
223 matriceus(:, :, j) = matmul(eignfnv, eignft)
224 END DO
225
226 ! Calcul de la matrice filtre 'matricev' pour les champs situes
227 ! sur la grille de V ou de Z
228
229 DO j = 1, jfiltnv
230 eignft(:modfrstnv(j) - 1, :) = 0.
231 forall (i = modfrstnv(j): iim) eignft(i, :) = eignfnu(:, i) &
232 * coefilnv(i, j)
233 matricevn(:, :, j) = matmul(eignfnu, eignft)
234 END DO
235
236 DO j = jfiltsv, jjm
237 eignft(:modfrstsv(j) - 1, :) = 0.
238 forall (i = modfrstsv(j):iim) eignft(i, :) = eignfnu(:, i) &
239 * coefilsv(i, j)
240 matricevs(:, :, j) = matmul(eignfnu, eignft)
241 END DO
242
243 ! Calcul de la matrice filtre 'matrinv' pour les champs situes
244 ! sur la grille scalaire , pour le filtre inverse
245
246 DO j = 2, jfiltnu
247 eignft(:modfrstnu(j) - 1, :) = 0.
248 forall (i = modfrstnu(j):iim) eignft(i, :) = eignfnv(:, i) &
249 * coefilnu(i, j) / (1. + coefilnu(i, j))
250 matrinvn(:, :, j) = matmul(eignfnv, eignft)
251 END DO
252
253 DO j = jfiltsu, jjm
254 eignft(:modfrstsu(j) - 1, :) = 0.
255 forall (i = modfrstsu(j):iim) eignft(i, :) = eignfnv(:, i) &
256 * coefilsu(i, j) / (1. + coefilsu(i, j))
257 matrinvs(:, :, j) = matmul(eignfnv, eignft)
258 END DO
259
260 END SUBROUTINE inifilr
261
262 end module inifilr_m

  ViewVC Help
Powered by ViewVC 1.1.21