Back to home page

MITgcm

 
 

    


File indexing completed on 2026-02-13 06:08:56 UTC

view on githubraw file Latest commit decd05a6 on 2026-02-12 15:55:36 UTC
feb7fa5d1e dngo*0001 #include "TAPENADE_OPTIONS.h"
                0002 #ifdef ALLOW_STREAMICE
                0003 # include "STREAMICE_OPTIONS.h"
6a9e386e2e dngo*0004 #endif
b4daa24319 Shre*0005 
6a9e386e2e dngo*0006 C--   File stubs_tap_adj.F: Subroutines with "manual" derivatives
                0007 C--  to be called by the tapenade adjoint and tangent linear code
                0008 C--    *_B: adjoint with checkpointing on
                0009 C--    *_D: tlm
                0010 C--    *_FWD, *_BWD: adjoint with no checkpointing
                0011 C--
                0012 C--    Contents:
                0013 C--    o  GLOBAL_MAX_R8_B
                0014 C--    o  GLOBAL_SUM_TILE_RL_B
                0015 C--    o  GLOBAL_SUM_TILE_RL_FWD
                0016 C--    o  GLOBAL_SUM_TILE_RL_BWD
                0017 C--    o  GLOBAL_SUM_R8_B
                0018 C--    o  CG2D_B0
                0019 C--    o  ADEXCH_3D_RL
                0020 C--    o  ADEXCH_UV_XY_RS
                0021 C--    o  ADEXCH_UV_3D_RL
                0022 C--    o  ADEXCH_XY_RS
                0023 C--    o  ADEXCH_XY_RL
                0024 C--    o  STREAMICE_CG_SOLVE_D
                0025 C--    o  STREAMICE_CG_SOLVE_B
                0026 C--    o  STREAMICE_CG_SOLVE_FWD
                0027 C--    o  STREAMICE_CG_SOLVE_BWD
                0028 
                0029 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
b4daa24319 Shre*0030       SUBROUTINE GLOBAL_MAX_R8_B(rhsmax, rhsmaxb, myThid)
                0031       IMPLICIT NONE
                0032 #include "SIZE.h"
                0033 #include "EEPARAMS.h"
                0034 #include "EESUPPORT.h"
                0035 #include "GLOBAL_MAX.h"
642de41482 dngo*0036       Real*8  rhsmax
                0037       Real*8  rhsmaxb
b4daa24319 Shre*0038       INTEGER myThid
                0039 
                0040       CALL GLOBAL_ADMAX_R8(rhsmaxb, myThid)
                0041       RETURN
                0042       END
                0043 
6a9e386e2e dngo*0044 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
b4daa24319 Shre*0045       SUBROUTINE GLOBAL_SUM_TILE_RL_B(phiTile, phiTileb,
642de41482 dngo*0046      &                                sumPhi, sumPhib, myThid)
b4daa24319 Shre*0047       IMPLICIT NONE
                0048 #include "SIZE.h"
                0049 #include "EEPARAMS.h"
                0050 #include "EESUPPORT.h"
                0051 #include "GLOBAL_SUM.h"
                0052       _RL     phiTile(nSx,nSy)
                0053       _RL     phiTileb(nSx,nSy)
                0054       _RL     sumPhib
                0055       _RL     sumPhi
                0056       INTEGER myThid
                0057 
                0058       CALL GLOBAL_ADSUM_TILE_RL(phiTileb, sumPhib, myThid)
642de41482 dngo*0059       RETURN
b4daa24319 Shre*0060       END
                0061 
6a9e386e2e dngo*0062 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0063       SUBROUTINE GLOBAL_SUM_TILE_RL_FWD(phiTile,
642de41482 dngo*0064      &                                  sumPhi, myThid)
6a9e386e2e dngo*0065       IMPLICIT NONE
                0066 #include "SIZE.h"
                0067 #include "EEPARAMS.h"
                0068 #include "EESUPPORT.h"
                0069 #include "GLOBAL_SUM.h"
                0070       _RL     phiTile(nSx,nSy)
                0071       _RL     sumPhi
                0072       INTEGER myThid
                0073 
                0074       CALL PUSHINTEGER4(myThid)
                0075       CALL GLOBAL_SUM_TILE_RL(phiTile, sumPhi, myThid)
                0076       CALL PUSHREAL8(sumPhi)
642de41482 dngo*0077       RETURN
6a9e386e2e dngo*0078       END
                0079 
                0080 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0081       SUBROUTINE GLOBAL_SUM_TILE_RL_BWD(phiTile, phiTileb,
642de41482 dngo*0082      &                                  sumPhi, sumPhib, myThid)
6a9e386e2e dngo*0083       IMPLICIT NONE
                0084 #include "SIZE.h"
                0085 #include "EEPARAMS.h"
                0086 #include "EESUPPORT.h"
                0087 #include "GLOBAL_SUM.h"
                0088       _RL     phiTile(nSx,nSy)
                0089       _RL     phiTileb(nSx,nSy)
                0090       _RL     sumPhib
                0091       _RL     sumPhi
                0092       INTEGER myThid
                0093 
                0094       CALL POPREAL8(sumPhi)
                0095       CALL POPINTEGER4(myThid)
                0096       CALL GLOBAL_ADSUM_TILE_RL(phiTileb, sumPhib, myThid)
642de41482 dngo*0097       RETURN
6a9e386e2e dngo*0098       END
                0099 
                0100 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
b4daa24319 Shre*0101 !     Adjoint of global_sum_r8 for arguments Arg1=(in;out)
                0102       SUBROUTINE GLOBAL_SUM_R8_B(sumPhi, sumPhib, myThid)
                0103       IMPLICIT NONE
                0104 #include "SIZE.h"
                0105 #include "EEPARAMS.h"
                0106 #include "EESUPPORT.h"
                0107 #include "GLOBAL_SUM.h"
642de41482 dngo*0108       Real*8 sumPhi
                0109       Real*8 sumPhib
b4daa24319 Shre*0110       INTEGER myThid
                0111 
                0112       CALL GLOBAL_ADSUM_R8(sumPhib, myThid)
642de41482 dngo*0113       RETURN
b4daa24319 Shre*0114       END
                0115 
6a9e386e2e dngo*0116 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
642de41482 dngo*0117       SUBROUTINE CG2D_B0(cg2d_b, cg2d_bb, cg2d_x, cg2d_xb,
                0118      &                   firstResidual, minResidualSq, lastResidual,
                0119      &                   numIters, nIterMin, myThid)
b4daa24319 Shre*0120       IMPLICIT NONE
                0121 #include "SIZE.h"
                0122 #include "EEPARAMS.h"
                0123 #include "PARAMS.h"
                0124 #include "CG2D.h"
                0125       _RL  cg2d_b(1-OLx:sNx+OLx, 1-OLy:sNy+OLy, nSx, nSy)
                0126       _RL  cg2d_bb(1-OLx:sNx+OLx, 1-OLy:sNy+OLy, nSx, nSy)
                0127       _RL  cg2d_x(1-OLx:sNx+OLx, 1-OLy:sNy+OLy, nSx, nSy)
                0128       _RL  cg2d_xb(1-OLx:sNx+OLx, 1-OLy:sNy+OLy, nSx, nSy)
                0129       _RL  firstResidual
                0130       _RL  minResidualSq
                0131       _RL  lastResidual
                0132       INTEGER numIters
                0133       INTEGER nIterMin
                0134       INTEGER myThid
                0135 
                0136 ! [llh] we assume the downstream cg2d_b is passive, which helps us
                0137 !    because it seems the input 2nd arg of cg2d() pollutes its output value:
642de41482 dngo*0138       cg2d_bb = 0. _d 0
                0139       CALL CG2D( cg2d_xb, cg2d_bb,
                0140      &           firstResidual, minResidualSq, lastResidual,
                0141      &           numIters, nIterMin, myThid)
b4daa24319 Shre*0142 ! [llh] the upstream cg2d_x is passive:
642de41482 dngo*0143       cg2d_xb = 0. _d 0
                0144       RETURN
b4daa24319 Shre*0145       END
                0146 
6a9e386e2e dngo*0147 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
b4daa24319 Shre*0148       SUBROUTINE ADEXCH_3D_RL( adVar, Nr, myThid )
                0149       IMPLICIT NONE
642de41482 dngo*0150       _RL    adVar
b4daa24319 Shre*0151       INTEGER Nr
                0152       INTEGER myThid
                0153       WRITE(*,*) "Called not yet defined"
642de41482 dngo*0154       RETURN
b4daa24319 Shre*0155       END
                0156 
6a9e386e2e dngo*0157 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
b4daa24319 Shre*0158       SUBROUTINE ADEXCH_UV_XY_RS( adU, adV, bool, myThid )
                0159       IMPLICIT NONE
642de41482 dngo*0160       _RS     adU, adV
b4daa24319 Shre*0161       LOGICAL bool
                0162       INTEGER myThid
                0163       WRITE(*,*) "Called not yet defined"
642de41482 dngo*0164       RETURN
b4daa24319 Shre*0165       END
                0166 
6a9e386e2e dngo*0167 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
642de41482 dngo*0168       SUBROUTINE ADEXCH_UV_3D_RL( aduVel, advVel, bool, Nr, myThid )
b4daa24319 Shre*0169       IMPLICIT NONE
642de41482 dngo*0170       _RL     aduVel, advVel
b4daa24319 Shre*0171       LOGICAL bool
                0172       INTEGER Nr
                0173       INTEGER myThid
                0174       WRITE(*,*) "Called not yet defined"
642de41482 dngo*0175       RETURN
b4daa24319 Shre*0176       END
                0177 
6a9e386e2e dngo*0178 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
b4daa24319 Shre*0179       SUBROUTINE ADEXCH_XY_RS( adVar, myThid )
                0180       IMPLICIT NONE
642de41482 dngo*0181       _RS     adVar
b4daa24319 Shre*0182       INTEGER myThid
                0183       WRITE(*,*) "Called not yet defined"
642de41482 dngo*0184       RETURN
b4daa24319 Shre*0185       END
                0186 
6a9e386e2e dngo*0187 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
b4daa24319 Shre*0188       SUBROUTINE ADEXCH_XY_RL( adVar, myThid )
                0189       IMPLICIT NONE
642de41482 dngo*0190       _RL     adVar
b4daa24319 Shre*0191       INTEGER myThid
                0192       WRITE(*,*) "Called not yet defined"
642de41482 dngo*0193       RETURN
b4daa24319 Shre*0194       END
                0195 
6a9e386e2e dngo*0196 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
642de41482 dngo*0197       SUBROUTINE STREAMICE_CG_SOLVE_D(
                0198      &                        cg_Uin, cg_Uind, cg_Vin, cg_Vind,
                0199      &                        cg_Bu, cg_Bud, cg_Bv, cg_Bvd,
                0200      &                        A_uu, A_uud, A_uv, A_uvd,
                0201      &                        A_vu, A_vud, A_vv, A_vvd,
                0202      &                        tolerance, iters, maxIter, myThid )
feb7fa5d1e dngo*0203       IMPLICIT NONE
                0204 
                0205 #include "SIZE.h"
                0206 #include "EEPARAMS.h"
                0207 #include "PARAMS.h"
                0208 
642de41482 dngo*0209       _RL cg_Uin(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0210       _RL cg_Vin(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
feb7fa5d1e dngo*0211       _RL cg_Bu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0212       _RL cg_Bv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
642de41482 dngo*0213       _RL A_uu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0214       _RL A_vu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0215       _RL A_uv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0216       _RL A_vv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0217       _RL cg_Uind(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0218       _RL cg_Vind(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0219       _RL cg_Bud (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0220       _RL cg_Bvd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0221       _RL A_uud(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0222       _RL A_vud(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0223       _RL A_uvd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0224       _RL A_vvd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0225       _RL tolerance
                0226       INTEGER iters
                0227       INTEGER maxIter
                0228       INTEGER myThid
                0229 
                0230 C     !LOCAL VARIABLES
                0231 #ifdef ALLOW_STREAMICE
                0232       INTEGER i, j, bi, bj
                0233       INTEGER colx, coly
                0234       _RL deltaAu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0235       _RL deltaAv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0236       _RL deltaAuMdBu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0237       _RL deltaAvMdBv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0238 
                0239       ! A x = b
                0240       ! delta A x + A delta x = delta b
                0241       ! A delta x = delta b - delta A x
                0242       ! delta x = A^{-1} (delta b - delta A x)
                0243 
                0244       CALL STREAMICE_CG_SOLVE( cg_Uin, cg_Vin,
                0245      &                         cg_Bu, cg_Bv,
                0246      &                         A_uu, A_uv, A_vu, A_vv,
                0247      &                         tolerance, iters, maxIter, myThid )
                0248 
                0249       _EXCH_XY_RL(cg_Uin, myThid)
                0250       _EXCH_XY_RL(cg_Vin, myThid)
                0251 
                0252       DO bj = myByLo(myThid), myByHi(myThid)
                0253        DO bi = myBxLo(myThid), myBxHi(myThid)
                0254         DO j=1,sNy
                0255          DO i=1,sNx
                0256 
                0257           deltaAu(i,j,bi,bj) = 0. _d 0
                0258           deltaAv(i,j,bi,bj) = 0. _d 0
                0259 
                0260           DO colx=-1,1
                0261            DO coly=-1,1
                0262             deltaAu(i,j,bi,bj) = deltaAu(i,j,bi,bj) +
                0263      &         A_uud(i,j,bi,bj,colx,coly)*
                0264      &         cg_Uin(i+colx,j+coly,bi,bj)+
                0265      &         A_uvd(i,j,bi,bj,colx,coly)*
                0266      &         cg_Vin(i+colx,j+coly,bi,bj)
                0267 
                0268             deltaAv(i,j,bi,bj) = deltaAv(i,j,bi,bj) +
                0269      &         A_vud(i,j,bi,bj,colx,coly)*
                0270      &         cg_Uin(i+colx,j+coly,bi,bj)+
                0271      &         A_vvd(i,j,bi,bj,colx,coly)*
                0272      &         cg_Vin(i+colx,j+coly,bi,bj)
                0273            ENDDO
                0274           ENDDO
                0275 
                0276          ENDDO
                0277         ENDDO
                0278        ENDDO
                0279       ENDDO
                0280 
                0281       _EXCH_XY_RL(deltaAu, myThid)
                0282       _EXCH_XY_RL(deltaAv, myThid)
                0283 
                0284       DO bj = myByLo(myThid), myByHi(myThid)
                0285        DO bi = myBxLo(myThid), myBxHi(myThid)
                0286         DO j=1-OLy,sNy+OLy
                0287          DO i=1-OLx,sNx+OLx
                0288           deltaAuMdBu(i,j,bi,bj) = cg_Bud(i,j,bi,bj)-deltaAu(i,j,bi,bj)
                0289           deltaAvMdBv(i,j,bi,bj) = cg_Bvd(i,j,bi,bj)-deltaAv(i,j,bi,bj)
                0290          ENDDO
                0291         ENDDO
                0292        ENDDO
                0293       ENDDO
                0294 
                0295       CALL STREAMICE_CG_SOLVE( cg_Uind, cg_Vind,
                0296      &                         deltaAuMdBu, deltaAvMdBv,
                0297      &                         A_uu, A_uv, A_vu, A_vv,
                0298      &                         tolerance, iters, maxIter, myThid )
                0299 
                0300 #endif
                0301       RETURN
feb7fa5d1e dngo*0302       END
                0303 
6a9e386e2e dngo*0304 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
642de41482 dngo*0305       SUBROUTINE STREAMICE_CG_SOLVE_B(
                0306      &                        cg_Uin, cg_Uinb, cg_Vin, cg_Vinb,
                0307      &                        cg_Bu, cg_Bub, cg_Bv, cg_Bvb,
                0308      &                        A_uu, A_uub, A_uv, A_uvb,
                0309      &                        A_vu, A_vub, A_vv, A_vvb,
                0310      &                        tolerance, iters, maxIter, myThid )
feb7fa5d1e dngo*0311 
                0312       IMPLICIT NONE
                0313 
                0314 #include "SIZE.h"
                0315 #include "EEPARAMS.h"
                0316 #include "PARAMS.h"
                0317 #ifdef ALLOW_STREAMICE
                0318 #include "STREAMICE.h"
                0319 #include "STREAMICE_CG.h"
                0320 #endif
                0321 
                0322       INTEGER myThid
                0323       _RL cg_Uin (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0324       _RL cg_Vin (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0325       _RL cg_Bu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0326       _RL cg_Bv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0327       _RL
                0328      & A_uu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1),
                0329      & A_vu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1),
                0330      & A_uv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1),
                0331      & A_vv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0332       _RL cg_Uinb (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0333       _RL cg_Vinb (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0334       _RL cg_Bub (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0335       _RL cg_Bvb (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0336       _RL
                0337      & A_uub (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1),
                0338      & A_vub (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1),
                0339      & A_uvb (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1),
                0340      & A_vvb (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
642de41482 dngo*0341       _RL tolerance
                0342       INTEGER iters
                0343       INTEGER maxIter
feb7fa5d1e dngo*0344 
                0345 #ifdef ALLOW_STREAMICE
                0346       INTEGER recalc_vel
                0347       recalc_vel = 1
                0348       CALL ADSTREAMICE_CG_SOLVE(
642de41482 dngo*0349      U  cg_Uin,    ! velocities - need to be recalc ed
                0350      I  cg_Uinb,      ! adjoint of vel (input)
                0351      U  cg_Vin,    ! velocities - need to be recalc ed
                0352      I  cg_Vinb,      ! adjoint of vel (input)
                0353      I  cg_Bu,   ! to recalc velocities
                0354      U  cg_Bub,     ! adjoint of RHS (output)
                0355      I  cg_Bv,   ! to recalc velocities
                0356      U  cg_Bvb,     ! adjoint of RHS (output)
                0357      I  A_uu,       ! section of matrix that multiplies u and projects on u
                0358      U  A_uub,     ! adjoint of matrix coeffs (output)
                0359      I  A_uv,       ! section of matrix that multiplies v and projects on u
                0360      U  A_uvb,     ! adjoint of matrix coeffs (output)
                0361      I  A_vu,       ! section of matrix that multiplies u and projects on v
                0362      U  A_vub,     ! adjoint of matrix coeffs (output)
                0363      I  A_vv,       ! section of matrix that multiplies v and projects on v
                0364      U  A_vvb,     ! adjoint of matrix coeffs (output)
                0365      I  tolerance,
                0366      I  maxIter,
                0367      I  recalc_vel,
                0368      I  myThid )
feb7fa5d1e dngo*0369 #endif
642de41482 dngo*0370       RETURN
feb7fa5d1e dngo*0371       END
                0372 
6a9e386e2e dngo*0373 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
feb7fa5d1e dngo*0374 C [llh] Hand-written forward-sweep adjoint code of STREAMICE_CG_SOLVE
                0375 C [llh] Adjoint of solving Ax=b,
                0376 C [llh]  where A is (Auu, Auv, Avu, Avv), x is (xu,xv), b is (bu,bv)
642de41482 dngo*0377       SUBROUTINE STREAMICE_CG_SOLVE_FWD(
                0378      &                        cg_Uin, cg_Vin, cg_Bu, cg_Bv,
                0379      &                        A_uu, A_uv, A_vu, A_vv,
                0380      &                        tolerance, iters, maxiters, myThid )
feb7fa5d1e dngo*0381 
                0382       IMPLICIT NONE
                0383 
                0384 #include "SIZE.h"
                0385 #include "EEPARAMS.h"
                0386 #include "PARAMS.h"
                0387 #ifdef ALLOW_STREAMICE
                0388 #include "STREAMICE.h"
                0389 #include "STREAMICE_CG.h"
                0390 #endif
                0391 
                0392       !INPUT/OUTPUT ARGUMENTS
                0393       _RL cg_Uin (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0394       _RL cg_Vin (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0395       _RL cg_Bu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0396       _RL cg_Bv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
642de41482 dngo*0397       _RL A_uu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0398       _RL A_vu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0399       _RL A_uv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0400       _RL A_vv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
feb7fa5d1e dngo*0401       _RL tolerance
                0402       INTEGER maxiters
                0403       INTEGER iters
                0404       INTEGER myThid
                0405 
                0406 #ifdef ALLOW_STREAMICE
                0407 C     !LOCAL VARIABLES
                0408       INTEGER i, j, bi, bj, colx, coly
                0409       INTEGER conv_flag, tmpiter
                0410       CHARACTER*(MAX_LEN_MBUF) msgBuf
                0411 
                0412       conv_flag = 0
                0413       CALL STREAMICE_CG_SOLVE(
642de41482 dngo*0414      &  cg_Uin,
                0415      &  cg_Vin,
                0416      &  cg_Bu,
                0417      &  cg_Bv,
feb7fa5d1e dngo*0418      &  A_uu,
                0419      &  A_uv,
                0420      &  A_vu,
                0421      &  A_vv,
                0422      &  tolerance,
                0423      &  iters,
                0424      &  maxiters,
                0425      &  myThid )
                0426 
642de41482 dngo*0427       CALL PUSHREAL8ARRAY(cg_Vin, (sNx+2*OLx)*(sNy+2*OLy)*nSx*nSy)
                0428       CALL PUSHREAL8ARRAY(cg_Uin, (sNx+2*OLx)*(sNy+2*OLy)*nSx*nSy)
feb7fa5d1e dngo*0429 #endif
642de41482 dngo*0430       RETURN
feb7fa5d1e dngo*0431       END
                0432 
6a9e386e2e dngo*0433 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
feb7fa5d1e dngo*0434 C [llh] Hand-written backward-sweep adjoint code of STREAMICE_CG_SOLVE
                0435 C [llh] Adjoint of solving Ax=b,
                0436 C [llh]  where A is (Auu, Auv, Avu, Avv), x is (xu,xv), b is (bu,bv)
                0437       SUBROUTINE STREAMICE_CG_SOLVE_BWD(
642de41482 dngo*0438      &                        cg_Uin, cg_Uinb, cg_Vin, cg_Vinb,
                0439      &                        cg_Bu, cg_Bub, cg_Bv, cg_Bvb,
                0440      &                        A_uu, A_uub, A_uv, A_uvb,
                0441      &                        A_vu, A_vub, A_vv, A_vvb,
decd05a68a dngo*0442      &                        tolerance, iters, maxIter, myThid )
feb7fa5d1e dngo*0443 
                0444       IMPLICIT NONE
                0445 
                0446 #include "SIZE.h"
                0447 #include "EEPARAMS.h"
                0448 #include "PARAMS.h"
                0449 #ifdef ALLOW_STREAMICE
                0450 #include "STREAMICE.h"
                0451 #include "STREAMICE_CG.h"
                0452 #endif
                0453 
                0454       _RL cg_Uin (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0455       _RL cg_Vin (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0456       _RL cg_Bu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0457       _RL cg_Bv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
642de41482 dngo*0458       _RL A_uu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0459       _RL A_vu (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0460       _RL A_uv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0461       _RL A_vv (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
feb7fa5d1e dngo*0462       _RL cg_Uinb (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0463       _RL cg_Vinb (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0464       _RL cg_Bub (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
                0465       _RL cg_Bvb (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
642de41482 dngo*0466       _RL A_uub (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0467       _RL A_vub (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0468       _RL A_uvb (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0469       _RL A_vvb (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy,-1:1,-1:1)
                0470       _RL tolerance
                0471 #ifdef ALLOW_PETSC
                0472       _RL toleranceb
                0473 #endif
                0474       INTEGER iters
                0475       INTEGER maxIter
                0476       INTEGER myThid
feb7fa5d1e dngo*0477 
                0478 #ifdef ALLOW_STREAMICE
6a9e386e2e dngo*0479       CALL POPREAL8ARRAY(cg_Uin, (sNx+2*OLx)*(sNy+2*OLy)*nSx*nSy)
                0480       CALL POPREAL8ARRAY(cg_Vin, (sNx+2*OLx)*(sNy+2*OLy)*nSx*nSy)
feb7fa5d1e dngo*0481 
                0482       CALL ADSTREAMICE_CG_SOLVE(
642de41482 dngo*0483      U  cg_Uin,    ! velocities - need to be recalc ed
                0484      I  cg_Uinb,      ! adjoint of vel (input)
                0485      U  cg_Vin,    ! velocities - need to be recalc ed
                0486      I  cg_Vinb,      ! adjoint of vel (input)
                0487      I  cg_Bu,   ! to recalc velocities
                0488      U  cg_Bub,     ! adjoint of RHS (output)
                0489      I  cg_Bv,   ! to recalc velocities
                0490      U  cg_Bvb,     ! adjoint of RHS (output)
                0491      I  A_uu,       ! section of matrix that multiplies u and projects on u
                0492      U  A_uub,     ! adjoint of matrix coeffs (output)
                0493      I  A_uv,       ! section of matrix that multiplies v and projects on u
                0494      U  A_uvb,     ! adjoint of matrix coeffs (output)
                0495      I  A_vu,       ! section of matrix that multiplies u and projects on v
                0496      U  A_vub,     ! adjoint of matrix coeffs (output)
                0497      I  A_vv,       ! section of matrix that multiplies v and projects on v
                0498      U  A_vvb,     ! adjoint of matrix coeffs (output)
                0499      I  tolerance,
                0500      I  maxIter,
                0501      I  myThid )
feb7fa5d1e dngo*0502 #endif
642de41482 dngo*0503       RETURN
feb7fa5d1e dngo*0504       END