1 | !!---------------------------------------------------------------------- |
---|
2 | !! *** solisl_fdir.h90 *** |
---|
3 | !!---------------------------------------------------------------------- |
---|
4 | |
---|
5 | !!---------------------------------------------------------------------- |
---|
6 | !! islmat : Compute the matrix associated with islands |
---|
7 | !! and save it in a direct access file |
---|
8 | !! islbsf : Compute the barotropic streamfunction of each island |
---|
9 | !! and save it in a direct access file |
---|
10 | !!---------------------------------------------------------------------- |
---|
11 | |
---|
12 | SUBROUTINE islmat |
---|
13 | !!---------------------------------------------------------------------- |
---|
14 | !! *** ROUTINE islmat *** |
---|
15 | !! |
---|
16 | !! ** Purpose : Compute and invert the island matrix |
---|
17 | !! |
---|
18 | !! ** Method : aisl(jni,jnj) matrix constituted by bsfisl circulation |
---|
19 | !! |
---|
20 | !! ** Action : - aisl : island matrix |
---|
21 | !! - aislm1 : invert of the island matrix |
---|
22 | !! |
---|
23 | !! ** file : islands: contain bsfisl, aisl and aislm1 |
---|
24 | !! |
---|
25 | !! History : |
---|
26 | !! 1.0 ! 88-10 (G. Madec) Original code |
---|
27 | !! 7.0 ! 96-01 (G. Madec) suppression of common work arrays |
---|
28 | !! 8.5 ! 02-08 (G. Madec) Free form, F90 |
---|
29 | !!---------------------------------------------------------------------- |
---|
30 | !! * Local declarations |
---|
31 | INTEGER :: ji, jj, jni, jnj, ios, jl |
---|
32 | REAL(wp) :: zwx(jpi,jpj,2) |
---|
33 | !!---------------------------------------------------------------------- |
---|
34 | |
---|
35 | |
---|
36 | ! I. Island matrix lecture in numisl (if it exists) |
---|
37 | ! ================================== |
---|
38 | ! file already open in islbsf routine |
---|
39 | ! Lecture |
---|
40 | READ( numisl, REC=jpisl+2, IOSTAT=ios ) aisl, aislm1 |
---|
41 | IF( ios /= 0 ) GOTO 110 |
---|
42 | |
---|
43 | ! Control print |
---|
44 | IF(lwp) THEN |
---|
45 | WRITE(numout,*) |
---|
46 | WRITE(numout,*) ' islmat: lecture aisl/aislm1 in numisl done' |
---|
47 | WRITE(numout,*) ' ~~~~~~' |
---|
48 | WRITE(numout,*) |
---|
49 | WRITE(numout,*) ' island matrix' |
---|
50 | WRITE(numout,*) ' -------------' |
---|
51 | WRITE(numout,*) |
---|
52 | |
---|
53 | DO jnj = 1, jpisl |
---|
54 | WRITE(numout,'(8e12.4)') ( aisl(jni,jnj), jni = 1, jpisl ) |
---|
55 | END DO |
---|
56 | |
---|
57 | WRITE(numout,*) |
---|
58 | WRITE(numout,*) ' inverse of the island matrix' |
---|
59 | WRITE(numout,*) ' ----------------------------' |
---|
60 | WRITE(numout,*) |
---|
61 | |
---|
62 | DO jnj = 1, jpisl |
---|
63 | WRITE(numout,'(12e11.3)') ( aislm1(jni,jnj), jni=1, jpisl ) |
---|
64 | END DO |
---|
65 | ENDIF |
---|
66 | |
---|
67 | CLOSE(numisl) |
---|
68 | |
---|
69 | RETURN |
---|
70 | |
---|
71 | 110 CONTINUE |
---|
72 | |
---|
73 | |
---|
74 | ! II. Island matrix computation |
---|
75 | ! ============================= |
---|
76 | |
---|
77 | DO jnj = 1, jpisl |
---|
78 | |
---|
79 | ! 1.1 Circulation of bsf(jnj) around island jni |
---|
80 | |
---|
81 | DO jj = 2, jpj |
---|
82 | zwx(:,jj,1) = -( bsfisl(:,jj,jnj) - bsfisl(:,jj-1,jnj) ) |
---|
83 | END DO |
---|
84 | DO ji = 1, jpi |
---|
85 | zwx(ji,1,1) = 0.e0 |
---|
86 | END DO |
---|
87 | |
---|
88 | DO jj = 1, jpj |
---|
89 | DO ji = 2, jpi |
---|
90 | zwx(ji,jj,2) = ( bsfisl(ji,jj,jnj) - bsfisl(ji-1,jj,jnj) ) |
---|
91 | END DO |
---|
92 | END DO |
---|
93 | DO jj = 1, jpj |
---|
94 | zwx(1,jj,2) = 0.e0 |
---|
95 | END DO |
---|
96 | |
---|
97 | ! 1.2 Island matrix |
---|
98 | |
---|
99 | DO jni = 1, jpisl |
---|
100 | aisl(jni,jnj) = 0. |
---|
101 | DO jl = 1, 2 |
---|
102 | DO jj = 1, jpj |
---|
103 | DO ji = 1, jpi |
---|
104 | aisl(jni,jnj) = aisl(jni,jnj) + acisl2(ji,jj,jl,jni) * zwx(ji,jj,jl) |
---|
105 | END DO |
---|
106 | END DO |
---|
107 | END DO |
---|
108 | END DO |
---|
109 | |
---|
110 | END DO |
---|
111 | IF( lk_mpp ) CALL mpp_sum( aisl, jpisl*jpisl ) ! sum over the global domain |
---|
112 | |
---|
113 | ! 1.3 Control print |
---|
114 | IF(lwp) THEN |
---|
115 | WRITE(numout,*) |
---|
116 | WRITE(numout,*) ' islmat: island matrix' |
---|
117 | WRITE(numout,*) ' ~~~~~~' |
---|
118 | WRITE(numout,*) |
---|
119 | |
---|
120 | DO jnj = 1, jpisl |
---|
121 | WRITE(numout,'(8e12.4)') ( aisl(jni,jnj), jni = 1, jpisl ) |
---|
122 | END DO |
---|
123 | ENDIF |
---|
124 | |
---|
125 | |
---|
126 | ! 2. Invertion of the island matrix |
---|
127 | ! --------------------------------- |
---|
128 | |
---|
129 | ! 2.1 Call of an imsl routine for the matrix invertion |
---|
130 | |
---|
131 | CALL linrg( jpisl, aisl, jpisl, aislm1, jpisl ) |
---|
132 | |
---|
133 | ! 2.2 Control print |
---|
134 | |
---|
135 | IF(lwp) THEN |
---|
136 | WRITE(numout,*) |
---|
137 | WRITE(numout,*) ' islmat: inverse of the island matrix' |
---|
138 | WRITE(numout,*) ' ~~~~~~' |
---|
139 | WRITE(numout,*) |
---|
140 | |
---|
141 | DO jnj = 1, jpisl |
---|
142 | WRITE(numout, '(12e11.3)') ( aislm1(jni,jnj), jni=1, jpisl ) |
---|
143 | END DO |
---|
144 | ENDIF |
---|
145 | |
---|
146 | |
---|
147 | ! 3. Output of aisl and aislm1 in numisl |
---|
148 | ! -------------------------------------- |
---|
149 | |
---|
150 | IF(lwp) WRITE(numisl,REC=jpisl+2) aisl, aislm1 |
---|
151 | |
---|
152 | CLOSE( numisl ) |
---|
153 | |
---|
154 | END SUBROUTINE islmat |
---|
155 | |
---|
156 | |
---|
157 | SUBROUTINE islbsf |
---|
158 | !!---------------------------------------------------------------------- |
---|
159 | !! *** ROUTINE islbsf *** |
---|
160 | !! |
---|
161 | !! ** Purpose : Compute the barotropic stream function associated with |
---|
162 | !! each island using FETI or preconditioned conjugate gradient |
---|
163 | !! method or read them in numisl. |
---|
164 | !! |
---|
165 | !! ** Method : Direct acces file |
---|
166 | !! |
---|
167 | !! ** file : numisl: barotropic stream function associated |
---|
168 | !! with each island of the domain |
---|
169 | !! |
---|
170 | !! ** Action : - bsfisl, the streamfunction which takes the value 1 |
---|
171 | !! over island ni, and 0 over the others islands |
---|
172 | !! - update 'numisl' file |
---|
173 | !! |
---|
174 | !! History : |
---|
175 | !! ! 87-10 (G. Madec) Original code |
---|
176 | !! ! 91-11 (G. Madec) |
---|
177 | !! ! 93-03 (G. Madec) release 7.1 |
---|
178 | !! ! 93-04 (M. Guyon) loops and suppress pointers |
---|
179 | !! ! 96-11 (A. Weaver) correction to preconditioning |
---|
180 | !! ! 98-02 (M. Guyon) FETI method |
---|
181 | !! ! 99-11 (M. Imbard) NetCDF FORMAT with IOIPSL |
---|
182 | !! 8.5 ! 02-08 (G. Madec) Free form, F90 |
---|
183 | !!---------------------------------------------------------------------- |
---|
184 | !! * Local declarations |
---|
185 | INTEGER :: ji, jj, jni, jii, jnp |
---|
186 | INTEGER :: iimlu, ijmlu, inmlu, ios, iju, iloc |
---|
187 | INTEGER :: ii, ij, icile, icut, inmax, indic |
---|
188 | REAL(wp) :: zepsr, zeplu, zgwgt |
---|
189 | REAL(wp) :: zep(jpisl) |
---|
190 | !!---------------------------------------------------------------------- |
---|
191 | |
---|
192 | |
---|
193 | ! 0. Initializations |
---|
194 | ! ================== |
---|
195 | |
---|
196 | ! set to zero the convergence indicator |
---|
197 | icile = 0 |
---|
198 | |
---|
199 | ! set to zero of bsfisl |
---|
200 | bsfisl(:,:,1:jpisl) = 0.e0 |
---|
201 | |
---|
202 | |
---|
203 | ! I. Lecture of bsfisl in numisl (if it exists) |
---|
204 | ! ============================================= |
---|
205 | |
---|
206 | ! open islands file |
---|
207 | ibloc = 4096 |
---|
208 | ilglo = ibloc*((jpiglo*jpjglo*jpbyt-1 )/ibloc+1) |
---|
209 | clname = 'islands' |
---|
210 | CALL ctlopn( numisl, 'islands', 'UNKNOWN', 'UNFORMATTED', 'DIRECT', & |
---|
211 | ilglo,numout,lwp,0) |
---|
212 | |
---|
213 | icut = 0 |
---|
214 | iimlu = 0 |
---|
215 | ijmlu = 0 |
---|
216 | inmlu = 0 |
---|
217 | zeplu = 0. |
---|
218 | READ (numisl, REC=1, IOSTAT=ios) iimlu, ijmlu, inmlu, zeplu |
---|
219 | IF( ios /= 0 ) GOTO 110 |
---|
220 | DO jni = 1, jpisl |
---|
221 | CALL read2( numisl, bsfisl(1,1,jni), 1, jni+1 ) |
---|
222 | END DO |
---|
223 | |
---|
224 | 110 CONTINUE |
---|
225 | |
---|
226 | ! the read precision is not the required one : icut=888 |
---|
227 | IF( zeplu > epsisl ) icut = 888 |
---|
228 | |
---|
229 | ! the read domain does not correspond to the model one : icut=999 |
---|
230 | IF( iimlu /= jpiglo .OR. ijmlu /= jpjglo .OR. inmlu /= jpisl ) icut = 999 |
---|
231 | |
---|
232 | ! Control print |
---|
233 | IF( icut == 999 ) THEN |
---|
234 | IF(lwp) THEN |
---|
235 | WRITE(numout,*) |
---|
236 | WRITE(numout,*) ' islbsf : lecture bsfisl in numisl failed' |
---|
237 | WRITE(numout,*) ' ~~~~~~' |
---|
238 | WRITE(numout,*) ' icut= ', icut |
---|
239 | WRITE(numout,*) ' imlu= ', iimlu, ' jmlu= ', ijmlu, & |
---|
240 | ' nilu= ', inmlu, ' epsisl lu= ', zeplu |
---|
241 | WRITE(numout,*) ' the bsfisl are computed from zero' |
---|
242 | ENDIF |
---|
243 | ELSEIF( icut == 888 ) THEN |
---|
244 | IF(lwp) THEN |
---|
245 | WRITE(numout,*) |
---|
246 | WRITE(numout,*) ' islbsf : lecture bsfisl in numisl done' |
---|
247 | WRITE(numout,*) ' ~~~~~~' |
---|
248 | WRITE(numout,*) ' the required accuracy is not reached' |
---|
249 | WRITE(numout,*) ' epsisl lu= ', zeplu,' epsisl required ', epsisl |
---|
250 | WRITE(numout,*) ' the bsfisl are computed from the read values' |
---|
251 | ENDIF |
---|
252 | ELSE |
---|
253 | IF(lwp) THEN |
---|
254 | WRITE(numout,*) |
---|
255 | WRITE(numout,*) ' islbsf : lecture bsfisl in numisl done' |
---|
256 | WRITE(numout,*) ' ~~~~~~' |
---|
257 | ENDIF |
---|
258 | RETURN |
---|
259 | ENDIF |
---|
260 | |
---|
261 | |
---|
262 | ! II. Compute the bsfisl (if icut=888 or 999) |
---|
263 | ! ====================== |
---|
264 | |
---|
265 | ! save nmax |
---|
266 | inmax = nmax |
---|
267 | |
---|
268 | ! set the number of iteration of island computation |
---|
269 | nmax = nmisl |
---|
270 | |
---|
271 | ! Loop over islands |
---|
272 | ! ----------------- |
---|
273 | |
---|
274 | DO jni = 1, jpisl |
---|
275 | |
---|
276 | |
---|
277 | ! 1. Initalizations of island computation |
---|
278 | ! --------------------------------------- |
---|
279 | |
---|
280 | ! 1.0 Set the pcg solution gcb to zero |
---|
281 | gcb(:,:) = 0.e0 |
---|
282 | |
---|
283 | ! 1.1 Set first guess gcx either to zero or to the read bsfisl |
---|
284 | |
---|
285 | IF( icut == 999 ) THEN |
---|
286 | gcx(:,:) = 0.e0 |
---|
287 | ELSEIF( icut == 888 ) THEN |
---|
288 | IF(lwp) WRITE(numout,*) ' islbsf: bsfisl read are used as first guess' |
---|
289 | ! C a u t i o n : bsfisl masked because it contains 1 along island seaside |
---|
290 | gcx(:,:) = bsfisl(:,:,jni) * bmask(:,:) |
---|
291 | ENDIF |
---|
292 | |
---|
293 | ! 1.2 Right hand side of the stream FUNCTION equation |
---|
294 | |
---|
295 | IF( lk_mpp ) THEN |
---|
296 | ! north fold treatment |
---|
297 | IF( npolj == 3 ) iloc = jpiglo -(nimpp-1+nimppt(nono+1)-1) |
---|
298 | IF( npolj == 4 ) iloc = jpiglo - 2*(nimpp-1) |
---|
299 | t2p1(:,1,1) = 0. |
---|
300 | ! north and south grid-points |
---|
301 | DO jii = 1, 2 |
---|
302 | DO jnp = 1, mnisl(jii,jni) |
---|
303 | ii = miisl(jnp,jii,jni) |
---|
304 | ij = mjisl(jnp,jii,jni) |
---|
305 | IF( ( npolj == 3 .OR. npolj == 4 ) .AND. ( ij == nlcj-1 .AND. jii == 1) ) THEN |
---|
306 | iju=iloc-ii+1 |
---|
307 | t2p1(iju,1,1) = t2p1(iju,1,1) + hur(ii,ij) * e1u(ii,ij) / e2u(ii,ij) |
---|
308 | ELSE |
---|
309 | gcb(ii,ij-jii+1) = gcb(ii,ij-jii+1) + hur(ii,ij) * e1u(ii,ij) / e2u(ii,ij) |
---|
310 | ENDIF |
---|
311 | END DO |
---|
312 | END DO |
---|
313 | |
---|
314 | ! east and west grid-points |
---|
315 | |
---|
316 | DO jii = 3, 4 |
---|
317 | DO jnp = 1, mnisl(jii,jni) |
---|
318 | ii = miisl(jnp,jii,jni) |
---|
319 | ij = mjisl(jnp,jii,jni) |
---|
320 | gcb(ii-jii+3,ij) = gcb(ii-jii+3,ij) + hvr(ii,ij) * e2v(ii,ij) / e1v(ii,ij) |
---|
321 | END DO |
---|
322 | END DO |
---|
323 | |
---|
324 | IF( lk_mpp ) CALL mpplnks( gcb ) !!bug ? should use an lbclnk ? is it possible??? |
---|
325 | |
---|
326 | ELSE |
---|
327 | ! north and south grid-points |
---|
328 | DO jii = 1, 2 |
---|
329 | DO jnp = 1, mnisl(jii,jni) |
---|
330 | ii = miisl(jnp,jii,jni) |
---|
331 | ij = mjisl(jnp,jii,jni) |
---|
332 | IF( ( nperio == 3 .OR. nperio == 4 ) .AND. ( ij == jpj-1 .AND. jii == 1) ) THEN |
---|
333 | gcb(jpi-ii+1,ij-1) = gcb(jpi-ii+1,ij-1) + hur(ii,ij) * e1u(ii,ij) / e2u(ii,ij) |
---|
334 | ELSE |
---|
335 | gcb(ii,ij-jii+1) = gcb(ii,ij-jii+1) + hur(ii,ij) * e1u(ii,ij) / e2u(ii,ij) |
---|
336 | ENDIF |
---|
337 | END DO |
---|
338 | END DO |
---|
339 | |
---|
340 | ! east and west grid-points |
---|
341 | DO jii = 3, 4 |
---|
342 | DO jnp = 1, mnisl(jii,jni) |
---|
343 | ii = miisl(jnp,jii,jni) |
---|
344 | ij = mjisl(jnp,jii,jni) |
---|
345 | IF( bmask(ii-jii+3,ij) /= 0. ) THEN |
---|
346 | gcb(ii-jii+3,ij) = gcb(ii-jii+3,ij) + hvr(ii,ij) * e2v(ii,ij) / e1v(ii,ij) |
---|
347 | ELSE |
---|
348 | ! east-west cyclic boundary conditions |
---|
349 | IF( ii-jii+3 == 1 ) THEN |
---|
350 | gcb(jpim1,ij) = gcb(jpim1,ij) + hvr(ii,ij) * e2v(ii,ij) / e1v(ii,ij) |
---|
351 | ENDIF |
---|
352 | ENDIF |
---|
353 | END DO |
---|
354 | END DO |
---|
355 | ENDIF |
---|
356 | |
---|
357 | ! 1.4 Preconditioned right hand side and absolute precision |
---|
358 | |
---|
359 | IF( nsolv == 3 ) THEN |
---|
360 | ! FETI method |
---|
361 | ncut = 0 |
---|
362 | rnorme=0. |
---|
363 | DO jj = 1, jpj |
---|
364 | DO ji = 1, jpi |
---|
365 | gcb (ji,jj) = bmask(ji,jj) * gcb(ji,jj) |
---|
366 | rnorme = rnorme + gcb(ji,jj) * gcb(ji,jj) |
---|
367 | END DO |
---|
368 | END DO |
---|
369 | |
---|
370 | CALL mpp_sum( rnorme ) |
---|
371 | |
---|
372 | IF(lwp) WRITE(numout,*) 'rnorme ', rnorme |
---|
373 | epsr = epsisl * epsisl * rnorme |
---|
374 | indic = 0 |
---|
375 | ELSE |
---|
376 | ncut = 0 |
---|
377 | rnorme = 0. |
---|
378 | DO jj = 1, jpj |
---|
379 | DO ji = 1, jpi |
---|
380 | gcb (ji,jj) = gcdprc(ji,jj) * gcb(ji,jj) |
---|
381 | zgwgt = gcdmat(ji,jj) * gcb(ji,jj) |
---|
382 | rnorme = rnorme + gcb(ji,jj) * zgwgt |
---|
383 | END DO |
---|
384 | END DO |
---|
385 | IF( lk_mpp ) CALL mpp_sum( rnorme ) ! sum over the global domain |
---|
386 | |
---|
387 | IF(lwp) WRITE(numout,*) 'rnorme ', rnorme |
---|
388 | epsr = epsisl * epsisl * rnorme |
---|
389 | indic = 0 |
---|
390 | ENDIF |
---|
391 | |
---|
392 | |
---|
393 | ! 3. PCG solver for gcp.gcx=gcb in monotask |
---|
394 | ! ------------------------------------------- |
---|
395 | |
---|
396 | IF(nsolv == 3) THEN |
---|
397 | ! FETI method |
---|
398 | ! precision to compute Islands matrix A |
---|
399 | |
---|
400 | epsilo = epsisl |
---|
401 | CALL solfet( indic ) |
---|
402 | ! ---> |
---|
403 | ! precision to compute grad PS |
---|
404 | |
---|
405 | epsilo = eps |
---|
406 | ELSE |
---|
407 | |
---|
408 | CALL solpcg( indic ) |
---|
409 | |
---|
410 | ENDIF |
---|
411 | |
---|
412 | |
---|
413 | ! 4. Save the solution in bsfisl |
---|
414 | ! ------------------------------ |
---|
415 | |
---|
416 | DO jj = 1, jpj |
---|
417 | DO ji = 1, jpi |
---|
418 | bsfisl(ji,jj,jni) = gcx(ji,jj) |
---|
419 | END DO |
---|
420 | END DO |
---|
421 | |
---|
422 | |
---|
423 | ! 5. Boundary conditions |
---|
424 | ! ---------------------- |
---|
425 | |
---|
426 | ! set to 1. coastal gridpoints of the island |
---|
427 | DO jnp = 1, mnisl(0,jni) |
---|
428 | ii = miisl(jnp,0,jni) |
---|
429 | ij = mjisl(jnp,0,jni) |
---|
430 | bsfisl(ii,ij,jni) = 1. |
---|
431 | END DO |
---|
432 | |
---|
433 | ! cyclic boundary conditions |
---|
434 | IF( nperio == 1 .OR. nperio == 4 ) THEN |
---|
435 | DO jj = 1, jpj |
---|
436 | bsfisl( 1 ,jj,jni) = bsfisl(jpim1,jj,jni) |
---|
437 | bsfisl(jpi,jj,jni) = bsfisl( 2 ,jj,jni) |
---|
438 | END DO |
---|
439 | ENDIF |
---|
440 | IF( nperio == 3 .OR. nperio == 4 ) THEN |
---|
441 | DO ji = 1, jpim1 |
---|
442 | iju = jpi-ji+1 |
---|
443 | bsfisl(ji,jpj-1,jni) = bsfisl(iju,jpj-2,jni) |
---|
444 | bsfisl(ji, jpj ,jni) = bsfisl(iju,jpj-3,jni) |
---|
445 | END DO |
---|
446 | ENDIF |
---|
447 | IF( lk_mpp ) CALL lbc_lnk( bsfisl(:,:,jni), 'G', 1. ) ! link at G-point |
---|
448 | |
---|
449 | |
---|
450 | ! 6. Control print |
---|
451 | ! ---------------- |
---|
452 | |
---|
453 | IF(lwp)WRITE(numout,*) |
---|
454 | IF(lwp)WRITE(numout,*) ' islbsf: island number: ', jni |
---|
455 | IF(lwp)WRITE (numout,9290) niter, res, SQRT(epsr)/epsisl |
---|
456 | zep(jni) = MAX(epsisl, res/(SQRT(epsr)/epsisl)) |
---|
457 | IF(indic < 0) THEN |
---|
458 | icile = icile-1 |
---|
459 | IF(lwp)WRITE(numout,*) ' pcg do not converge for island: ', jni |
---|
460 | IF(lwp)WRITE(numout,*) ' Precision reached: ',zep(jni) |
---|
461 | ENDIF |
---|
462 | |
---|
463 | 9290 FORMAT(' niter :',i4,' , res :',e20.10,' , gcb :',e20.10) |
---|
464 | |
---|
465 | |
---|
466 | ! End loop for islands |
---|
467 | ! -------------------- |
---|
468 | |
---|
469 | END DO |
---|
470 | |
---|
471 | |
---|
472 | ! 7. Reset PCG |
---|
473 | ! ------------ |
---|
474 | |
---|
475 | ! reset the number of iteration for pcg |
---|
476 | nmax = inmax |
---|
477 | |
---|
478 | ! reset to zero pcg arrays |
---|
479 | gcx (:,:) = 0.e0 |
---|
480 | gcxb (:,:) = 0.e0 |
---|
481 | gcb (:,:) = 0.e0 |
---|
482 | gcr (:,:) = 0.e0 |
---|
483 | gcdes(:,:) = 0.e0 |
---|
484 | gccd (:,:) = 0.e0 |
---|
485 | |
---|
486 | |
---|
487 | ! III. Output of bsfisl in numisl |
---|
488 | ! =============================== |
---|
489 | |
---|
490 | IF( icile == 0 .AND. icut /= 0 ) THEN |
---|
491 | IF(lwp) THEN |
---|
492 | WRITE(numout,*) |
---|
493 | WRITE(numout,*)' islbsf : write bsfisl in numisl ', numisl |
---|
494 | WRITE(numout,*)' ~~~~~~' |
---|
495 | ENDIF |
---|
496 | |
---|
497 | IF(lwp) WRITE(numisl,REC=1) jpiglo, jpjglo, jpisl, epsisl |
---|
498 | DO jni = 1, jpisl |
---|
499 | CALL write2( numisl, bsfisl(1,1,jni), 1, jni+1 ) |
---|
500 | END DO |
---|
501 | ENDIF |
---|
502 | |
---|
503 | IF( icile < 0 ) THEN |
---|
504 | IF(lwp) THEN |
---|
505 | WRITE(numout,*) |
---|
506 | WRITE(numout,*) ' islbsf : number of island without convergence: ',ABS(icile) |
---|
507 | WRITE(numout,*) ' ~~~~~~' |
---|
508 | ENDIF |
---|
509 | zepsr = epsisl |
---|
510 | DO jni = 1, jpisl |
---|
511 | IF(lwp)WRITE(numout,*) ' isl ',jni,' precision reached ', zep(jni) |
---|
512 | zepsr = MAX( zep(jni), zepsr ) |
---|
513 | END DO |
---|
514 | IF( zepsr == 0. ) zepsr = epsisl |
---|
515 | IF(lwp) THEN |
---|
516 | WRITE(numout,*) ' save value of precision reached: ',zepsr |
---|
517 | WRITE(numout,*) |
---|
518 | WRITE(numout,*)' islbsf : save bsfisl in numisl ',numisl |
---|
519 | WRITE(numout,*)' ~~~~~~' |
---|
520 | ENDIF |
---|
521 | |
---|
522 | IF(lwp) WRITE(numisl,REC=1) jpiglo, jpjglo, jpisl, zepsr |
---|
523 | DO jni = 1, jpisl |
---|
524 | CALL write2( numisl, bsfisl(1,1,jni), 1, jni+1 ) |
---|
525 | END DO |
---|
526 | nstop = nstop + 1 |
---|
527 | ENDIF |
---|
528 | |
---|
529 | END SUBROUTINE islbsf |
---|