| 1 | dgoldberg | 1.9 | C $Header: /u/gcmpack/MITgcm_contrib/dgoldberg/streamice/streamice_cg_solve.F,v 1.8 2013/05/28 22:32:39 dgoldberg Exp $ | 
| 2 | heimbach | 1.1 | 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 | dgoldberg | 1.4 | 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 | heimbach | 1.1 | I                               tolerance, | 
| 19 | heimbach | 1.2 | O                               iters, | 
| 20 |  |  | I                               myThid ) | 
| 21 | heimbach | 1.1 | C     /============================================================\ | 
| 22 |  |  | C     | SUBROUTINE                                                 | | 
| 23 |  |  | C     | o                                                          | | 
| 24 |  |  | C     |============================================================| | 
| 25 |  |  | C     |                                                            | | 
| 26 |  |  | C     \============================================================/ | 
| 27 |  |  | IMPLICIT NONE | 
| 28 |  |  |  | 
| 29 |  |  | #include "SIZE.h" | 
| 30 |  |  | #include "EEPARAMS.h" | 
| 31 |  |  | #include "PARAMS.h" | 
| 32 |  |  | #include "STREAMICE.h" | 
| 33 |  |  | #include "STREAMICE_CG.h" | 
| 34 | dgoldberg | 1.8 |  | 
| 35 |  |  |  | 
| 36 |  |  |  | 
| 37 | dgoldberg | 1.7 | #ifdef ALLOW_PETSC | 
| 38 |  |  | #include "finclude/petsc.h" | 
| 39 | dgoldberg | 1.9 | ! UNCOMMENT IF V3.0 | 
| 40 |  |  | !#include "finclude/petscvec.h" | 
| 41 |  |  | !#include "finclude/petscmat.h" | 
| 42 |  |  | !#include "finclude/petscksp.h" | 
| 43 |  |  | !#include "finclude/petscpc.h" | 
| 44 | dgoldberg | 1.7 | #endif | 
| 45 |  |  | C     === Global variables === | 
| 46 | heimbach | 1.1 |  | 
| 47 |  |  |  | 
| 48 |  |  | C     !INPUT/OUTPUT ARGUMENTS | 
| 49 |  |  | C     cg_Uin, cg_Vin - input and output velocities | 
| 50 |  |  | C     cg_Bu, cg_Bv - driving stress | 
| 51 |  |  | INTEGER myThid | 
| 52 |  |  | INTEGER iters | 
| 53 |  |  | _RL tolerance | 
| 54 |  |  | _RL cg_Uin (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) | 
| 55 |  |  | _RL cg_Vin (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) | 
| 56 |  |  | _RL cg_Bu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) | 
| 57 |  |  | _RL cg_Bv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) | 
| 58 | dgoldberg | 1.4 | _RL | 
| 59 |  |  | & A_uu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1), | 
| 60 |  |  | & A_vu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1), | 
| 61 |  |  | & A_uv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1), | 
| 62 |  |  | & A_vv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1) | 
| 63 | heimbach | 1.1 |  | 
| 64 |  |  | C     LOCAL VARIABLES | 
| 65 |  |  | INTEGER i, j, bi, bj, cg_halo, conv_flag | 
| 66 |  |  | INTEGER iter, is, js, ie, je, colx, coly, k | 
| 67 | dgoldberg | 1.4 | _RL dot_p1, dot_p2, alpha_k, beta_k, resid, resid_0 | 
| 68 | heimbach | 1.1 | _RL dot_p1_tile (nSx,nSy) | 
| 69 |  |  | _RL dot_p2_tile (nSx,nSy) | 
| 70 | dgoldberg | 1.5 | CHARACTER*(MAX_LEN_MBUF) msgBuf | 
| 71 | heimbach | 1.1 |  | 
| 72 | dgoldberg | 1.7 |  | 
| 73 |  |  | #ifdef ALLOW_PETSC | 
| 74 |  |  | INTEGER indices(2*(snx*nsx*sny*nsy)) | 
| 75 |  |  | INTEGER n_dofs_cum_sum (0:nPx*nPy-1), idx(1) | 
| 76 |  |  | _RL rhs_values(2*(snx*nsx*sny*nsy)) | 
| 77 |  |  | _RL solution_values(2*(snx*nsx*sny*nsy)) | 
| 78 |  |  | !      _RL mat_values (2*Nx*Ny,2*(snx*nsx*sny*nsy)) | 
| 79 |  |  | _RL mat_values (18,1), mat_val_return(1) | 
| 80 |  |  | INTEGER indices_col(18) | 
| 81 |  |  | INTEGER local_dofs, global_dofs, dof_index, dof_index_col | 
| 82 |  |  | INTEGER local_offset | 
| 83 |  |  | Mat matrix | 
| 84 |  |  | KSP ksp | 
| 85 |  |  | PC  pc | 
| 86 |  |  | Vec rhs | 
| 87 |  |  | Vec solution | 
| 88 |  |  | PetscErrorCode ierr | 
| 89 |  |  | #ifdef ALLOW_USE_MPI | 
| 90 |  |  | integer mpiRC, mpiMyWid | 
| 91 |  |  | #endif | 
| 92 |  |  | #endif | 
| 93 |  |  |  | 
| 94 | heimbach | 1.1 |  | 
| 95 |  |  | #ifdef ALLOW_STREAMICE | 
| 96 |  |  |  | 
| 97 | dgoldberg | 1.8 |  | 
| 98 |  |  |  | 
| 99 | dgoldberg | 1.7 | CALL TIMER_START ('STREAMICE_CG_SOLVE',myThid) | 
| 100 | dgoldberg | 1.8 | #ifndef STREAMICE_SERIAL_TRISOLVE | 
| 101 | dgoldberg | 1.7 |  | 
| 102 |  |  | #ifdef ALLOW_PETSC | 
| 103 |  |  |  | 
| 104 |  |  | #ifdef ALLOW_USE_MPI | 
| 105 |  |  |  | 
| 106 |  |  |  | 
| 107 |  |  | CALL MPI_COMM_RANK( MPI_COMM_WORLD, mpiMyWId, mpiRC ) | 
| 108 |  |  | local_dofs = n_dofs_process (mpiMyWid) | 
| 109 |  |  | global_dofs = 0 | 
| 110 |  |  |  | 
| 111 |  |  | n_dofs_cum_sum(0) = 0 | 
| 112 |  |  | DO i=0,nPx*nPy-1 | 
| 113 |  |  | global_dofs = global_dofs + n_dofs_process (i) | 
| 114 |  |  | if (i.ge.1) THEN | 
| 115 |  |  | n_dofs_cum_sum(i) = n_dofs_cum_sum(i-1)+ | 
| 116 |  |  | &                     n_dofs_process(i-1) | 
| 117 |  |  | endif | 
| 118 |  |  | ENDDO | 
| 119 |  |  | local_offset = n_dofs_cum_sum(mpimywid) | 
| 120 |  |  |  | 
| 121 |  |  | #else | 
| 122 |  |  |  | 
| 123 |  |  | local_dofs = n_dofs_process (0) | 
| 124 |  |  | global_dofs = local_dofs | 
| 125 |  |  | local_offset = 0 | 
| 126 |  |  |  | 
| 127 |  |  | #endif | 
| 128 |  |  |  | 
| 129 |  |  | !      call petscInitialize(PETSC_NULL_CHARACTER,ierr) | 
| 130 |  |  |  | 
| 131 |  |  | !---------------------- | 
| 132 |  |  |  | 
| 133 |  |  | call VecCreate(PETSC_COMM_WORLD, rhs, ierr) | 
| 134 |  |  | call VecSetSizes(rhs, local_dofs, global_dofs, ierr) | 
| 135 |  |  | call VecSetType(rhs, VECMPI, ierr) | 
| 136 |  |  |  | 
| 137 |  |  | call VecCreate(PETSC_COMM_WORLD, solution, ierr) | 
| 138 |  |  | call VecSetSizes(solution, local_dofs, global_dofs, ierr) | 
| 139 |  |  | call VecSetType(solution, VECMPI, ierr) | 
| 140 |  |  |  | 
| 141 |  |  | do i=1,local_dofs | 
| 142 |  |  | indices(i) = i-1 + local_offset | 
| 143 |  |  | end do | 
| 144 |  |  | do i=1,2*nSx*nSy*sNx*sNy | 
| 145 |  |  | rhs_values (i) = 0. _d 0 | 
| 146 |  |  | solution_values (i) = 0. _d 0 | 
| 147 |  |  | enddo | 
| 148 |  |  |  | 
| 149 |  |  | ! gather rhs and initial guess values to populate petsc vectors | 
| 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 |  |  |  | 
| 156 |  |  | dof_index = INT(streamice_petsc_dofs_u(i,j,bi,bj)) | 
| 157 |  |  | &                - local_offset | 
| 158 |  |  |  | 
| 159 |  |  | if (dof_index.ge.0) THEN | 
| 160 |  |  |  | 
| 161 |  |  | rhs_values(dof_index+1) = cg_Bu(i,j,bi,bj) | 
| 162 |  |  | solution_values(dof_index+1) = cg_Uin(i,j,bi,bj) | 
| 163 |  |  |  | 
| 164 |  |  | endif | 
| 165 |  |  |  | 
| 166 |  |  | !--------------- | 
| 167 |  |  |  | 
| 168 |  |  | dof_index = INT(streamice_petsc_dofs_v(i,j,bi,bj)) | 
| 169 |  |  | &                - local_offset | 
| 170 |  |  |  | 
| 171 |  |  | if (dof_index.ge.0) THEN | 
| 172 |  |  |  | 
| 173 |  |  | rhs_values(dof_index+1) = cg_Bv(i,j,bi,bj) | 
| 174 |  |  | solution_values(dof_index+1) = cg_Vin(i,j,bi,bj) | 
| 175 |  |  |  | 
| 176 |  |  | endif | 
| 177 |  |  |  | 
| 178 |  |  | ENDDO | 
| 179 |  |  | ENDDO | 
| 180 |  |  | ENDDO | 
| 181 |  |  | ENDDO | 
| 182 |  |  |  | 
| 183 |  |  |  | 
| 184 |  |  | call VecSetValues(rhs, local_dofs, indices, rhs_values, | 
| 185 |  |  | &                  INSERT_VALUES, ierr) | 
| 186 |  |  | call VecAssemblyBegin(rhs, ierr) | 
| 187 |  |  | call VecAssemblyEnd(rhs, ierr) | 
| 188 |  |  |  | 
| 189 |  |  |  | 
| 190 |  |  | call VecSetValues(solution, local_dofs, indices, | 
| 191 |  |  | &                  solution_values, INSERT_VALUES, ierr) | 
| 192 |  |  | call VecAssemblyBegin(solution, ierr) | 
| 193 |  |  | call VecAssemblyEnd(solution, ierr) | 
| 194 |  |  |  | 
| 195 | dgoldberg | 1.9 | !     IF USING v3.0 THEN | 
| 196 |  |  | !     call MatCreateMPIAIJ (PETSC_COMM_WORLD, | 
| 197 |  |  | call MatCreateAIJ (PETSC_COMM_WORLD, | 
| 198 | dgoldberg | 1.7 | &                      local_dofs, local_dofs, | 
| 199 |  |  | &                      global_dofs, global_dofs, | 
| 200 |  |  | &                      18, PETSC_NULL_INTEGER, | 
| 201 |  |  | &                      18, PETSC_NULL_INTEGER, | 
| 202 |  |  | &                      matrix, ierr) | 
| 203 |  |  |  | 
| 204 | dgoldberg | 1.8 |  | 
| 205 | dgoldberg | 1.7 | ! populate petsc matrix | 
| 206 |  |  |  | 
| 207 |  |  | DO bj = myByLo(myThid), myByHi(myThid) | 
| 208 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 209 |  |  | DO j=1,sNy | 
| 210 |  |  | DO i=1,sNx | 
| 211 |  |  |  | 
| 212 |  |  | dof_index = INT(streamice_petsc_dofs_u(i,j,bi,bj)) | 
| 213 |  |  | !     &                - local_offset | 
| 214 |  |  |  | 
| 215 |  |  | IF (dof_index .ge. 0) THEN | 
| 216 |  |  |  | 
| 217 |  |  | DO k=1,18 | 
| 218 |  |  | indices_col(k) = 0 | 
| 219 |  |  | mat_values(k,1) = 0. _d 0 | 
| 220 |  |  | ENDDO | 
| 221 |  |  | k=0 | 
| 222 |  |  |  | 
| 223 |  |  | DO coly=-1,1 | 
| 224 |  |  | DO colx=-1,1 | 
| 225 |  |  |  | 
| 226 |  |  | dof_index_col = streamice_petsc_dofs_u(i+colx,j+coly,bi,bj) | 
| 227 | dgoldberg | 1.4 |  | 
| 228 | dgoldberg | 1.7 | if (dof_index_col.ge.0) THEN | 
| 229 |  |  | !               pscal = A_uu(i,j,bi,bj,colx,coly) | 
| 230 |  |  | !               CALL MatSetValue (matrix,dof_index, dof_index_col, | 
| 231 |  |  | !     &              pscal,INSERT_VALUES,ierr) | 
| 232 |  |  | k=k+1 | 
| 233 |  |  | mat_values (k,1) = A_uu(i,j,bi,bj,colx,coly) | 
| 234 |  |  | indices_col (k) = dof_index_col | 
| 235 |  |  | endif | 
| 236 |  |  |  | 
| 237 |  |  | dof_index_col = streamice_petsc_dofs_v(i+colx,j+coly,bi,bj) | 
| 238 |  |  |  | 
| 239 |  |  | if (dof_index_col.ge.0) THEN | 
| 240 |  |  | !               CALL MatSetValue (matrix,dof_index, dof_index_col, | 
| 241 |  |  | !     &              A_uv(i,j,bi,bj,colx,coly),INSERT_VALUES,ierr) | 
| 242 |  |  | k=k+1 | 
| 243 |  |  | mat_values (k,1) = A_uv(i,j,bi,bj,colx,coly) | 
| 244 |  |  | indices_col (k) = dof_index_col | 
| 245 |  |  | endif | 
| 246 |  |  |  | 
| 247 |  |  | ENDDO | 
| 248 |  |  | ENDDO | 
| 249 |  |  |  | 
| 250 |  |  | call matSetValues (matrix, 1, dof_index, k, indices_col, | 
| 251 |  |  | &                        mat_values,INSERT_VALUES,ierr) | 
| 252 |  |  |  | 
| 253 |  |  |  | 
| 254 |  |  | ENDIF | 
| 255 |  |  |  | 
| 256 |  |  | ! ---------------------------------------------- | 
| 257 |  |  |  | 
| 258 |  |  | dof_index = INT(streamice_petsc_dofs_v(i,j,bi,bj)) | 
| 259 |  |  | !     &                - local_offset | 
| 260 |  |  |  | 
| 261 |  |  | IF (dof_index .ge. 0) THEN | 
| 262 |  |  |  | 
| 263 |  |  | DO k=1,18 | 
| 264 |  |  | indices_col(k) = 0 | 
| 265 |  |  | mat_values(k,1) = 0. _d 0 | 
| 266 |  |  | ENDDO | 
| 267 |  |  | k=0 | 
| 268 |  |  |  | 
| 269 |  |  | DO coly=-1,1 | 
| 270 |  |  | DO colx=-1,1 | 
| 271 |  |  |  | 
| 272 |  |  | dof_index_col = streamice_petsc_dofs_u(i+colx,j+coly,bi,bj) | 
| 273 |  |  |  | 
| 274 |  |  | if (dof_index_col.ge.0) THEN | 
| 275 |  |  | !               CALL MatSetValue (matrix,dof_index, dof_index_col, | 
| 276 |  |  | !     &              A_vu(i,j,bi,bj,colx,coly),INSERT_VALUES,ierr) | 
| 277 |  |  | k=k+1 | 
| 278 |  |  | mat_values (k,1) = A_vu(i,j,bi,bj,colx,coly) | 
| 279 |  |  | indices_col (k) = dof_index_col | 
| 280 |  |  | endif | 
| 281 |  |  |  | 
| 282 |  |  | dof_index_col = streamice_petsc_dofs_v(i+colx,j+coly,bi,bj) | 
| 283 |  |  |  | 
| 284 |  |  | if (dof_index_col.ge.0) THEN | 
| 285 |  |  | !               CALL MatSetValue (matrix,dof_index, dof_index_col, | 
| 286 |  |  | !     &              A_vv(i,j,bi,bj,colx,coly),INSERT_VALUES,ierr) | 
| 287 |  |  | k=k+1 | 
| 288 |  |  | mat_values (k,1) = A_vv(i,j,bi,bj,colx,coly) | 
| 289 |  |  | indices_col (k) = dof_index_col | 
| 290 |  |  | endif | 
| 291 |  |  |  | 
| 292 |  |  | ENDDO | 
| 293 |  |  | ENDDO | 
| 294 |  |  |  | 
| 295 |  |  | call matSetValues (matrix, 1, dof_index, k, indices_col, | 
| 296 |  |  | &                        mat_values,INSERT_VALUES,ierr) | 
| 297 |  |  | ENDIF | 
| 298 |  |  |  | 
| 299 |  |  | ENDDO | 
| 300 |  |  | ENDDO | 
| 301 |  |  | ENDDO | 
| 302 |  |  | ENDDO | 
| 303 |  |  |  | 
| 304 |  |  | call MatAssemblyBegin(matrix,MAT_FINAL_ASSEMBLY,ierr) | 
| 305 |  |  | call MatAssemblyEnd(matrix,MAT_FINAL_ASSEMBLY,ierr) | 
| 306 |  |  |  | 
| 307 |  |  |  | 
| 308 |  |  | call KSPCreate(PETSC_COMM_WORLD, ksp, ierr) | 
| 309 |  |  | call KSPSetOperators(ksp, matrix, matrix, | 
| 310 |  |  | &                     DIFFERENT_NONZERO_PATTERN, ierr) | 
| 311 |  |  |  | 
| 312 |  |  | SELECT CASE (PETSC_SOLVER_TYPE) | 
| 313 |  |  | CASE ('CG') | 
| 314 |  |  | PRINT *, "PETSC SOLVER: SELECTED CG" | 
| 315 |  |  | call KSPSetType(ksp, KSPCG, ierr) | 
| 316 |  |  | CASE ('GMRES') | 
| 317 |  |  | PRINT *, "PETSC SOLVER: SELECTED GMRES" | 
| 318 |  |  | call KSPSetType(ksp, KSPGMRES, ierr) | 
| 319 |  |  | CASE ('BICG') | 
| 320 |  |  | PRINT *, "PETSC SOLVER: SELECTED BICG" | 
| 321 |  |  | call KSPSetType(ksp, KSPBICG, ierr) | 
| 322 |  |  | CASE DEFAULT | 
| 323 |  |  | PRINT *, "PETSC SOLVER: SELECTED DEFAULT" | 
| 324 |  |  | call KSPSetType(ksp, KSPCG, ierr) | 
| 325 |  |  | END SELECT | 
| 326 |  |  |  | 
| 327 |  |  | call KSPGetPC(ksp, pc, ierr) | 
| 328 |  |  | call KSPSetTolerances(ksp,tolerance, | 
| 329 |  |  | &     PETSC_DEFAULT_DOUBLE_PRECISION, | 
| 330 |  |  | &     PETSC_DEFAULT_DOUBLE_PRECISION, | 
| 331 |  |  | &     streamice_max_cg_iter,ierr) | 
| 332 |  |  |  | 
| 333 |  |  | SELECT CASE (PETSC_PRECOND_TYPE) | 
| 334 |  |  | CASE ('BLOCKJACOBI') | 
| 335 |  |  | PRINT *, "PETSC PRECOND: SELECTED BJACOBI" | 
| 336 |  |  | call PCSetType(pc, PCBJACOBI, ierr) | 
| 337 |  |  | CASE ('JACOBI') | 
| 338 |  |  | PRINT *, "PETSC PRECOND: SELECTED JACOBI" | 
| 339 |  |  | call PCSetType(pc, PCJACOBI, ierr) | 
| 340 |  |  | CASE ('ILU') | 
| 341 |  |  | PRINT *, "PETSC PRECOND: SELECTED ILU" | 
| 342 |  |  | call PCSetType(pc, PCILU, ierr) | 
| 343 |  |  | CASE DEFAULT | 
| 344 |  |  | PRINT *, "PETSC PRECOND: SELECTED DEFAULT" | 
| 345 |  |  | call PCSetType(pc, PCBJACOBI, ierr) | 
| 346 |  |  | END SELECT | 
| 347 |  |  |  | 
| 348 |  |  |  | 
| 349 |  |  | call KSPSolve(ksp, rhs, solution, ierr) | 
| 350 |  |  | call KSPGetIterationNumber(ksp,iters,ierr) | 
| 351 |  |  |  | 
| 352 |  |  | call VecGetValues(solution,local_dofs,indices, | 
| 353 |  |  | &      solution_values,ierr) | 
| 354 |  |  |  | 
| 355 |  |  | DO bj = myByLo(myThid), myByHi(myThid) | 
| 356 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 357 |  |  | DO j=1,sNy | 
| 358 |  |  | DO i=1,sNx | 
| 359 |  |  |  | 
| 360 |  |  | dof_index = INT(streamice_petsc_dofs_u(i,j,bi,bj)) | 
| 361 |  |  | &                - local_offset | 
| 362 |  |  | if (dof_index.ge.0) THEN | 
| 363 |  |  | cg_Uin(i,j,bi,bj) = solution_values(dof_index+1) | 
| 364 |  |  | endif | 
| 365 |  |  |  | 
| 366 |  |  | dof_index = INT(streamice_petsc_dofs_v(i,j,bi,bj)) | 
| 367 |  |  | &                - local_offset | 
| 368 |  |  | if (dof_index.ge.0) THEN | 
| 369 |  |  | cg_Vin(i,j,bi,bj) = solution_values(dof_index+1) | 
| 370 |  |  | endif | 
| 371 |  |  |  | 
| 372 |  |  | ENDDO | 
| 373 |  |  | ENDDO | 
| 374 |  |  | ENDDO | 
| 375 |  |  | ENDDO | 
| 376 |  |  |  | 
| 377 |  |  |  | 
| 378 |  |  |  | 
| 379 | dgoldberg | 1.8 | #else  /* ALLOW_PETSC */ | 
| 380 | dgoldberg | 1.7 |  | 
| 381 |  |  |  | 
| 382 |  |  | iters = streamice_max_cg_iter | 
| 383 | heimbach | 1.1 | conv_flag = 0 | 
| 384 |  |  |  | 
| 385 | dgoldberg | 1.4 |  | 
| 386 | heimbach | 1.1 | DO bj = myByLo(myThid), myByHi(myThid) | 
| 387 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 388 |  |  | DO j=1,sNy | 
| 389 |  |  | DO i=1,sNx | 
| 390 |  |  | Zu_SI (i,j,bi,bj) = 0. _d 0 | 
| 391 |  |  | Zv_SI (i,j,bi,bj) = 0. _d 0 | 
| 392 |  |  | Ru_SI (i,j,bi,bj) = 0. _d 0 | 
| 393 |  |  | Rv_SI (i,j,bi,bj) = 0. _d 0 | 
| 394 |  |  | Au_SI (i,j,bi,bj) = 0. _d 0 | 
| 395 |  |  | Av_SI (i,j,bi,bj) = 0. _d 0 | 
| 396 |  |  | Du_SI (i,j,bi,bj) = 0. _d 0 | 
| 397 |  |  | Dv_SI (i,j,bi,bj) = 0. _d 0 | 
| 398 |  |  | ENDDO | 
| 399 |  |  | ENDDO | 
| 400 |  |  | ENDDO | 
| 401 |  |  | ENDDO | 
| 402 | dgoldberg | 1.4 |  | 
| 403 |  |  | C     FIND INITIAL RESIDUAL, and initialize r | 
| 404 | heimbach | 1.1 |  | 
| 405 | dgoldberg | 1.4 | ! #ifdef STREAMICE_CONSTRUCT_MATRIX | 
| 406 | heimbach | 1.1 |  | 
| 407 | dgoldberg | 1.4 |  | 
| 408 |  |  |  | 
| 409 |  |  | DO bj = myByLo(myThid), myByHi(myThid) | 
| 410 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 411 |  |  | DO j=1,sNy | 
| 412 |  |  | DO i=1,sNx | 
| 413 |  |  | DO colx=-1,1 | 
| 414 |  |  | DO coly=-1,1 | 
| 415 |  |  | Au_SI(i,j,bi,bj) = Au_SI(i,j,bi,bj) + | 
| 416 |  |  | &         A_uu(i,j,bi,bj,colx,coly)* | 
| 417 |  |  | &         cg_Uin(i+colx,j+coly,bi,bj)+ | 
| 418 |  |  | &         A_uv(i,j,bi,bj,colx,coly)* | 
| 419 |  |  | &         cg_Vin(i+colx,j+coly,bi,bj) | 
| 420 | heimbach | 1.1 |  | 
| 421 | dgoldberg | 1.4 |  | 
| 422 |  |  | Av_SI(i,j,bi,bj) = Av_SI(i,j,bi,bj) + | 
| 423 |  |  | &         A_vu(i,j,bi,bj,colx,coly)* | 
| 424 |  |  | &         cg_Uin(i+colx,j+coly,bi,bj)+ | 
| 425 |  |  | &         A_vv(i,j,bi,bj,colx,coly)* | 
| 426 |  |  | &         cg_Vin(i+colx,j+coly,bi,bj) | 
| 427 |  |  | ENDDO | 
| 428 |  |  | ENDDO | 
| 429 |  |  | ENDDO | 
| 430 |  |  | ENDDO | 
| 431 | heimbach | 1.1 | ENDDO | 
| 432 |  |  | ENDDO | 
| 433 |  |  |  | 
| 434 | dgoldberg | 1.4 |  | 
| 435 | heimbach | 1.1 | _EXCH_XY_RL( Au_SI, myThid ) | 
| 436 |  |  | _EXCH_XY_RL( Av_SI, myThid ) | 
| 437 |  |  |  | 
| 438 | dgoldberg | 1.4 |  | 
| 439 | heimbach | 1.1 | DO bj = myByLo(myThid), myByHi(myThid) | 
| 440 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 441 |  |  | DO j=1-OLy,sNy+OLy | 
| 442 |  |  | DO i=1-OLx,sNx+OLx | 
| 443 | dgoldberg | 1.4 | Ru_SI(i,j,bi,bj)=cg_Bu(i,j,bi,bj)- | 
| 444 | heimbach | 1.1 | &     Au_SI(i,j,bi,bj) | 
| 445 | dgoldberg | 1.4 | Rv_SI(i,j,bi,bj)=cg_Bv(i,j,bi,bj)- | 
| 446 | heimbach | 1.1 | &     Av_SI(i,j,bi,bj) | 
| 447 |  |  | ENDDO | 
| 448 |  |  | ENDDO | 
| 449 |  |  | dot_p1_tile(bi,bj) = 0. _d 0 | 
| 450 |  |  | dot_p2_tile(bi,bj) = 0. _d 0 | 
| 451 |  |  | ENDDO | 
| 452 |  |  | ENDDO | 
| 453 |  |  |  | 
| 454 | dgoldberg | 1.4 |  | 
| 455 |  |  |  | 
| 456 | heimbach | 1.1 | DO bj = myByLo(myThid), myByHi(myThid) | 
| 457 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 458 |  |  | DO j=1,sNy | 
| 459 |  |  | DO i=1,sNx | 
| 460 |  |  | IF (STREAMICE_umask(i,j,bi,bj).eq.1.0) | 
| 461 |  |  | &      dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Ru_SI(i,j,bi,bj)**2 | 
| 462 |  |  | IF (STREAMICE_vmask(i,j,bi,bj).eq.1.0) | 
| 463 |  |  | &      dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Rv_SI(i,j,bi,bj)**2 | 
| 464 |  |  | ENDDO | 
| 465 |  |  | ENDDO | 
| 466 |  |  | ENDDO | 
| 467 |  |  | ENDDO | 
| 468 |  |  |  | 
| 469 |  |  | CALL GLOBAL_SUM_TILE_RL( dot_p1_tile, dot_p1, myThid ) | 
| 470 |  |  | resid_0 = sqrt(dot_p1) | 
| 471 |  |  |  | 
| 472 | dgoldberg | 1.5 | WRITE(msgBuf,'(A,E14.7)') 'CONJ GRAD INIT RESID, ', | 
| 473 |  |  | &         resid_0 | 
| 474 |  |  | CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, | 
| 475 |  |  | &                     SQUEEZE_RIGHT , 1) | 
| 476 |  |  |  | 
| 477 | dgoldberg | 1.4 | C    CCCCCCCCCCCCCCCCCCCC | 
| 478 |  |  |  | 
| 479 | heimbach | 1.1 | DO bj = myByLo(myThid), myByHi(myThid) | 
| 480 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 481 |  |  | DO j=1-OLy,sNy+OLy | 
| 482 |  |  | DO i=1-OLx,sNx+OLx | 
| 483 |  |  | IF (STREAMICE_umask(i,j,bi,bj).eq.1.0) | 
| 484 |  |  | &      Zu_SI(i,j,bi,bj)=Ru_SI(i,j,bi,bj) / DIAGu_SI(i,j,bi,bj) | 
| 485 |  |  | IF (STREAMICE_vmask(i,j,bi,bj).eq.1.0) | 
| 486 |  |  | &      Zv_SI(i,j,bi,bj)=Rv_SI(i,j,bi,bj) / DIAGv_SI(i,j,bi,bj) | 
| 487 |  |  | ENDDO | 
| 488 |  |  | ENDDO | 
| 489 |  |  | ENDDO | 
| 490 |  |  | ENDDO | 
| 491 |  |  |  | 
| 492 |  |  | cg_halo = min(OLx-1,OLy-1) | 
| 493 |  |  | conv_flag = 0 | 
| 494 |  |  |  | 
| 495 |  |  | DO bj = myByLo(myThid), myByHi(myThid) | 
| 496 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 497 |  |  | DO j=1-OLy,sNy+OLy | 
| 498 |  |  | DO i=1-OLx,sNx+OLx | 
| 499 |  |  | Du_SI(i,j,bi,bj)=Zu_SI(i,j,bi,bj) | 
| 500 |  |  | Dv_SI(i,j,bi,bj)=Zv_SI(i,j,bi,bj) | 
| 501 |  |  | ENDDO | 
| 502 |  |  | ENDDO | 
| 503 |  |  | ENDDO | 
| 504 |  |  | ENDDO | 
| 505 |  |  |  | 
| 506 |  |  | resid = resid_0 | 
| 507 |  |  | iters = 0 | 
| 508 |  |  |  | 
| 509 |  |  | c  !!!!!!!!!!!!!!!!!! | 
| 510 |  |  | c  !!              !! | 
| 511 |  |  | c  !! MAIN CG LOOP !! | 
| 512 |  |  | c  !!              !! | 
| 513 |  |  | c  !!!!!!!!!!!!!!!!!! | 
| 514 | dgoldberg | 1.4 |  | 
| 515 |  |  |  | 
| 516 |  |  |  | 
| 517 | heimbach | 1.1 |  | 
| 518 |  |  | c  ! initially, b-grid data is valid up to 3 halo nodes out -- right? (check for MITgcm!!) | 
| 519 |  |  |  | 
| 520 | dgoldberg | 1.6 | WRITE(msgBuf,'(A)') 'BEGINNING MAIN CG LOOP' | 
| 521 |  |  | CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, | 
| 522 |  |  | &                     SQUEEZE_RIGHT , 1) | 
| 523 | heimbach | 1.1 |  | 
| 524 | dgoldberg | 1.3 | !       IF(STREAMICE_construct_matrix) CALL STREAMICE_CG_MAKE_A(myThid) | 
| 525 |  |  |  | 
| 526 | heimbach | 1.1 |  | 
| 527 |  |  | do iter = 1, streamice_max_cg_iter | 
| 528 |  |  | if (resid .gt. tolerance*resid_0) then | 
| 529 |  |  |  | 
| 530 |  |  | c      to avoid using "exit" | 
| 531 |  |  | iters = iters + 1 | 
| 532 |  |  |  | 
| 533 |  |  | is = 1 - cg_halo | 
| 534 |  |  | ie = sNx + cg_halo | 
| 535 |  |  | js = 1 - cg_halo | 
| 536 |  |  | je = sNy + cg_halo | 
| 537 |  |  |  | 
| 538 |  |  | DO bj = myByLo(myThid), myByHi(myThid) | 
| 539 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 540 |  |  | DO j=1-OLy,sNy+OLy | 
| 541 |  |  | DO i=1-OLx,sNx+OLx | 
| 542 |  |  | Au_SI(i,j,bi,bj) = 0. _d 0 | 
| 543 |  |  | Av_SI(i,j,bi,bj) = 0. _d 0 | 
| 544 |  |  | ENDDO | 
| 545 |  |  | ENDDO | 
| 546 |  |  | ENDDO | 
| 547 |  |  | ENDDO | 
| 548 |  |  |  | 
| 549 | dgoldberg | 1.3 | !        IF (STREAMICE_construct_matrix) THEN | 
| 550 |  |  |  | 
| 551 | dgoldberg | 1.4 | ! #ifdef STREAMICE_CONSTRUCT_MATRIX | 
| 552 | dgoldberg | 1.3 |  | 
| 553 | heimbach | 1.1 | DO bj = myByLo(myThid), myByHi(myThid) | 
| 554 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 555 |  |  | DO j=js,je | 
| 556 |  |  | DO i=is,ie | 
| 557 |  |  | DO colx=-1,1 | 
| 558 |  |  | DO coly=-1,1 | 
| 559 |  |  | Au_SI(i,j,bi,bj) = Au_SI(i,j,bi,bj) + | 
| 560 | dgoldberg | 1.4 | &         A_uu(i,j,bi,bj,colx,coly)* | 
| 561 | heimbach | 1.1 | &         Du_SI(i+colx,j+coly,bi,bj)+ | 
| 562 | dgoldberg | 1.4 | &         A_uv(i,j,bi,bj,colx,coly)* | 
| 563 | heimbach | 1.1 | &         Dv_SI(i+colx,j+coly,bi,bj) | 
| 564 |  |  | Av_SI(i,j,bi,bj) = Av_SI(i,j,bi,bj) + | 
| 565 | dgoldberg | 1.4 | &         A_vu(i,j,bi,bj,colx,coly)* | 
| 566 | heimbach | 1.1 | &         Du_SI(i+colx,j+coly,bi,bj)+ | 
| 567 | dgoldberg | 1.4 | &         A_vv(i,j,bi,bj,colx,coly)* | 
| 568 | heimbach | 1.1 | &         Dv_SI(i+colx,j+coly,bi,bj) | 
| 569 |  |  | ENDDO | 
| 570 |  |  | ENDDO | 
| 571 |  |  | ENDDO | 
| 572 |  |  | ENDDO | 
| 573 |  |  | ENDDO | 
| 574 |  |  | ENDDO | 
| 575 |  |  |  | 
| 576 | dgoldberg | 1.3 | !        else | 
| 577 | dgoldberg | 1.4 | ! #else | 
| 578 |  |  | ! | 
| 579 |  |  | !         CALL STREAMICE_CG_ACTION( myThid, | 
| 580 |  |  | !      O     Au_SI, | 
| 581 |  |  | !      O     Av_SI, | 
| 582 |  |  | !      I     Du_SI, | 
| 583 |  |  | !      I     Dv_SI, | 
| 584 |  |  | !      I     is,ie,js,je) | 
| 585 |  |  | ! | 
| 586 |  |  | ! !        ENDIF | 
| 587 |  |  | ! | 
| 588 |  |  | ! #endif | 
| 589 | heimbach | 1.1 |  | 
| 590 |  |  |  | 
| 591 |  |  | DO bj = myByLo(myThid), myByHi(myThid) | 
| 592 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 593 |  |  | dot_p1_tile(bi,bj) = 0. _d 0 | 
| 594 |  |  | dot_p2_tile(bi,bj) = 0. _d 0 | 
| 595 |  |  | ENDDO | 
| 596 |  |  | ENDDO | 
| 597 |  |  |  | 
| 598 |  |  | DO bj = myByLo(myThid), myByHi(myThid) | 
| 599 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 600 |  |  | DO j=1,sNy | 
| 601 |  |  | DO i=1,sNx | 
| 602 |  |  | IF (STREAMICE_umask(i,j,bi,bj).eq.1.0) THEN | 
| 603 |  |  | dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Zu_SI(i,j,bi,bj)* | 
| 604 |  |  | &        Ru_SI(i,j,bi,bj) | 
| 605 |  |  | dot_p2_tile(bi,bj)=dot_p2_tile(bi,bj)+Du_SI(i,j,bi,bj)* | 
| 606 |  |  | &        Au_SI(i,j,bi,bj) | 
| 607 |  |  | ENDIF | 
| 608 |  |  | IF (STREAMICE_vmask(i,j,bi,bj).eq.1.0) THEN | 
| 609 |  |  | dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Zv_SI(i,j,bi,bj)* | 
| 610 |  |  | &        Rv_SI(i,j,bi,bj) | 
| 611 |  |  | dot_p2_tile(bi,bj)=dot_p2_tile(bi,bj)+Dv_SI(i,j,bi,bj)* | 
| 612 |  |  | &        Av_SI(i,j,bi,bj) | 
| 613 |  |  | ENDIF | 
| 614 |  |  | ENDDO | 
| 615 |  |  | ENDDO | 
| 616 |  |  | ENDDO | 
| 617 |  |  | ENDDO | 
| 618 |  |  |  | 
| 619 |  |  | CALL GLOBAL_SUM_TILE_RL( dot_p1_tile, dot_p1, myThid ) | 
| 620 |  |  | CALL GLOBAL_SUM_TILE_RL( dot_p2_tile, dot_p2, myThid ) | 
| 621 |  |  | alpha_k = dot_p1/dot_p2 | 
| 622 |  |  |  | 
| 623 |  |  | DO bj = myByLo(myThid), myByHi(myThid) | 
| 624 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 625 |  |  | DO j=1-OLy,sNy+OLy | 
| 626 |  |  | DO i=1-OLx,sNx+OLx | 
| 627 |  |  |  | 
| 628 |  |  | IF (STREAMICE_umask(i,j,bi,bj).eq.1.0) THEN | 
| 629 |  |  | cg_Uin(i,j,bi,bj)=cg_Uin(i,j,bi,bj)+ | 
| 630 |  |  | &       alpha_k*Du_SI(i,j,bi,bj) | 
| 631 |  |  | Ru_old_SI(i,j,bi,bj) = Ru_SI(i,j,bi,bj) | 
| 632 |  |  | Zu_old_SI(i,j,bi,bj) = Zu_SI(i,j,bi,bj) | 
| 633 |  |  | Ru_SI(i,j,bi,bj) = Ru_SI(i,j,bi,bj)- | 
| 634 |  |  | &       alpha_k*Au_SI(i,j,bi,bj) | 
| 635 |  |  | Zu_SI(i,j,bi,bj) = Ru_SI(i,j,bi,bj) / | 
| 636 |  |  | &       DIAGu_SI(i,j,bi,bj) | 
| 637 |  |  | ENDIF | 
| 638 |  |  |  | 
| 639 |  |  | IF (STREAMICE_vmask(i,j,bi,bj).eq.1.0) THEN | 
| 640 |  |  | cg_Vin(i,j,bi,bj)=cg_Vin(i,j,bi,bj)+ | 
| 641 |  |  | &       alpha_k*Dv_SI(i,j,bi,bj) | 
| 642 |  |  | Rv_old_SI(i,j,bi,bj) = Rv_SI(i,j,bi,bj) | 
| 643 |  |  | Zv_old_SI(i,j,bi,bj) = Zv_SI(i,j,bi,bj) | 
| 644 |  |  | Rv_SI(i,j,bi,bj) = Rv_SI(i,j,bi,bj)- | 
| 645 |  |  | &       alpha_k*Av_SI(i,j,bi,bj) | 
| 646 |  |  | Zv_SI(i,j,bi,bj) = Rv_SI(i,j,bi,bj) / | 
| 647 |  |  | &       DIAGv_SI(i,j,bi,bj) | 
| 648 |  |  |  | 
| 649 |  |  | ENDIF | 
| 650 |  |  | ENDDO | 
| 651 |  |  | ENDDO | 
| 652 |  |  | ENDDO | 
| 653 |  |  | ENDDO | 
| 654 |  |  |  | 
| 655 |  |  | DO bj = myByLo(myThid), myByHi(myThid) | 
| 656 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 657 |  |  | dot_p1_tile(bi,bj) = 0. _d 0 | 
| 658 |  |  | dot_p2_tile(bi,bj) = 0. _d 0 | 
| 659 |  |  | ENDDO | 
| 660 |  |  | ENDDO | 
| 661 |  |  |  | 
| 662 |  |  | DO bj = myByLo(myThid), myByHi(myThid) | 
| 663 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 664 |  |  | DO j=1,sNy | 
| 665 |  |  | DO i=1,sNx | 
| 666 |  |  |  | 
| 667 |  |  | IF (STREAMICE_umask(i,j,bi,bj).eq.1.0) THEN | 
| 668 |  |  | dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Zu_SI(i,j,bi,bj)* | 
| 669 |  |  | &        Ru_SI(i,j,bi,bj) | 
| 670 |  |  | dot_p2_tile(bi,bj)=dot_p2_tile(bi,bj)+Zu_old_SI(i,j,bi,bj)* | 
| 671 |  |  | &        Ru_old_SI(i,j,bi,bj) | 
| 672 |  |  | ENDIF | 
| 673 |  |  |  | 
| 674 |  |  | IF (STREAMICE_vmask(i,j,bi,bj).eq.1.0) THEN | 
| 675 |  |  | dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Zv_SI(i,j,bi,bj)* | 
| 676 |  |  | &        Rv_SI(i,j,bi,bj) | 
| 677 |  |  | dot_p2_tile(bi,bj)=dot_p2_tile(bi,bj)+Zv_old_SI(i,j,bi,bj)* | 
| 678 |  |  | &        Rv_old_SI(i,j,bi,bj) | 
| 679 |  |  | ENDIF | 
| 680 |  |  |  | 
| 681 |  |  | ENDDO | 
| 682 |  |  | ENDDO | 
| 683 |  |  | ENDDO | 
| 684 |  |  | ENDDO | 
| 685 |  |  |  | 
| 686 |  |  | CALL GLOBAL_SUM_TILE_RL( dot_p1_tile, dot_p1, myThid ) | 
| 687 |  |  | CALL GLOBAL_SUM_TILE_RL( dot_p2_tile, dot_p2, myThid ) | 
| 688 |  |  |  | 
| 689 |  |  | beta_k = dot_p1/dot_p2 | 
| 690 |  |  |  | 
| 691 |  |  | DO bj = myByLo(myThid), myByHi(myThid) | 
| 692 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 693 |  |  | DO j=1-OLy,sNy+OLy | 
| 694 |  |  | DO i=1-OLx,sNx+OLx | 
| 695 |  |  | IF (STREAMICE_umask(i,j,bi,bj).eq.1.0) | 
| 696 |  |  | &      Du_SI(i,j,bi,bj)=beta_k*Du_SI(i,j,bi,bj)+ | 
| 697 |  |  | &      Zu_SI(i,j,bi,bj) | 
| 698 |  |  | IF (STREAMICE_vmask(i,j,bi,bj).eq.1.0) | 
| 699 |  |  | &      Dv_SI(i,j,bi,bj)=beta_k*Dv_SI(i,j,bi,bj)+ | 
| 700 |  |  | &      Zv_SI(i,j,bi,bj) | 
| 701 |  |  | ENDDO | 
| 702 |  |  | ENDDO | 
| 703 |  |  | ENDDO | 
| 704 |  |  | ENDDO | 
| 705 |  |  |  | 
| 706 |  |  | DO bj = myByLo(myThid), myByHi(myThid) | 
| 707 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 708 |  |  | dot_p1_tile(bi,bj) = 0. _d 0 | 
| 709 |  |  | ENDDO | 
| 710 |  |  | ENDDO | 
| 711 |  |  |  | 
| 712 |  |  | DO bj = myByLo(myThid), myByHi(myThid) | 
| 713 |  |  | DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 714 |  |  | DO j=1,sNy | 
| 715 |  |  | DO i=1,sNx | 
| 716 |  |  | IF (STREAMICE_umask(i,j,bi,bj).eq.1.0) | 
| 717 |  |  | &      dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Ru_SI(i,j,bi,bj)**2 | 
| 718 |  |  | IF (STREAMICE_vmask(i,j,bi,bj).eq.1.0) | 
| 719 |  |  | &      dot_p1_tile(bi,bj)=dot_p1_tile(bi,bj)+Rv_SI(i,j,bi,bj)**2 | 
| 720 |  |  | ENDDO | 
| 721 |  |  | ENDDO | 
| 722 |  |  | ENDDO | 
| 723 |  |  | ENDDO | 
| 724 |  |  |  | 
| 725 |  |  | CALL GLOBAL_SUM_TILE_RL( dot_p1_tile, dot_p1, myThid ) | 
| 726 |  |  | resid = sqrt(dot_p1) | 
| 727 |  |  |  | 
| 728 |  |  | !        IF (iter .eq. 1) then | 
| 729 |  |  | !         print *, alpha_k, beta_k, resid | 
| 730 |  |  | !        ENDIF | 
| 731 |  |  |  | 
| 732 |  |  | cg_halo = cg_halo - 1 | 
| 733 |  |  |  | 
| 734 |  |  | if (cg_halo .eq. 0) then | 
| 735 |  |  | cg_halo = min(OLx-1,OLy-1) | 
| 736 |  |  | _EXCH_XY_RL( Du_SI, myThid ) | 
| 737 |  |  | _EXCH_XY_RL( Dv_SI, myThid ) | 
| 738 |  |  | _EXCH_XY_RL( Ru_SI, myThid ) | 
| 739 |  |  | _EXCH_XY_RL( Rv_SI, myThid ) | 
| 740 |  |  | _EXCH_XY_RL( cg_Uin, myThid ) | 
| 741 |  |  | _EXCH_XY_RL( cg_Vin, myThid ) | 
| 742 |  |  | endif | 
| 743 |  |  |  | 
| 744 | dgoldberg | 1.4 |  | 
| 745 | heimbach | 1.1 | endif | 
| 746 |  |  | enddo ! end of CG loop | 
| 747 |  |  |  | 
| 748 |  |  | c     to avoid using "exit" | 
| 749 |  |  | c     if iters has reached max_iters there is no convergence | 
| 750 |  |  |  | 
| 751 |  |  | IF (iters .lt. streamice_max_cg_iter) THEN | 
| 752 |  |  | conv_flag = 1 | 
| 753 |  |  | ENDIF | 
| 754 |  |  |  | 
| 755 | dgoldberg | 1.4 | !       DO bj = myByLo(myThid), myByHi(myThid) | 
| 756 |  |  | !        DO bi = myBxLo(myThid), myBxHi(myThid) | 
| 757 |  |  | !         DO j=1-OLy,sNy+OLy | 
| 758 |  |  | !          DO i=1-OLy,sNx+OLy | 
| 759 |  |  | !           IF (STREAMICE_umask(i,j,bi,bj).eq.3.0) | 
| 760 |  |  | !      &     cg_Uin(i,j,bi,bj)=u_bdry_values_SI(i,j,bi,bj) | 
| 761 |  |  | !           IF (STREAMICE_vmask(i,j,bi,bj).eq.3.0) | 
| 762 |  |  | !      &     cg_Vin(i,j,bi,bj)=v_bdry_values_SI(i,j,bi,bj) | 
| 763 |  |  | !          ENDDO | 
| 764 |  |  | !         ENDDO | 
| 765 |  |  | !        ENDDO | 
| 766 |  |  | !       ENDDO | 
| 767 |  |  | ! | 
| 768 |  |  | !       _EXCH_XY_RL( cg_Uin, myThid ) | 
| 769 |  |  | !       _EXCH_XY_RL( cg_Vin, myThid ) | 
| 770 | heimbach | 1.1 |  | 
| 771 | dgoldberg | 1.8 | #endif /* ifndef ALLOW_PETSC */ | 
| 772 |  |  |  | 
| 773 |  |  | #else /* STREAMICE_SERIAL_TRISOLVE */ | 
| 774 |  |  |  | 
| 775 |  |  | CALL STREAMICE_TRIDIAG_SOLVE( | 
| 776 |  |  | U                               cg_Uin,     ! x-velocities | 
| 777 |  |  | U                               cg_Vin, | 
| 778 |  |  | U                               cg_Bu,      ! force in x dir | 
| 779 |  |  | I                               A_uu,       ! section of matrix that multiplies u and projects on u | 
| 780 |  |  | I                               STREAMICE_umask, | 
| 781 |  |  | I                               myThid ) | 
| 782 |  |  |  | 
| 783 | heimbach | 1.1 | #endif | 
| 784 |  |  |  | 
| 785 | dgoldberg | 1.7 | CALL TIMER_STOP ('STREAMICE_CG_SOLVE',myThid) | 
| 786 | heimbach | 1.1 |  | 
| 787 |  |  |  | 
| 788 | dgoldberg | 1.7 | #endif | 
| 789 |  |  | RETURN | 
| 790 |  |  | END |