/[MITgcm]/MITgcm_contrib/dgoldberg/streamice/streamice_cg_solve.F
ViewVC logotype

Contents of /MITgcm_contrib/dgoldberg/streamice/streamice_cg_solve.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.10 - (show annotations) (download)
Wed Aug 27 19:29:13 2014 UTC (10 years, 10 months ago) by dgoldberg
Branch: MAIN
CVS Tags: HEAD
Changes since 1.9: +55 -301 lines
updating contrib streamice repo with latest files, and separated out convergence checks; and parameterised maximum iteration counts and interface w shelfice for coupling

1 C $Header: /u/gcmpack/MITgcm/pkg/streamice/streamice_cg_solve.F,v 1.5 2014/07/09 16:09:23 dgoldberg Exp $
2 C $Name: $
3
4 #include "STREAMICE_OPTIONS.h"
5
6 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
7
8 CBOP
9 SUBROUTINE STREAMICE_CG_SOLVE(
10 U cg_Uin, ! x-velocities
11 U cg_Vin, ! y-velocities
12 I cg_Bu, ! force in x dir
13 I cg_Bv, ! force in y dir
14 I A_uu, ! section of matrix that multiplies u and projects on u
15 I A_uv, ! section of matrix that multiplies v and projects on u
16 I A_vu, ! section of matrix that multiplies u and projects on v
17 I A_vv, ! section of matrix that multiplies v and projects on v
18 I tolerance,
19 O iters,
20 I maxIter,
21 I myThid )
22 C /============================================================\
23 C | SUBROUTINE |
24 C | o |
25 C |============================================================|
26 C | |
27 C \============================================================/
28 IMPLICIT NONE
29
30 #include "SIZE.h"
31 #include "EEPARAMS.h"
32 #include "PARAMS.h"
33 #include "STREAMICE.h"
34 #include "STREAMICE_CG.h"
35
36
37
38 !#ifdef ALLOW_PETSC
39 !#include "finclude/petsc.h"
40 ! UNCOMMENT IF V3.0
41 !#include "finclude/petscvec.h"
42 !#include "finclude/petscmat.h"
43 !#include "finclude/petscksp.h"
44 !#include "finclude/petscpc.h"
45 !#endif
46 C === Global variables ===
47
48
49 C !INPUT/OUTPUT ARGUMENTS
50 C cg_Uin, cg_Vin - input and output velocities
51 C cg_Bu, cg_Bv - driving stress
52 INTEGER myThid
53 INTEGER iters
54 INTEGER maxIter
55 _RL tolerance
56 _RL cg_Uin (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
57 _RL cg_Vin (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
58 _RL cg_Bu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
59 _RL cg_Bv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
60 _RL
61 & A_uu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1),
62 & A_vu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1),
63 & A_uv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1),
64 & A_vv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
65
66 C LOCAL VARIABLES
67 INTEGER i, j, bi, bj, cg_halo, conv_flag
68 INTEGER iter, is, js, ie, je, colx, coly, k
69 _RL dot_p1, dot_p2, alpha_k, beta_k, resid, resid_0
70 _RL dot_p1_tile (nSx,nSy)
71 _RL dot_p2_tile (nSx,nSy)
72 CHARACTER*(MAX_LEN_MBUF) msgBuf
73
74
75 !#ifdef ALLOW_PETSC
76 ! INTEGER indices(2*(snx*nsx*sny*nsy))
77 ! INTEGER n_dofs_cum_sum (0:nPx*nPy-1), idx(1)
78 ! _RL rhs_values(2*(snx*nsx*sny*nsy))
79 ! _RL solution_values(2*(snx*nsx*sny*nsy))
80 ! _RL mat_values (2*Nx*Ny,2*(snx*nsx*sny*nsy))
81 ! _RL mat_values (18,1), mat_val_return(1)
82 ! INTEGER indices_col(18)
83 ! INTEGER local_dofs, global_dofs, dof_index, dof_index_col
84 ! INTEGER local_offset
85 ! Mat matrix
86 ! KSP ksp
87 ! PC pc
88 ! Vec rhs
89 ! Vec solution
90 ! PetscErrorCode ierr
91 !#ifdef ALLOW_USE_MPI
92 ! integer mpiRC, mpiMyWid
93 !#endif
94 !#endif
95
96
97 #ifdef ALLOW_STREAMICE
98
99
100
101 CALL TIMER_START ('STREAMICE_CG_SOLVE',myThid)
102 #ifndef STREAMICE_SERIAL_TRISOLVE
103
104 #ifdef ALLOW_PETSC
105
106 CALL STREAMICE_CG_SOLVE_PETSC(
107 U cg_Uin, ! x-velocities
108 U cg_Vin, ! y-velocities
109 I cg_Bu, ! force in x dir
110 I cg_Bv, ! force in y dir
111 I A_uu, ! section of matrix that multiplies u and projects on u
112 I A_uv, ! section of matrix that multiplies v and projects on u
113 I A_vu, ! section of matrix that multiplies u and projects on v
114 I A_vv, ! section of matrix that multiplies v and projects on v
115 I tolerance,
116 I maxIter,
117 O iters,
118 I myThid )
119
120
121 #else /* ALLOW_PETSC */
122
123
124 iters = maxIter
125 conv_flag = 0
126
127
128 DO bj = myByLo(myThid), myByHi(myThid)
129 DO bi = myBxLo(myThid), myBxHi(myThid)
130 DO j=1,sNy
131 DO i=1,sNx
132 Zu_SI (i,j,bi,bj) = 0. _d 0
133 Zv_SI (i,j,bi,bj) = 0. _d 0
134 Ru_SI (i,j,bi,bj) = 0. _d 0
135 Rv_SI (i,j,bi,bj) = 0. _d 0
136 Au_SI (i,j,bi,bj) = 0. _d 0
137 Av_SI (i,j,bi,bj) = 0. _d 0
138 Du_SI (i,j,bi,bj) = 0. _d 0
139 Dv_SI (i,j,bi,bj) = 0. _d 0
140 ENDDO
141 ENDDO
142 ENDDO
143 ENDDO
144
145 C FIND INITIAL RESIDUAL, and initialize r
146
147 ! #ifdef STREAMICE_CONSTRUCT_MATRIX
148
149
150
151 DO bj = myByLo(myThid), myByHi(myThid)
152 DO bi = myBxLo(myThid), myBxHi(myThid)
153 DO j=1,sNy
154 DO i=1,sNx
155 DO colx=-1,1
156 DO coly=-1,1
157 Au_SI(i,j,bi,bj) = Au_SI(i,j,bi,bj) +
158 & A_uu(i,j,bi,bj,colx,coly)*
159 & cg_Uin(i+colx,j+coly,bi,bj)+
160 & A_uv(i,j,bi,bj,colx,coly)*
161 & cg_Vin(i+colx,j+coly,bi,bj)
162
163
164 Av_SI(i,j,bi,bj) = Av_SI(i,j,bi,bj) +
165 & A_vu(i,j,bi,bj,colx,coly)*
166 & cg_Uin(i+colx,j+coly,bi,bj)+
167 & A_vv(i,j,bi,bj,colx,coly)*
168 & cg_Vin(i+colx,j+coly,bi,bj)
169 ENDDO
170 ENDDO
171 ENDDO
172 ENDDO
173 ENDDO
174 ENDDO
175
176
177 _EXCH_XY_RL( Au_SI, myThid )
178 _EXCH_XY_RL( Av_SI, myThid )
179
180
181 DO bj = myByLo(myThid), myByHi(myThid)
182 DO bi = myBxLo(myThid), myBxHi(myThid)
183 DO j=1-OLy,sNy+OLy
184 DO i=1-OLx,sNx+OLx
185 Ru_SI(i,j,bi,bj)=cg_Bu(i,j,bi,bj)-
186 & Au_SI(i,j,bi,bj)
187 Rv_SI(i,j,bi,bj)=cg_Bv(i,j,bi,bj)-
188 & Av_SI(i,j,bi,bj)
189 ENDDO
190 ENDDO
191 dot_p1_tile(bi,bj) = 0. _d 0
192 dot_p2_tile(bi,bj) = 0. _d 0
193 ENDDO
194 ENDDO
195
196 DO bj = myByLo(myThid), myByHi(myThid)
197 DO bi = myBxLo(myThid), myBxHi(myThid)
198 DO j=1,sNy
199 DO i=1,sNx
200 IF (STREAMICE_umask(i,j,bi,bj).eq.1.0)
201 & dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Ru_SI(i,j,bi,bj)**2
202 IF (STREAMICE_vmask(i,j,bi,bj).eq.1.0)
203 & dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Rv_SI(i,j,bi,bj)**2
204 ENDDO
205 ENDDO
206 ENDDO
207 ENDDO
208
209
210 CALL GLOBAL_SUM_TILE_RL( dot_p1_tile, dot_p1, myThid )
211 resid_0 = sqrt(dot_p1)
212
213 DO bj = myByLo(myThid), myByHi(myThid)
214 DO bi = myBxLo(myThid), myBxHi(myThid)
215
216 WRITE(msgBuf,'(A,I1,I1,E14.7)') 'CONJ GRAD INIT RESID LOCAL, ',
217 & bi,bj, dot_p1_tile(bi,bj)
218 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
219 & SQUEEZE_RIGHT , 1)
220
221 enddo
222 enddo
223
224 WRITE(msgBuf,'(A,E14.7)') 'CONJ GRAD INIT RESID, ',
225 & resid_0
226 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
227 & SQUEEZE_RIGHT , 1)
228
229 C CCCCCCCCCCCCCCCCCCCC
230
231 DO bj = myByLo(myThid), myByHi(myThid)
232 DO bi = myBxLo(myThid), myBxHi(myThid)
233 DO j=1-OLy,sNy+OLy
234 DO i=1-OLx,sNx+OLx
235 IF (STREAMICE_umask(i,j,bi,bj).eq.1.0)
236 & Zu_SI(i,j,bi,bj)=Ru_SI(i,j,bi,bj) / DIAGu_SI(i,j,bi,bj)
237 IF (STREAMICE_vmask(i,j,bi,bj).eq.1.0)
238 & Zv_SI(i,j,bi,bj)=Rv_SI(i,j,bi,bj) / DIAGv_SI(i,j,bi,bj)
239 ENDDO
240 ENDDO
241 ENDDO
242 ENDDO
243
244 cg_halo = min(OLx-1,OLy-1)
245 conv_flag = 0
246
247 DO bj = myByLo(myThid), myByHi(myThid)
248 DO bi = myBxLo(myThid), myBxHi(myThid)
249 DO j=1-OLy,sNy+OLy
250 DO i=1-OLx,sNx+OLx
251 Du_SI(i,j,bi,bj)=Zu_SI(i,j,bi,bj)
252 Dv_SI(i,j,bi,bj)=Zv_SI(i,j,bi,bj)
253 ENDDO
254 ENDDO
255 ENDDO
256 ENDDO
257
258 resid = resid_0
259 iters = 0
260
261 c !!!!!!!!!!!!!!!!!!
262 c !! !!
263 c !! MAIN CG LOOP !!
264 c !! !!
265 c !!!!!!!!!!!!!!!!!!
266
267
268
269
270 c ! initially, b-grid data is valid up to 3 halo nodes out -- right? (check for MITgcm!!)
271
272 WRITE(msgBuf,'(A)') 'BEGINNING MAIN CG LOOP'
273 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
274 & SQUEEZE_RIGHT , 1)
275
276 ! IF(STREAMICE_construct_matrix) CALL STREAMICE_CG_MAKE_A(myThid)
277
278
279 do iter = 1, maxIter
280 if (resid .gt. tolerance*resid_0) then
281
282 c to avoid using "exit"
283 iters = iters + 1
284
285 is = 1 - cg_halo
286 ie = sNx + cg_halo
287 js = 1 - cg_halo
288 je = sNy + cg_halo
289
290 DO bj = myByLo(myThid), myByHi(myThid)
291 DO bi = myBxLo(myThid), myBxHi(myThid)
292 DO j=1-OLy,sNy+OLy
293 DO i=1-OLx,sNx+OLx
294 Au_SI(i,j,bi,bj) = 0. _d 0
295 Av_SI(i,j,bi,bj) = 0. _d 0
296 ENDDO
297 ENDDO
298 ENDDO
299 ENDDO
300
301 ! IF (STREAMICE_construct_matrix) THEN
302
303 ! #ifdef STREAMICE_CONSTRUCT_MATRIX
304
305 DO bj = myByLo(myThid), myByHi(myThid)
306 DO bi = myBxLo(myThid), myBxHi(myThid)
307 DO j=js,je
308 DO i=is,ie
309 DO colx=-1,1
310 DO coly=-1,1
311 Au_SI(i,j,bi,bj) = Au_SI(i,j,bi,bj) +
312 & A_uu(i,j,bi,bj,colx,coly)*
313 & Du_SI(i+colx,j+coly,bi,bj)+
314 & A_uv(i,j,bi,bj,colx,coly)*
315 & Dv_SI(i+colx,j+coly,bi,bj)
316 Av_SI(i,j,bi,bj) = Av_SI(i,j,bi,bj) +
317 & A_vu(i,j,bi,bj,colx,coly)*
318 & Du_SI(i+colx,j+coly,bi,bj)+
319 & A_vv(i,j,bi,bj,colx,coly)*
320 & Dv_SI(i+colx,j+coly,bi,bj)
321 ENDDO
322 ENDDO
323 ENDDO
324 ENDDO
325 ENDDO
326 ENDDO
327
328 ! else
329 ! #else
330 !
331 ! CALL STREAMICE_CG_ACTION( myThid,
332 ! O Au_SI,
333 ! O Av_SI,
334 ! I Du_SI,
335 ! I Dv_SI,
336 ! I is,ie,js,je)
337 !
338 ! ! ENDIF
339 !
340 ! #endif
341
342
343 DO bj = myByLo(myThid), myByHi(myThid)
344 DO bi = myBxLo(myThid), myBxHi(myThid)
345 dot_p1_tile(bi,bj) = 0. _d 0
346 dot_p2_tile(bi,bj) = 0. _d 0
347 ENDDO
348 ENDDO
349
350 DO bj = myByLo(myThid), myByHi(myThid)
351 DO bi = myBxLo(myThid), myBxHi(myThid)
352 DO j=1,sNy
353 DO i=1,sNx
354 IF (STREAMICE_umask(i,j,bi,bj).eq.1.0) THEN
355 dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Zu_SI(i,j,bi,bj)*
356 & Ru_SI(i,j,bi,bj)
357 dot_p2_tile(bi,bj)=dot_p2_tile(bi,bj)+Du_SI(i,j,bi,bj)*
358 & Au_SI(i,j,bi,bj)
359 ENDIF
360 IF (STREAMICE_vmask(i,j,bi,bj).eq.1.0) THEN
361 dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Zv_SI(i,j,bi,bj)*
362 & Rv_SI(i,j,bi,bj)
363 dot_p2_tile(bi,bj)=dot_p2_tile(bi,bj)+Dv_SI(i,j,bi,bj)*
364 & Av_SI(i,j,bi,bj)
365 ENDIF
366 ENDDO
367 ENDDO
368 ENDDO
369 ENDDO
370
371 CALL GLOBAL_SUM_TILE_RL( dot_p1_tile, dot_p1, myThid )
372 CALL GLOBAL_SUM_TILE_RL( dot_p2_tile, dot_p2, myThid )
373 alpha_k = dot_p1/dot_p2
374
375 DO bj = myByLo(myThid), myByHi(myThid)
376 DO bi = myBxLo(myThid), myBxHi(myThid)
377 DO j=1-OLy,sNy+OLy
378 DO i=1-OLx,sNx+OLx
379
380 IF (STREAMICE_umask(i,j,bi,bj).eq.1.0) THEN
381 cg_Uin(i,j,bi,bj)=cg_Uin(i,j,bi,bj)+
382 & alpha_k*Du_SI(i,j,bi,bj)
383 Ru_old_SI(i,j,bi,bj) = Ru_SI(i,j,bi,bj)
384 Zu_old_SI(i,j,bi,bj) = Zu_SI(i,j,bi,bj)
385 Ru_SI(i,j,bi,bj) = Ru_SI(i,j,bi,bj)-
386 & alpha_k*Au_SI(i,j,bi,bj)
387 Zu_SI(i,j,bi,bj) = Ru_SI(i,j,bi,bj) /
388 & DIAGu_SI(i,j,bi,bj)
389 ENDIF
390
391 IF (STREAMICE_vmask(i,j,bi,bj).eq.1.0) THEN
392 cg_Vin(i,j,bi,bj)=cg_Vin(i,j,bi,bj)+
393 & alpha_k*Dv_SI(i,j,bi,bj)
394 Rv_old_SI(i,j,bi,bj) = Rv_SI(i,j,bi,bj)
395 Zv_old_SI(i,j,bi,bj) = Zv_SI(i,j,bi,bj)
396 Rv_SI(i,j,bi,bj) = Rv_SI(i,j,bi,bj)-
397 & alpha_k*Av_SI(i,j,bi,bj)
398 Zv_SI(i,j,bi,bj) = Rv_SI(i,j,bi,bj) /
399 & DIAGv_SI(i,j,bi,bj)
400
401 ENDIF
402 ENDDO
403 ENDDO
404 ENDDO
405 ENDDO
406
407 DO bj = myByLo(myThid), myByHi(myThid)
408 DO bi = myBxLo(myThid), myBxHi(myThid)
409 dot_p1_tile(bi,bj) = 0. _d 0
410 dot_p2_tile(bi,bj) = 0. _d 0
411 ENDDO
412 ENDDO
413
414 DO bj = myByLo(myThid), myByHi(myThid)
415 DO bi = myBxLo(myThid), myBxHi(myThid)
416 DO j=1,sNy
417 DO i=1,sNx
418
419 IF (STREAMICE_umask(i,j,bi,bj).eq.1.0) THEN
420 dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Zu_SI(i,j,bi,bj)*
421 & Ru_SI(i,j,bi,bj)
422 dot_p2_tile(bi,bj)=dot_p2_tile(bi,bj)+Zu_old_SI(i,j,bi,bj)*
423 & Ru_old_SI(i,j,bi,bj)
424 ENDIF
425
426 IF (STREAMICE_vmask(i,j,bi,bj).eq.1.0) THEN
427 dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Zv_SI(i,j,bi,bj)*
428 & Rv_SI(i,j,bi,bj)
429 dot_p2_tile(bi,bj)=dot_p2_tile(bi,bj)+Zv_old_SI(i,j,bi,bj)*
430 & Rv_old_SI(i,j,bi,bj)
431 ENDIF
432
433 ENDDO
434 ENDDO
435 ENDDO
436 ENDDO
437
438 CALL GLOBAL_SUM_TILE_RL( dot_p1_tile, dot_p1, myThid )
439 CALL GLOBAL_SUM_TILE_RL( dot_p2_tile, dot_p2, myThid )
440
441 beta_k = dot_p1/dot_p2
442
443 DO bj = myByLo(myThid), myByHi(myThid)
444 DO bi = myBxLo(myThid), myBxHi(myThid)
445 DO j=1-OLy,sNy+OLy
446 DO i=1-OLx,sNx+OLx
447 IF (STREAMICE_umask(i,j,bi,bj).eq.1.0)
448 & Du_SI(i,j,bi,bj)=beta_k*Du_SI(i,j,bi,bj)+
449 & Zu_SI(i,j,bi,bj)
450 IF (STREAMICE_vmask(i,j,bi,bj).eq.1.0)
451 & Dv_SI(i,j,bi,bj)=beta_k*Dv_SI(i,j,bi,bj)+
452 & Zv_SI(i,j,bi,bj)
453 ENDDO
454 ENDDO
455 ENDDO
456 ENDDO
457
458 DO bj = myByLo(myThid), myByHi(myThid)
459 DO bi = myBxLo(myThid), myBxHi(myThid)
460 dot_p1_tile(bi,bj) = 0. _d 0
461 ENDDO
462 ENDDO
463
464 DO bj = myByLo(myThid), myByHi(myThid)
465 DO bi = myBxLo(myThid), myBxHi(myThid)
466 DO j=1,sNy
467 DO i=1,sNx
468 IF (STREAMICE_umask(i,j,bi,bj).eq.1.0)
469 & dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Ru_SI(i,j,bi,bj)**2
470 IF (STREAMICE_vmask(i,j,bi,bj).eq.1.0)
471 & dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Rv_SI(i,j,bi,bj)**2
472 ENDDO
473 ENDDO
474 ENDDO
475 ENDDO
476
477 CALL GLOBAL_SUM_TILE_RL( dot_p1_tile, dot_p1, myThid )
478 resid = sqrt(dot_p1)
479
480 ! IF (iter .eq. 1) then
481 ! print *, alpha_k, beta_k, resid
482 ! ENDIF
483
484 cg_halo = cg_halo - 1
485
486 if (cg_halo .eq. 0) then
487 cg_halo = min(OLx-1,OLy-1)
488 _EXCH_XY_RL( Du_SI, myThid )
489 _EXCH_XY_RL( Dv_SI, myThid )
490 _EXCH_XY_RL( Ru_SI, myThid )
491 _EXCH_XY_RL( Rv_SI, myThid )
492 _EXCH_XY_RL( cg_Uin, myThid )
493 _EXCH_XY_RL( cg_Vin, myThid )
494 endif
495
496
497 endif
498 enddo ! end of CG loop
499
500 c to avoid using "exit"
501 c if iters has reached max_iters there is no convergence
502
503 IF (iters .lt. maxIter) THEN
504 conv_flag = 1
505 ENDIF
506
507 ! DO bj = myByLo(myThid), myByHi(myThid)
508 ! DO bi = myBxLo(myThid), myBxHi(myThid)
509 ! DO j=1-OLy,sNy+OLy
510 ! DO i=1-OLy,sNx+OLy
511 ! IF (STREAMICE_umask(i,j,bi,bj).eq.3.0)
512 ! & cg_Uin(i,j,bi,bj)=u_bdry_values_SI(i,j,bi,bj)
513 ! IF (STREAMICE_vmask(i,j,bi,bj).eq.3.0)
514 ! & cg_Vin(i,j,bi,bj)=v_bdry_values_SI(i,j,bi,bj)
515 ! ENDDO
516 ! ENDDO
517 ! ENDDO
518 ! ENDDO
519 !
520 ! _EXCH_XY_RL( cg_Uin, myThid )
521 ! _EXCH_XY_RL( cg_Vin, myThid )
522
523 #endif /* ifndef ALLOW_PETSC */
524
525 #else /* STREAMICE_SERIAL_TRISOLVE */
526
527 iters = 0
528
529 CALL STREAMICE_TRIDIAG_SOLVE(
530 U cg_Uin, ! x-velocities
531 U cg_Vin,
532 U cg_Bu, ! force in x dir
533 I A_uu, ! section of matrix that multiplies u and projects on u
534 I STREAMICE_umask,
535 I myThid )
536
537 #endif
538
539 CALL TIMER_STOP ('STREAMICE_CG_SOLVE',myThid)
540
541
542 #endif
543 RETURN
544 END

  ViewVC Help
Powered by ViewVC 1.1.22