| 1 |
C $Header: /u/gcmpack/MITgcm/pkg/streamice/streamice_cg_solve_petsc.F,v 1.1 2013/08/24 20:32:03 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_PETSC( |
| 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, maxiter |
| 54 |
_RL tolerance |
| 55 |
_RL cg_Uin (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
| 56 |
_RL cg_Vin (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
| 57 |
_RL cg_Bu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
| 58 |
_RL cg_Bv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
| 59 |
_RL |
| 60 |
& A_uu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1), |
| 61 |
& A_vu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1), |
| 62 |
& A_uv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1), |
| 63 |
& A_vv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1) |
| 64 |
|
| 65 |
C LOCAL VARIABLES |
| 66 |
INTEGER i, j, bi, bj, cg_halo, conv_flag |
| 67 |
INTEGER iter, is, js, ie, je, colx, coly, k |
| 68 |
_RL dot_p1, dot_p2, alpha_k, beta_k, resid, resid_0 |
| 69 |
_RL dot_p1_tile (nSx,nSy) |
| 70 |
_RL dot_p2_tile (nSx,nSy) |
| 71 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
| 72 |
|
| 73 |
|
| 74 |
#ifdef ALLOW_PETSC |
| 75 |
INTEGER indices(2*(snx*nsx*sny*nsy)) |
| 76 |
INTEGER n_dofs_cum_sum (0:nPx*nPy-1), idx(1) |
| 77 |
_RL rhs_values(2*(snx*nsx*sny*nsy)) |
| 78 |
_RL solution_values(2*(snx*nsx*sny*nsy)) |
| 79 |
! _RL mat_values (2*Nx*Ny,2*(snx*nsx*sny*nsy)) |
| 80 |
_RL mat_values (18,1), mat_val_return(1) |
| 81 |
INTEGER indices_col(18) |
| 82 |
INTEGER local_dofs, global_dofs, dof_index, dof_index_col |
| 83 |
INTEGER local_offset |
| 84 |
Mat matrix |
| 85 |
KSP ksp |
| 86 |
PC pc |
| 87 |
Vec rhs |
| 88 |
Vec solution |
| 89 |
PetscErrorCode ierr |
| 90 |
#ifdef ALLOW_USE_MPI |
| 91 |
integer mpiRC, mpiMyWid |
| 92 |
#endif |
| 93 |
#endif |
| 94 |
|
| 95 |
|
| 96 |
#ifdef ALLOW_STREAMICE |
| 97 |
|
| 98 |
|
| 99 |
|
| 100 |
|
| 101 |
#ifdef ALLOW_PETSC |
| 102 |
|
| 103 |
#ifdef ALLOW_USE_MPI |
| 104 |
|
| 105 |
|
| 106 |
CALL MPI_COMM_RANK( MPI_COMM_WORLD, mpiMyWId, mpiRC ) |
| 107 |
local_dofs = n_dofs_process (mpiMyWid) |
| 108 |
global_dofs = 0 |
| 109 |
|
| 110 |
n_dofs_cum_sum(0) = 0 |
| 111 |
DO i=0,nPx*nPy-1 |
| 112 |
global_dofs = global_dofs + n_dofs_process (i) |
| 113 |
if (i.ge.1) THEN |
| 114 |
n_dofs_cum_sum(i) = n_dofs_cum_sum(i-1)+ |
| 115 |
& n_dofs_process(i-1) |
| 116 |
endif |
| 117 |
ENDDO |
| 118 |
local_offset = n_dofs_cum_sum(mpimywid) |
| 119 |
|
| 120 |
#else |
| 121 |
|
| 122 |
local_dofs = n_dofs_process (0) |
| 123 |
global_dofs = local_dofs |
| 124 |
local_offset = 0 |
| 125 |
|
| 126 |
#endif |
| 127 |
|
| 128 |
! call petscInitialize(PETSC_NULL_CHARACTER,ierr) |
| 129 |
|
| 130 |
!---------------------- |
| 131 |
|
| 132 |
call VecCreate(PETSC_COMM_WORLD, rhs, ierr) |
| 133 |
call VecSetSizes(rhs, local_dofs, global_dofs, ierr) |
| 134 |
call VecSetType(rhs, VECMPI, ierr) |
| 135 |
|
| 136 |
call VecCreate(PETSC_COMM_WORLD, solution, ierr) |
| 137 |
call VecSetSizes(solution, local_dofs, global_dofs, ierr) |
| 138 |
call VecSetType(solution, VECMPI, ierr) |
| 139 |
|
| 140 |
do i=1,local_dofs |
| 141 |
indices(i) = i-1 + local_offset |
| 142 |
end do |
| 143 |
do i=1,2*nSx*nSy*sNx*sNy |
| 144 |
rhs_values (i) = 0. _d 0 |
| 145 |
solution_values (i) = 0. _d 0 |
| 146 |
enddo |
| 147 |
|
| 148 |
! gather rhs and initial guess values to populate petsc vectors |
| 149 |
|
| 150 |
DO bj = myByLo(myThid), myByHi(myThid) |
| 151 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
| 152 |
DO j=1,sNy |
| 153 |
DO i=1,sNx |
| 154 |
|
| 155 |
dof_index = INT(streamice_petsc_dofs_u(i,j,bi,bj)) |
| 156 |
& - local_offset |
| 157 |
|
| 158 |
if (dof_index.ge.0) THEN |
| 159 |
|
| 160 |
rhs_values(dof_index+1) = cg_Bu(i,j,bi,bj) |
| 161 |
solution_values(dof_index+1) = cg_Uin(i,j,bi,bj) |
| 162 |
|
| 163 |
endif |
| 164 |
|
| 165 |
!--------------- |
| 166 |
|
| 167 |
dof_index = INT(streamice_petsc_dofs_v(i,j,bi,bj)) |
| 168 |
& - local_offset |
| 169 |
|
| 170 |
if (dof_index.ge.0) THEN |
| 171 |
|
| 172 |
rhs_values(dof_index+1) = cg_Bv(i,j,bi,bj) |
| 173 |
solution_values(dof_index+1) = cg_Vin(i,j,bi,bj) |
| 174 |
|
| 175 |
endif |
| 176 |
|
| 177 |
ENDDO |
| 178 |
ENDDO |
| 179 |
ENDDO |
| 180 |
ENDDO |
| 181 |
|
| 182 |
|
| 183 |
call VecSetValues(rhs, local_dofs, indices, rhs_values, |
| 184 |
& INSERT_VALUES, ierr) |
| 185 |
call VecAssemblyBegin(rhs, ierr) |
| 186 |
call VecAssemblyEnd(rhs, ierr) |
| 187 |
|
| 188 |
|
| 189 |
call VecSetValues(solution, local_dofs, indices, |
| 190 |
& solution_values, INSERT_VALUES, ierr) |
| 191 |
call VecAssemblyBegin(solution, ierr) |
| 192 |
call VecAssemblyEnd(solution, ierr) |
| 193 |
|
| 194 |
! IF USING v3.0 THEN |
| 195 |
! call MatCreateMPIAIJ (PETSC_COMM_WORLD, |
| 196 |
call MatCreateAIJ (PETSC_COMM_WORLD, |
| 197 |
& local_dofs, local_dofs, |
| 198 |
& global_dofs, global_dofs, |
| 199 |
& 18, PETSC_NULL_INTEGER, |
| 200 |
& 18, PETSC_NULL_INTEGER, |
| 201 |
& matrix, ierr) |
| 202 |
|
| 203 |
|
| 204 |
! populate petsc matrix |
| 205 |
|
| 206 |
DO bj = myByLo(myThid), myByHi(myThid) |
| 207 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
| 208 |
DO j=1,sNy |
| 209 |
DO i=1,sNx |
| 210 |
|
| 211 |
dof_index = INT(streamice_petsc_dofs_u(i,j,bi,bj)) |
| 212 |
! & - local_offset |
| 213 |
|
| 214 |
IF (dof_index .ge. 0) THEN |
| 215 |
|
| 216 |
DO k=1,18 |
| 217 |
indices_col(k) = 0 |
| 218 |
mat_values(k,1) = 0. _d 0 |
| 219 |
ENDDO |
| 220 |
k=0 |
| 221 |
|
| 222 |
DO coly=-1,1 |
| 223 |
DO colx=-1,1 |
| 224 |
|
| 225 |
dof_index_col = streamice_petsc_dofs_u(i+colx,j+coly,bi,bj) |
| 226 |
|
| 227 |
if (dof_index_col.ge.0) THEN |
| 228 |
! pscal = A_uu(i,j,bi,bj,colx,coly) |
| 229 |
! CALL MatSetValue (matrix,dof_index, dof_index_col, |
| 230 |
! & pscal,INSERT_VALUES,ierr) |
| 231 |
k=k+1 |
| 232 |
mat_values (k,1) = A_uu(i,j,bi,bj,colx,coly) |
| 233 |
indices_col (k) = dof_index_col |
| 234 |
endif |
| 235 |
|
| 236 |
dof_index_col = streamice_petsc_dofs_v(i+colx,j+coly,bi,bj) |
| 237 |
|
| 238 |
if (dof_index_col.ge.0) THEN |
| 239 |
! CALL MatSetValue (matrix,dof_index, dof_index_col, |
| 240 |
! & A_uv(i,j,bi,bj,colx,coly),INSERT_VALUES,ierr) |
| 241 |
k=k+1 |
| 242 |
mat_values (k,1) = A_uv(i,j,bi,bj,colx,coly) |
| 243 |
indices_col (k) = dof_index_col |
| 244 |
endif |
| 245 |
|
| 246 |
ENDDO |
| 247 |
ENDDO |
| 248 |
|
| 249 |
call matSetValues (matrix, 1, dof_index, k, indices_col, |
| 250 |
& mat_values,INSERT_VALUES,ierr) |
| 251 |
|
| 252 |
|
| 253 |
ENDIF |
| 254 |
|
| 255 |
! ---------------------------------------------- |
| 256 |
|
| 257 |
dof_index = INT(streamice_petsc_dofs_v(i,j,bi,bj)) |
| 258 |
! & - local_offset |
| 259 |
|
| 260 |
IF (dof_index .ge. 0) THEN |
| 261 |
|
| 262 |
DO k=1,18 |
| 263 |
indices_col(k) = 0 |
| 264 |
mat_values(k,1) = 0. _d 0 |
| 265 |
ENDDO |
| 266 |
k=0 |
| 267 |
|
| 268 |
DO coly=-1,1 |
| 269 |
DO colx=-1,1 |
| 270 |
|
| 271 |
dof_index_col = streamice_petsc_dofs_u(i+colx,j+coly,bi,bj) |
| 272 |
|
| 273 |
if (dof_index_col.ge.0) THEN |
| 274 |
! CALL MatSetValue (matrix,dof_index, dof_index_col, |
| 275 |
! & A_vu(i,j,bi,bj,colx,coly),INSERT_VALUES,ierr) |
| 276 |
k=k+1 |
| 277 |
mat_values (k,1) = A_vu(i,j,bi,bj,colx,coly) |
| 278 |
indices_col (k) = dof_index_col |
| 279 |
endif |
| 280 |
|
| 281 |
dof_index_col = streamice_petsc_dofs_v(i+colx,j+coly,bi,bj) |
| 282 |
|
| 283 |
if (dof_index_col.ge.0) THEN |
| 284 |
! CALL MatSetValue (matrix,dof_index, dof_index_col, |
| 285 |
! & A_vv(i,j,bi,bj,colx,coly),INSERT_VALUES,ierr) |
| 286 |
k=k+1 |
| 287 |
mat_values (k,1) = A_vv(i,j,bi,bj,colx,coly) |
| 288 |
indices_col (k) = dof_index_col |
| 289 |
endif |
| 290 |
|
| 291 |
ENDDO |
| 292 |
ENDDO |
| 293 |
|
| 294 |
call matSetValues (matrix, 1, dof_index, k, indices_col, |
| 295 |
& mat_values,INSERT_VALUES,ierr) |
| 296 |
ENDIF |
| 297 |
|
| 298 |
ENDDO |
| 299 |
ENDDO |
| 300 |
ENDDO |
| 301 |
ENDDO |
| 302 |
|
| 303 |
call MatAssemblyBegin(matrix,MAT_FINAL_ASSEMBLY,ierr) |
| 304 |
call MatAssemblyEnd(matrix,MAT_FINAL_ASSEMBLY,ierr) |
| 305 |
|
| 306 |
|
| 307 |
call KSPCreate(PETSC_COMM_WORLD, ksp, ierr) |
| 308 |
call KSPSetOperators(ksp, matrix, matrix, |
| 309 |
& DIFFERENT_NONZERO_PATTERN, ierr) |
| 310 |
|
| 311 |
SELECT CASE (PETSC_SOLVER_TYPE) |
| 312 |
CASE ('CG') |
| 313 |
PRINT *, "PETSC SOLVER: SELECTED CG" |
| 314 |
call KSPSetType(ksp, KSPCG, ierr) |
| 315 |
CASE ('GMRES') |
| 316 |
PRINT *, "PETSC SOLVER: SELECTED GMRES" |
| 317 |
call KSPSetType(ksp, KSPGMRES, ierr) |
| 318 |
CASE ('BICG') |
| 319 |
PRINT *, "PETSC SOLVER: SELECTED BICG" |
| 320 |
call KSPSetType(ksp, KSPBICG, ierr) |
| 321 |
CASE DEFAULT |
| 322 |
PRINT *, "PETSC SOLVER: SELECTED DEFAULT" |
| 323 |
call KSPSetType(ksp, KSPCG, ierr) |
| 324 |
END SELECT |
| 325 |
|
| 326 |
call KSPGetPC(ksp, pc, ierr) |
| 327 |
call KSPSetTolerances(ksp,tolerance, |
| 328 |
& PETSC_DEFAULT_DOUBLE_PRECISION, |
| 329 |
& PETSC_DEFAULT_DOUBLE_PRECISION, |
| 330 |
& maxiter,ierr) |
| 331 |
|
| 332 |
SELECT CASE (PETSC_PRECOND_TYPE) |
| 333 |
CASE ('BLOCKJACOBI') |
| 334 |
PRINT *, "PETSC PRECOND: SELECTED BJACOBI" |
| 335 |
call PCSetType(pc, PCBJACOBI, ierr) |
| 336 |
CASE ('JACOBI') |
| 337 |
PRINT *, "PETSC PRECOND: SELECTED JACOBI" |
| 338 |
call PCSetType(pc, PCJACOBI, ierr) |
| 339 |
CASE ('ILU') |
| 340 |
PRINT *, "PETSC PRECOND: SELECTED ILU" |
| 341 |
call PCSetType(pc, PCILU, ierr) |
| 342 |
CASE DEFAULT |
| 343 |
PRINT *, "PETSC PRECOND: SELECTED DEFAULT" |
| 344 |
call PCSetType(pc, PCBJACOBI, ierr) |
| 345 |
END SELECT |
| 346 |
|
| 347 |
|
| 348 |
call KSPSolve(ksp, rhs, solution, ierr) |
| 349 |
call KSPGetIterationNumber(ksp,iters,ierr) |
| 350 |
|
| 351 |
call VecGetValues(solution,local_dofs,indices, |
| 352 |
& solution_values,ierr) |
| 353 |
|
| 354 |
DO bj = myByLo(myThid), myByHi(myThid) |
| 355 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
| 356 |
DO j=1,sNy |
| 357 |
DO i=1,sNx |
| 358 |
|
| 359 |
dof_index = INT(streamice_petsc_dofs_u(i,j,bi,bj)) |
| 360 |
& - local_offset |
| 361 |
if (dof_index.ge.0) THEN |
| 362 |
cg_Uin(i,j,bi,bj) = solution_values(dof_index+1) |
| 363 |
endif |
| 364 |
|
| 365 |
dof_index = INT(streamice_petsc_dofs_v(i,j,bi,bj)) |
| 366 |
& - local_offset |
| 367 |
if (dof_index.ge.0) THEN |
| 368 |
cg_Vin(i,j,bi,bj) = solution_values(dof_index+1) |
| 369 |
endif |
| 370 |
|
| 371 |
ENDDO |
| 372 |
ENDDO |
| 373 |
ENDDO |
| 374 |
ENDDO |
| 375 |
|
| 376 |
call KSPDestroy (ksp, ierr) |
| 377 |
call VecDestroy (rhs, ierr) |
| 378 |
call VecDestroy (solution, ierr) |
| 379 |
call MatDestroy (matrix, ierr) |
| 380 |
|
| 381 |
|
| 382 |
! call PetscFinalize(ierr) |
| 383 |
|
| 384 |
|
| 385 |
#endif |
| 386 |
|
| 387 |
|
| 388 |
|
| 389 |
#endif |
| 390 |
RETURN |
| 391 |
END |