1 |
C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_fgmres.F,v 1.16 2013/02/13 09:14:58 mlosch Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "SEAICE_OPTIONS.h" |
5 |
|
6 |
C-- File seaice_fgmres.F: seaice fgmres dynamical (linear) solver S/R: |
7 |
C-- Contents |
8 |
C-- o SEAICE_FGMRES_DRIVER |
9 |
C-- o SEAICE_MAP2VEC |
10 |
C-- o SEAICE_MAP_RS2VEC |
11 |
C-- o SEAICE_FGMRES |
12 |
C-- o SEAICE_SCALPROD |
13 |
|
14 |
CBOP |
15 |
C !ROUTINE: SEAICE_FGMRES_DRIVER |
16 |
C !INTERFACE: |
17 |
|
18 |
SUBROUTINE SEAICE_FGMRES_DRIVER( |
19 |
I uIceRes, vIceRes, |
20 |
U duIce, dvIce, |
21 |
U iCode, |
22 |
I FGMRESeps, iOutFGMRES, |
23 |
I newtonIter, |
24 |
U krylovIter, |
25 |
I myTime, myIter, myThid ) |
26 |
|
27 |
C !DESCRIPTION: \bv |
28 |
C *==========================================================* |
29 |
C | SUBROUTINE SEAICE_FGMRES_DRIVER |
30 |
C | o driver routine for fgmres |
31 |
C | o does the conversion between 2D fields and 1D vector |
32 |
C | back and forth |
33 |
C *==========================================================* |
34 |
C | written by Martin Losch, Oct 2012 |
35 |
C *==========================================================* |
36 |
C \ev |
37 |
|
38 |
C !USES: |
39 |
IMPLICIT NONE |
40 |
|
41 |
C === Global variables === |
42 |
#include "SIZE.h" |
43 |
#include "EEPARAMS.h" |
44 |
#include "SEAICE_SIZE.h" |
45 |
#include "SEAICE_PARAMS.h" |
46 |
|
47 |
C !INPUT/OUTPUT PARAMETERS: |
48 |
C === Routine arguments === |
49 |
C myTime :: Simulation time |
50 |
C myIter :: Simulation timestep number |
51 |
C myThid :: my Thread Id. number |
52 |
C newtonIter :: current iterate of Newton iteration (for diagnostics) |
53 |
C krylovIter :: current iterate of Newton iteration (updated) |
54 |
C iCode :: FGMRES parameter to determine next step |
55 |
C iOutFGMRES :: control output of fgmres |
56 |
_RL myTime |
57 |
INTEGER myIter |
58 |
INTEGER myThid |
59 |
INTEGER newtonIter |
60 |
INTEGER krylovIter |
61 |
INTEGER iOutFGMRES |
62 |
INTEGER iCode |
63 |
C FGMRESeps :: tolerance for FGMRES |
64 |
_RL FGMRESeps |
65 |
C du/vIce :: solution vector |
66 |
_RL duIce(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
67 |
_RL dvIce(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
68 |
C u/vIceRes :: residual F(u) |
69 |
_RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
70 |
_RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
71 |
|
72 |
#ifdef SEAICE_ALLOW_JFNK |
73 |
C Local variables: |
74 |
C k :: loop indices |
75 |
INTEGER k, bi, bj |
76 |
C FGMRES parameters |
77 |
C nVec :: size of the input vector(s) |
78 |
C im :: size of Krylov space |
79 |
C ifgmres :: interation counter |
80 |
INTEGER nVec |
81 |
PARAMETER ( nVec = 2*sNx*sNy ) |
82 |
INTEGER im |
83 |
PARAMETER ( im = 50 ) |
84 |
INTEGER ifgmres |
85 |
C work arrays |
86 |
_RL rhs(nVec,nSx,nSy), sol(nVec,nSx,nSy) |
87 |
_RL vv(nVec,im+1,nSx,nSy), w(nVec,im,nSx,nSy) |
88 |
_RL wk1(nVec,nSx,nSy), wk2(nVec,nSx,nSy) |
89 |
C need to store some of the fgmres parameters and fields so that |
90 |
C they are not forgotten between Krylov iterations |
91 |
COMMON /FGMRES_I/ ifgmres |
92 |
COMMON /FGMRES_RL/ sol, rhs, vv, w |
93 |
CEOP |
94 |
|
95 |
IF ( iCode .EQ. 0 ) THEN |
96 |
C The first guess is zero because it is a correction, but this |
97 |
C is implemented by setting du/vIce=0 outside of this routine; |
98 |
C this make it possible to restart FGMRES with a nonzero sol |
99 |
CALL SEAICE_MAP2VEC(nVec,duIce,dvIce,sol,.TRUE.,myThid) |
100 |
C wk2 needs to be reset for iCode = 0, because it may contain |
101 |
C remains of the previous Krylov iteration |
102 |
DO bj=myByLo(myThid),myByHi(myThid) |
103 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
104 |
DO k=1,nVec |
105 |
wk2(k,bi,bj) = 0. _d 0 |
106 |
ENDDO |
107 |
ENDDO |
108 |
ENDDO |
109 |
ELSEIF ( iCode .EQ. 3 ) THEN |
110 |
CALL SEAICE_MAP2VEC(nVec,uIceRes,vIceRes,rhs,.TRUE.,myThid) |
111 |
C change sign of rhs because we are solving J*u = -F |
112 |
C wk2 needs to be initialised for iCode = 3, because it may contain |
113 |
C garbage |
114 |
DO bj=myByLo(myThid),myByHi(myThid) |
115 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
116 |
DO k=1,nVec |
117 |
rhs(k,bi,bj) = -rhs(k,bi,bj) |
118 |
wk2(k,bi,bj) = 0. _d 0 |
119 |
ENDDO |
120 |
ENDDO |
121 |
ENDDO |
122 |
ELSE |
123 |
C map preconditioner results or Jacobian times vector, |
124 |
C stored in du/vIce to wk2 |
125 |
CALL SEAICE_MAP2VEC(nVec,duIce,dvIce,wk2,.TRUE.,myThid) |
126 |
ENDIF |
127 |
C |
128 |
CALL SEAICE_FGMRES (nVec,im,rhs,sol,ifgmres,krylovIter, |
129 |
U vv,w,wk1,wk2, |
130 |
I FGMRESeps,SEAICEkrylovIterMax,iOutFGMRES, |
131 |
U iCode, |
132 |
I myThid) |
133 |
C |
134 |
IF ( iCode .EQ. 0 ) THEN |
135 |
C map sol(ution) vector to du/vIce |
136 |
CALL SEAICE_MAP2VEC(nVec,duIce,dvIce,sol,.FALSE.,myThid) |
137 |
ELSE |
138 |
C map work vector to du/vIce to either compute a preconditioner |
139 |
C solution (wk1=rhs) or a Jacobian times wk1 |
140 |
CALL SEAICE_MAP2VEC(nVec,duIce,dvIce,wk1,.FALSE.,myThid) |
141 |
ENDIF |
142 |
|
143 |
C Fill overlaps in updated fields |
144 |
CALL EXCH_UV_XY_RL( duIce, dvIce,.TRUE.,myThid) |
145 |
|
146 |
RETURN |
147 |
END |
148 |
|
149 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
150 |
CBOP |
151 |
C !ROUTINE: SEAICE_MAP2VEC |
152 |
C !INTERFACE: |
153 |
|
154 |
SUBROUTINE SEAICE_MAP2VEC( |
155 |
I n, |
156 |
O xfld2d, yfld2d, |
157 |
U vector, |
158 |
I map2vec, myThid ) |
159 |
|
160 |
C !DESCRIPTION: \bv |
161 |
C *==========================================================* |
162 |
C | SUBROUTINE SEAICE_MAP2VEC |
163 |
C | o maps 2 2D-fields to vector and back |
164 |
C *==========================================================* |
165 |
C | written by Martin Losch, Oct 2012 |
166 |
C *==========================================================* |
167 |
C \ev |
168 |
|
169 |
C !USES: |
170 |
IMPLICIT NONE |
171 |
|
172 |
C === Global variables === |
173 |
#include "SIZE.h" |
174 |
#include "EEPARAMS.h" |
175 |
C === Routine arguments === |
176 |
INTEGER n |
177 |
LOGICAL map2vec |
178 |
INTEGER myThid |
179 |
_RL xfld2d (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
180 |
_RL yfld2d (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
181 |
_RL vector (n,nSx,nSy) |
182 |
C === local variables === |
183 |
INTEGER I, J, bi, bj |
184 |
INTEGER ii, jj, m |
185 |
CEOP |
186 |
|
187 |
m = n/2 |
188 |
DO bj=myByLo(myThid),myByHi(myThid) |
189 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
190 |
#ifdef SEAICE_JFNK_MAP_REORDER |
191 |
ii = 0 |
192 |
IF ( map2vec ) THEN |
193 |
DO J=1,sNy |
194 |
jj = 2*sNx*(J-1) |
195 |
DO I=1,sNx |
196 |
ii = jj + 2*I |
197 |
vector(ii-1,bi,bj) = xfld2d(I,J,bi,bj) |
198 |
vector(ii, bi,bj) = yfld2d(I,J,bi,bj) |
199 |
ENDDO |
200 |
ENDDO |
201 |
ELSE |
202 |
DO J=1,sNy |
203 |
jj = 2*sNx*(J-1) |
204 |
DO I=1,sNx |
205 |
ii = jj + 2*I |
206 |
xfld2d(I,J,bi,bj) = vector(ii-1,bi,bj) |
207 |
yfld2d(I,J,bi,bj) = vector(ii, bi,bj) |
208 |
ENDDO |
209 |
ENDDO |
210 |
ENDIF |
211 |
#else |
212 |
IF ( map2vec ) THEN |
213 |
DO J=1,sNy |
214 |
jj = sNx*(J-1) |
215 |
DO I=1,sNx |
216 |
ii = jj + I |
217 |
vector(ii, bi,bj) = xfld2d(I,J,bi,bj) |
218 |
vector(ii+m,bi,bj) = yfld2d(I,J,bi,bj) |
219 |
ENDDO |
220 |
ENDDO |
221 |
ELSE |
222 |
DO J=1,sNy |
223 |
jj = sNx*(J-1) |
224 |
DO I=1,sNx |
225 |
ii = jj + I |
226 |
xfld2d(I,J,bi,bj) = vector(ii, bi,bj) |
227 |
yfld2d(I,J,bi,bj) = vector(ii+m,bi,bj) |
228 |
ENDDO |
229 |
ENDDO |
230 |
ENDIF |
231 |
#endif /* SEAICE_JFNK_MAP_REORDER */ |
232 |
C bi,bj-loops |
233 |
ENDDO |
234 |
ENDDO |
235 |
|
236 |
RETURN |
237 |
END |
238 |
|
239 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
240 |
CBOP |
241 |
C !ROUTINE: SEAICE_MAP_RS2VEC |
242 |
C !INTERFACE: |
243 |
|
244 |
SUBROUTINE SEAICE_MAP_RS2VEC( |
245 |
I n, |
246 |
O xfld2d, yfld2d, |
247 |
U vector, |
248 |
I map2vec, myThid ) |
249 |
|
250 |
C !DESCRIPTION: \bv |
251 |
C *==========================================================* |
252 |
C | SUBROUTINE SEAICE_MAP_RS2VEC |
253 |
C | o maps 2 2D-RS-fields to vector and back |
254 |
C *==========================================================* |
255 |
C | written by Martin Losch, Oct 2012 |
256 |
C *==========================================================* |
257 |
C \ev |
258 |
|
259 |
C !USES: |
260 |
IMPLICIT NONE |
261 |
|
262 |
C === Global variables === |
263 |
#include "SIZE.h" |
264 |
#include "EEPARAMS.h" |
265 |
C === Routine arguments === |
266 |
INTEGER n |
267 |
LOGICAL map2vec |
268 |
INTEGER myThid |
269 |
_RS xfld2d (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
270 |
_RS yfld2d (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
271 |
_RL vector (n,nSx,nSy) |
272 |
C === local variables === |
273 |
INTEGER I, J, bi, bj |
274 |
INTEGER ii, jj, m |
275 |
CEOP |
276 |
|
277 |
m = n/2 |
278 |
DO bj=myByLo(myThid),myByHi(myThid) |
279 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
280 |
#ifdef SEAICE_JFNK_MAP_REORDER |
281 |
ii = 0 |
282 |
IF ( map2vec ) THEN |
283 |
DO J=1,sNy |
284 |
jj = 2*sNx*(J-1) |
285 |
DO I=1,sNx |
286 |
ii = jj + 2*I |
287 |
vector(ii-1,bi,bj) = xfld2d(I,J,bi,bj) |
288 |
vector(ii, bi,bj) = yfld2d(I,J,bi,bj) |
289 |
ENDDO |
290 |
ENDDO |
291 |
ELSE |
292 |
DO J=1,sNy |
293 |
jj = 2*sNx*(J-1) |
294 |
DO I=1,sNx |
295 |
ii = jj + 2*I |
296 |
xfld2d(I,J,bi,bj) = vector(ii-1,bi,bj) |
297 |
yfld2d(I,J,bi,bj) = vector(ii, bi,bj) |
298 |
ENDDO |
299 |
ENDDO |
300 |
ENDIF |
301 |
#else |
302 |
IF ( map2vec ) THEN |
303 |
DO J=1,sNy |
304 |
jj = sNx*(J-1) |
305 |
DO I=1,sNx |
306 |
ii = jj + I |
307 |
vector(ii, bi,bj) = xfld2d(I,J,bi,bj) |
308 |
vector(ii+m,bi,bj) = yfld2d(I,J,bi,bj) |
309 |
ENDDO |
310 |
ENDDO |
311 |
ELSE |
312 |
DO J=1,sNy |
313 |
jj = sNx*(J-1) |
314 |
DO I=1,sNx |
315 |
ii = jj + I |
316 |
xfld2d(I,J,bi,bj) = vector(ii, bi,bj) |
317 |
yfld2d(I,J,bi,bj) = vector(ii+m,bi,bj) |
318 |
ENDDO |
319 |
ENDDO |
320 |
ENDIF |
321 |
#endif /* SEAICE_JFNK_MAP_REORDER */ |
322 |
C bi,bj-loops |
323 |
ENDDO |
324 |
ENDDO |
325 |
|
326 |
RETURN |
327 |
END |
328 |
|
329 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
330 |
CBOP |
331 |
C !ROUTINE: SEAICE_FGMRES |
332 |
C !INTERFACE: |
333 |
SUBROUTINE SEAICE_FGMRES ( |
334 |
I n,im,rhs, |
335 |
U sol,i,its,vv,w,wk1,wk2, |
336 |
I eps,maxits,iout, |
337 |
U icode, |
338 |
I myThid ) |
339 |
|
340 |
C----------------------------------------------------------------------- |
341 |
C mlosch Oct 2012: modified the routine further to be compliant with |
342 |
C MITgcm standards: |
343 |
C f90 -> F |
344 |
C !-comment -> C-comment |
345 |
C add its to list of arguments |
346 |
C double precision -> _RL |
347 |
C implicit none |
348 |
C |
349 |
C jfl Dec 1st 2006. We modified the routine so that it is double precison. |
350 |
C Here are the modifications: |
351 |
C 1) implicit real (a-h,o-z) becomes implicit real*8 (a-h,o-z) |
352 |
C 2) real bocomes real*8 |
353 |
C 3) subroutine scopy.f has been changed for dcopy.f |
354 |
C 4) subroutine saxpy.f has been changed for daxpy.f |
355 |
C 5) function sdot.f has been changed for ddot.f |
356 |
C 6) 1e-08 becomes 1d-08 |
357 |
C |
358 |
C Be careful with the dcopy, daxpy and ddot code...there is a slight |
359 |
C difference with the single precision versions (scopy, saxpy and sdot). |
360 |
C In the single precision versions, the array are declared sightly differently. |
361 |
C It is written for single precision: |
362 |
C |
363 |
C modified 12/3/93, array(1) declarations changed to array(*) |
364 |
C----------------------------------------------------------------------- |
365 |
|
366 |
implicit none |
367 |
C === Global variables === |
368 |
#include "SIZE.h" |
369 |
#include "EEPARAMS.h" |
370 |
CML implicit double precision (a-h,o-z) !jfl modification |
371 |
integer myThid |
372 |
integer n, im, its, maxits, iout, icode |
373 |
_RL rhs(n,nSx,nSy), sol(n,nSx,nSy) |
374 |
_RL vv(n,im+1,nSx,nSy), w(n,im,nSx,nSy) |
375 |
_RL wk1(n,nSx,nSy), wk2(n,nSx,nSy), eps |
376 |
C----------------------------------------------------------------------- |
377 |
C flexible GMRES routine. This is a version of GMRES which allows a |
378 |
C a variable preconditioner. Implemented with a reverse communication |
379 |
C protocole for flexibility - |
380 |
C DISTRIBUTED VERSION (USES DISTDOT FOR DDOT) |
381 |
C explicit (exact) residual norms for restarts |
382 |
C written by Y. Saad, modified by A. Malevsky, version February 1, 1995 |
383 |
C----------------------------------------------------------------------- |
384 |
C This Is A Reverse Communication Implementation. |
385 |
C------------------------------------------------- |
386 |
C USAGE: (see also comments for icode below). FGMRES |
387 |
C should be put in a loop and the loop should be active for as |
388 |
C long as icode is not equal to 0. On return fgmres will |
389 |
C 1) either be requesting the new preconditioned vector applied |
390 |
C to wk1 in case icode.eq.1 (result should be put in wk2) |
391 |
C 2) or be requesting the product of A applied to the vector wk1 |
392 |
C in case icode.eq.2 (result should be put in wk2) |
393 |
C 3) or be terminated in case icode .eq. 0. |
394 |
C on entry always set icode = 0. So icode should be set back to zero |
395 |
C upon convergence. |
396 |
C----------------------------------------------------------------------- |
397 |
C Here is a typical way of running fgmres: |
398 |
C |
399 |
C icode = 0 |
400 |
C 1 continue |
401 |
C call fgmres (n,im,rhs,sol,i,vv,w,wk1,wk2,eps,maxits,iout, |
402 |
C & icode,its,mythid) |
403 |
C |
404 |
C if (icode .eq. 1) then |
405 |
C call precon(n, wk1, wk2) <--- user variable preconditioning |
406 |
C goto 1 |
407 |
C else if (icode .ge. 2) then |
408 |
C call matvec (n,wk1, wk2) <--- user matrix vector product. |
409 |
C goto 1 |
410 |
C else |
411 |
C ----- done ---- |
412 |
C ......... |
413 |
C----------------------------------------------------------------------- |
414 |
C list of parameters |
415 |
C------------------- |
416 |
C |
417 |
C n == integer. the dimension of the problem |
418 |
C im == size of Krylov subspace: should not exceed 50 in this |
419 |
C version (can be reset in code. looking at comment below) |
420 |
C rhs == vector of length n containing the right hand side |
421 |
C sol == initial guess on input, approximate solution on output |
422 |
C vv == work space of size n x (im+1) |
423 |
C w == work space of length n x im |
424 |
C wk1, |
425 |
C wk2, == two work vectors of length n each used for the reverse |
426 |
C communication protocole. When on return (icode .ne. 1) |
427 |
C the user should call fgmres again with wk2 = precon * wk1 |
428 |
C and icode untouched. When icode.eq.1 then it means that |
429 |
C convergence has taken place. |
430 |
C |
431 |
C eps == tolerance for stopping criterion. process is stopped |
432 |
C as soon as ( ||.|| is the euclidean norm): |
433 |
C || current residual||/||initial residual|| <= eps |
434 |
C |
435 |
C maxits== maximum number of iterations allowed |
436 |
C |
437 |
C i == internal iteration counter, updated in this routine |
438 |
C its == current (Krylov) iteration counter, updated in this routine |
439 |
C |
440 |
C iout == output unit number number for printing intermediate results |
441 |
C if (iout .le. 0) no statistics are printed. |
442 |
C |
443 |
C icode = integer. indicator for the reverse communication protocole. |
444 |
C ON ENTRY : icode should be set to icode = 0. |
445 |
C ON RETURN: |
446 |
C * icode .eq. 1 value means that fgmres has not finished |
447 |
C and that it is requesting a preconditioned vector before |
448 |
C continuing. The user must compute M**(-1) wk1, where M is |
449 |
C the preconditioing matrix (may vary at each call) and wk1 is |
450 |
C the vector as provided by fgmres upun return, and put the |
451 |
C result in wk2. Then fgmres must be called again without |
452 |
C changing any other argument. |
453 |
C * icode .eq. 2 value means that fgmres has not finished |
454 |
C and that it is requesting a matrix vector product before |
455 |
C continuing. The user must compute A * wk1, where A is the |
456 |
C coefficient matrix and wk1 is the vector provided by |
457 |
C upon return. The result of the operation is to be put in |
458 |
C the vector wk2. Then fgmres must be called again without |
459 |
C changing any other argument. |
460 |
C * icode .eq. 0 means that fgmres has finished and sol contains |
461 |
C the approximate solution. |
462 |
C comment: typically fgmres must be implemented in a loop |
463 |
C with fgmres being called as long icode is returned with |
464 |
C a value .ne. 0. |
465 |
C----------------------------------------------------------------------- |
466 |
C local variables -- !jfl modif |
467 |
integer imax |
468 |
parameter ( imax = 50 ) |
469 |
_RL hh(4*imax+1,4*imax),c(4*imax),s(4*imax) |
470 |
_RL rs(4*imax+1),t,ro |
471 |
C------------------------------------------------------------- |
472 |
C arnoldi size should not exceed 50 in this version.. |
473 |
C------------------------------------------------------------- |
474 |
integer i, i1, ii, j, jj, k, k1!, n1 |
475 |
integer bi, bj |
476 |
_RL r0, gam, epsmac, eps1 |
477 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
478 |
|
479 |
CEOP |
480 |
CML save |
481 |
C local common block to replace the save statement |
482 |
COMMON /SEAICE_FMRES_LOC_I/ i1 |
483 |
COMMON /SEAICE_FMRES_LOC_RL/ |
484 |
& hh, c, s, rs, t, ro, r0, gam, epsmac, eps1 |
485 |
data epsmac/1.d-16/ |
486 |
C |
487 |
C computed goto |
488 |
C |
489 |
if ( im .gt. imax ) stop 'size of krylov space > 50' |
490 |
goto (100,200,300,11) icode +1 |
491 |
100 continue |
492 |
CML n1 = n + 1 |
493 |
its = 0 |
494 |
C------------------------------------------------------------- |
495 |
C ** outer loop starts here.. |
496 |
C--------------compute initial residual vector -------------- |
497 |
C 10 continue |
498 |
CML call dcopy (n, sol, 1, wk1, 1) !jfl modification |
499 |
do bj=myByLo(myThid),myByHi(myThid) |
500 |
do bi=myBxLo(myThid),myBxHi(myThid) |
501 |
do j=1,n |
502 |
wk1(j,bi,bj)=sol(j,bi,bj) |
503 |
enddo |
504 |
enddo |
505 |
enddo |
506 |
icode = 3 |
507 |
RETURN |
508 |
11 continue |
509 |
do bj=myByLo(myThid),myByHi(myThid) |
510 |
do bi=myBxLo(myThid),myBxHi(myThid) |
511 |
do j=1,n |
512 |
vv(j,1,bi,bj) = rhs(j,bi,bj) - wk2(j,bi,bj) |
513 |
enddo |
514 |
enddo |
515 |
enddo |
516 |
20 continue |
517 |
CML ro = ddot(n, vv, 1, vv,1) !jfl modification |
518 |
call SEAICE_SCALPROD(n, im+1, 1, 1, vv, vv, ro, myThid) |
519 |
ro = sqrt(ro) |
520 |
if (ro .eq. 0.0 _d 0) goto 999 |
521 |
t = 1.0 _d 0/ ro |
522 |
do bj=myByLo(myThid),myByHi(myThid) |
523 |
do bi=myBxLo(myThid),myBxHi(myThid) |
524 |
do j=1, n |
525 |
vv(j,1,bi,bj) = vv(j,1,bi,bj)*t |
526 |
enddo |
527 |
enddo |
528 |
enddo |
529 |
if (its .eq. 0) eps1=eps |
530 |
C not sure what this is, r0 is never used again |
531 |
if (its .eq. 0) r0 = ro |
532 |
if (iout .gt. 0) then |
533 |
_BEGIN_MASTER( myThid ) |
534 |
write(msgBuf, 199) its, ro |
535 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
536 |
& SQUEEZE_RIGHT, myThid ) |
537 |
C print *,'chau',its, ro !write(iout, 199) its, ro |
538 |
_END_MASTER( myThid ) |
539 |
endif |
540 |
C |
541 |
C initialize 1-st term of rhs of hessenberg system.. |
542 |
C |
543 |
rs(1) = ro |
544 |
i = 0 |
545 |
4 continue |
546 |
i=i+1 |
547 |
its = its + 1 |
548 |
i1 = i + 1 |
549 |
do bj=myByLo(myThid),myByHi(myThid) |
550 |
do bi=myBxLo(myThid),myBxHi(myThid) |
551 |
do k=1, n |
552 |
wk1(k,bi,bj) = vv(k,i,bi,bj) |
553 |
enddo |
554 |
enddo |
555 |
enddo |
556 |
C |
557 |
C return |
558 |
C |
559 |
icode = 1 |
560 |
RETURN |
561 |
200 continue |
562 |
do bj=myByLo(myThid),myByHi(myThid) |
563 |
do bi=myBxLo(myThid),myBxHi(myThid) |
564 |
do k=1, n |
565 |
w(k,i,bi,bj) = wk2(k,bi,bj) |
566 |
enddo |
567 |
enddo |
568 |
enddo |
569 |
C |
570 |
C call matvec operation |
571 |
C |
572 |
CML call dcopy(n, wk2, 1, wk1, 1) !jfl modification |
573 |
do bj=myByLo(myThid),myByHi(myThid) |
574 |
do bi=myBxLo(myThid),myBxHi(myThid) |
575 |
do k=1,n |
576 |
wk1(k,bi,bj)=wk2(k,bi,bj) |
577 |
enddo |
578 |
enddo |
579 |
enddo |
580 |
C |
581 |
C return |
582 |
C |
583 |
icode = 2 |
584 |
RETURN |
585 |
300 continue |
586 |
C |
587 |
C first call to ope corresponds to intialization goto back to 11. |
588 |
C |
589 |
C if (icode .eq. 3) goto 11 |
590 |
CML call dcopy (n, wk2, 1, vv(1,i1), 1) !jfl modification |
591 |
do bj=myByLo(myThid),myByHi(myThid) |
592 |
do bi=myBxLo(myThid),myBxHi(myThid) |
593 |
do k=1,n |
594 |
vv(k,i1,bi,bj)=wk2(k,bi,bj) |
595 |
enddo |
596 |
enddo |
597 |
enddo |
598 |
C |
599 |
C modified gram - schmidt... |
600 |
C |
601 |
do j=1, i |
602 |
CML t = ddot(n, vv(1,j), 1, vv(1,i1), 1) !jfl modification |
603 |
call SEAICE_SCALPROD(n, im+1, j, i1, vv, vv, t, myThid) |
604 |
hh(j,i) = t |
605 |
CML call daxpy(n, -t, vv(1,j), 1, vv(1,i1), 1) !jfl modification |
606 |
CML enddo |
607 |
CML do j=1, i |
608 |
CML t = hh(j,i) |
609 |
do bj=myByLo(myThid),myByHi(myThid) |
610 |
do bi=myBxLo(myThid),myBxHi(myThid) |
611 |
do k=1,n |
612 |
vv(k,i1,bi,bj) = vv(k,i1,bi,bj) - t*vv(k,j,bi,bj) |
613 |
enddo |
614 |
enddo |
615 |
enddo |
616 |
enddo |
617 |
CML t = sqrt(ddot(n, vv(1,i1), 1, vv(1,i1), 1)) !jfl modification |
618 |
call SEAICE_SCALPROD(n, im+1, i1, i1, vv, vv, t, myThid) |
619 |
t = sqrt(t) |
620 |
hh(i1,i) = t |
621 |
if (t .ne. 0.0 _d 0) then |
622 |
t = 1.0 _d 0 / t |
623 |
do bj=myByLo(myThid),myByHi(myThid) |
624 |
do bi=myBxLo(myThid),myBxHi(myThid) |
625 |
do k=1,n |
626 |
vv(k,i1,bi,bj) = vv(k,i1,bi,bj)*t |
627 |
enddo |
628 |
enddo |
629 |
enddo |
630 |
endif |
631 |
C |
632 |
C done with modified gram schimd and arnoldi step. |
633 |
C now update factorization of hh |
634 |
C |
635 |
if (i .ne. 1) then |
636 |
C |
637 |
C perfrom previous transformations on i-th column of h |
638 |
C |
639 |
do k=2,i |
640 |
k1 = k-1 |
641 |
t = hh(k1,i) |
642 |
hh(k1,i) = c(k1)*t + s(k1)*hh(k,i) |
643 |
hh(k,i) = -s(k1)*t + c(k1)*hh(k,i) |
644 |
enddo |
645 |
endif |
646 |
gam = sqrt(hh(i,i)**2 + hh(i1,i)**2) |
647 |
if (gam .eq. 0.0 _d 0) gam = epsmac |
648 |
C-----------#determine next plane rotation #------------------- |
649 |
c(i) = hh(i,i)/gam |
650 |
s(i) = hh(i1,i)/gam |
651 |
C numerically more stable Givens rotation, but the results |
652 |
C are not better |
653 |
CML c(i)=1. _d 0 |
654 |
CML s(i)=0. _d 0 |
655 |
CML if ( abs(hh(i1,i)) .gt. 0.0 _d 0) then |
656 |
CML if ( abs(hh(i1,i)) .gt. abs(hh(i,i)) ) then |
657 |
CML gam = hh(i,i)/hh(i1,i) |
658 |
CML s(i) = 1./sqrt(1.+gam*gam) |
659 |
CML c(i) = s(i)*gam |
660 |
CML else |
661 |
CML gam = hh(i1,i)/hh(i,i) |
662 |
CML c(i) = 1./sqrt(1.+gam*gam) |
663 |
CML s(i) = c(i)*gam |
664 |
CML endif |
665 |
CML endif |
666 |
rs(i1) = -s(i)*rs(i) |
667 |
rs(i) = c(i)*rs(i) |
668 |
C |
669 |
C determine res. norm. and test for convergence |
670 |
C |
671 |
hh(i,i) = c(i)*hh(i,i) + s(i)*hh(i1,i) |
672 |
ro = abs(rs(i1)) |
673 |
if (iout .gt. 0) then |
674 |
_BEGIN_MASTER( myThid ) |
675 |
write(msgBuf, 199) its, ro |
676 |
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
677 |
& SQUEEZE_RIGHT, myThid ) |
678 |
_END_MASTER( myThid ) |
679 |
endif |
680 |
if (i .lt. im .and. (ro .gt. eps1)) goto 4 |
681 |
C |
682 |
C now compute solution. first solve upper triangular system. |
683 |
C |
684 |
rs(i) = rs(i)/hh(i,i) |
685 |
do ii=2,i |
686 |
k=i-ii+1 |
687 |
k1 = k+1 |
688 |
t=rs(k) |
689 |
do j=k1,i |
690 |
t = t-hh(k,j)*rs(j) |
691 |
enddo |
692 |
rs(k) = t/hh(k,k) |
693 |
enddo |
694 |
C |
695 |
C done with back substitution.. |
696 |
C now form linear combination to get solution |
697 |
C |
698 |
do j=1, i |
699 |
t = rs(j) |
700 |
CML call daxpy(n, t, w(1,j), 1, sol,1) !jfl modification |
701 |
do bj=myByLo(myThid),myByHi(myThid) |
702 |
do bi=myBxLo(myThid),myBxHi(myThid) |
703 |
do k=1,n |
704 |
sol(k,bi,bj) = sol(k,bi,bj) + t*w(k,j,bi,bj) |
705 |
enddo |
706 |
enddo |
707 |
enddo |
708 |
enddo |
709 |
C |
710 |
C test for return |
711 |
C |
712 |
if (ro .le. eps1 .or. its .ge. maxits) goto 999 |
713 |
C |
714 |
C else compute residual vector and continue.. |
715 |
C |
716 |
C goto 10 |
717 |
|
718 |
do j=1,i |
719 |
jj = i1-j+1 |
720 |
rs(jj-1) = -s(jj-1)*rs(jj) |
721 |
rs(jj) = c(jj-1)*rs(jj) |
722 |
enddo |
723 |
do j=1,i1 |
724 |
t = rs(j) |
725 |
if (j .eq. 1) t = t-1.0 _d 0 |
726 |
CML call daxpy (n, t, vv(1,j), 1, vv, 1) |
727 |
do bj=myByLo(myThid),myByHi(myThid) |
728 |
do bi=myBxLo(myThid),myBxHi(myThid) |
729 |
do k=1,n |
730 |
vv(k,1,bi,bj) = vv(k,1,bi,bj) + t*vv(k,j,bi,bj) |
731 |
enddo |
732 |
enddo |
733 |
enddo |
734 |
enddo |
735 |
C |
736 |
C restart outer loop. |
737 |
C |
738 |
goto 20 |
739 |
999 icode = 0 |
740 |
|
741 |
199 format(' SEAICE_FGMRES: its =', i4, ' res. norm =', d26.16) |
742 |
C |
743 |
RETURN |
744 |
C-----end-of-fgmres----------------------------------------------------- |
745 |
C----------------------------------------------------------------------- |
746 |
END |
747 |
|
748 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
749 |
CBOP |
750 |
C !ROUTINE: SEAICE_SCALPROD |
751 |
C !INTERFACE: |
752 |
|
753 |
subroutine SEAICE_SCALPROD(n,im,i1,i2,dx,dy,t,myThid) |
754 |
|
755 |
C forms the dot product of two vectors. |
756 |
C uses unrolled loops for increments equal to one. |
757 |
C jack dongarra, linpack, 3/11/78. |
758 |
C ML: code stolen from BLAS-ddot and adapted for parallel applications |
759 |
|
760 |
implicit none |
761 |
#include "SIZE.h" |
762 |
#include "EEPARAMS.h" |
763 |
#include "EESUPPORT.h" |
764 |
#include "SEAICE_SIZE.h" |
765 |
#include "SEAICE.h" |
766 |
integer n, im, i1, i2 |
767 |
_RL dx(n,im,nSx,nSy),dy(n,im,nSx,nSy) |
768 |
_RL t |
769 |
integer myThid |
770 |
C local arrays |
771 |
_RL dtemp(nSx,nSy) |
772 |
integer i,m,mp1,bi,bj |
773 |
CEOP |
774 |
|
775 |
m = mod(n,5) |
776 |
mp1 = m + 1 |
777 |
t = 0. _d 0 |
778 |
c if( m .eq. 0 ) go to 40 |
779 |
do bj=myByLo(myThid),myByHi(myThid) |
780 |
do bi=myBxLo(myThid),myBxHi(myThid) |
781 |
dtemp(bi,bj) = 0. _d 0 |
782 |
if ( m .ne. 0 ) then |
783 |
do i = 1,m |
784 |
dtemp(bi,bj) = dtemp(bi,bj) + dx(i,i1,bi,bj)*dy(i,i2,bi,bj) |
785 |
& * scalarProductMetric(i,1,bi,bj) |
786 |
enddo |
787 |
endif |
788 |
if ( n .ge. 5 ) then |
789 |
c if( n .lt. 5 ) go to 60 |
790 |
c40 mp1 = m + 1 |
791 |
do i = mp1,n,5 |
792 |
dtemp(bi,bj) = dtemp(bi,bj) + |
793 |
& dx(i, i1,bi,bj)*dy(i, i2,bi,bj) |
794 |
& * scalarProductMetric(i, 1, bi,bj) + |
795 |
& dx(i + 1,i1,bi,bj)*dy(i + 1,i2,bi,bj) |
796 |
& * scalarProductMetric(i + 1,1, bi,bj) + |
797 |
& dx(i + 2,i1,bi,bj)*dy(i + 2,i2,bi,bj) |
798 |
& * scalarProductMetric(i + 2,1, bi,bj) + |
799 |
& dx(i + 3,i1,bi,bj)*dy(i + 3,i2,bi,bj) |
800 |
& * scalarProductMetric(i + 3,1, bi,bj) + |
801 |
& dx(i + 4,i1,bi,bj)*dy(i + 4,i2,bi,bj) |
802 |
& * scalarProductMetric(i + 4,1, bi,bj) |
803 |
enddo |
804 |
c60 continue |
805 |
endif |
806 |
enddo |
807 |
enddo |
808 |
CALL GLOBAL_SUM_TILE_RL( dtemp,t,myThid ) |
809 |
|
810 |
#endif /* SEAICE_ALLOW_JFNK */ |
811 |
|
812 |
RETURN |
813 |
END |