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