/[MITgcm]/MITgcm/model/src/impldiff.F
ViewVC logotype

Diff of /MITgcm/model/src/impldiff.F

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

revision 1.7 by cnh, Fri Nov 6 22:44:46 1998 UTC revision 1.8 by adcroft, Tue May 18 18:01:13 1999 UTC
# Line 8  C     | o Solve implicit diffusion equat Line 8  C     | o Solve implicit diffusion equat
8  C     |   diffusivity.                                           |  C     |   diffusivity.                                           |
9  C     \==========================================================/  C     \==========================================================/
10        SUBROUTINE IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,        SUBROUTINE IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,
11       I                     KappaRT,KappaRS,       I                     deltaTX,KappaRX,recip_hFac,
12         U                     gXNm1,
13       I                     myThid )       I                     myThid )
14        IMPLICIT NONE        IMPLICIT NONE
15  C     == Global data ==  C     == Global data ==
# Line 20  C     == Global data == Line 21  C     == Global data ==
21    
22  C     == Routine Arguments ==  C     == Routine Arguments ==
23        INTEGER bi,bj,iMin,iMax,jMin,jMax        INTEGER bi,bj,iMin,iMax,jMin,jMax
24        _RL KappaRT(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)        _RL deltaTX
25        _RL KappaRS(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)        _RL KappaRX(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
26          _RS recip_hFac(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr,nSx,nSy)
27          _RL gXnm1(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr,nSx,nSy)
28        INTEGER myThid        INTEGER myThid
29    
30  C     == Local variables ==  C     == Local variables ==
# Line 33  C     == Local variables == Line 36  C     == Local variables ==
36        _RL bet(1-Olx:sNx+Olx,1-Oly:sNy+Oly)        _RL bet(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
37        _RL gam(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)        _RL gam(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
38    
 C ***************************************************************  
 C *****************                              ****************  
 C *****************           N O T E            ****************  
 C *****************                              ****************  
 C ***************************************************************  
 C  
 C The implicit diffusion of SALT currently uses the diffusivity  
 C of the THETA.  
 C                  ie. KappaRS is ignored.  
 C  
 C ***************************************************************  
 C ***************************************************************  
   
39        IF (Nr.GT.1) THEN         ! Only need do anything if Nr>1        IF (Nr.GT.1) THEN         ! Only need do anything if Nr>1
40  C--    Beginning of forward sweep (top level)  C--    Beginning of forward sweep (top level)
41         DO j=jMin,jMax         DO j=jMin,jMax
42          DO i=iMin,iMax          DO i=iMin,iMax
43           c(i,j)=-deltaTtracer*recip_hFacC(i,j,1,bi,bj)*recip_drF(1)           c(i,j)=-deltaTX*recip_hFac(i,j,1,bi,bj)*recip_drF(1)
44       &         *KappaRT(i,j,2)*recip_drC(2)       &         *KappaRX(i,j,2)*recip_drC(2)
45           b(i,j)=1.-c(i,j)           b(i,j)=1.-c(i,j)
46           bet(i,j)=0.           bet(i,j)=0.
47           IF (b(i,j).NE.0.) bet(i,j)=1. / b(i,j)           IF (b(i,j).NE.0.) bet(i,j)=1. / b(i,j)
48           gTnm1(i,j,1,bi,bj)=gTnm1(i,j,1,bi,bj)*bet(i,j)           gXnm1(i,j,1,bi,bj)=gXnm1(i,j,1,bi,bj)*bet(i,j)
          gSnm1(i,j,1,bi,bj)=gSnm1(i,j,1,bi,bj)*bet(i,j)  
49          ENDDO          ENDDO
50         ENDDO         ENDDO
51        ENDIF        ENDIF
# Line 66  C--   Middle of forward sweep Line 55  C--   Middle of forward sweep
55          DO j=jMin,jMax          DO j=jMin,jMax
56           DO i=iMin,iMax           DO i=iMin,iMax
57            ckm1(i,j)=c(i,j)            ckm1(i,j)=c(i,j)
58            a(i,j)=-deltaTtracer*recip_hFacC(i,j,k,bi,bj)*recip_drF(k)            a(i,j)=-deltaTX*recip_hFac(i,j,k,bi,bj)*recip_drF(k)
59       &         *KappaRT(i,j, k )*recip_drC( k )       &         *KappaRX(i,j, k )*recip_drC( k )
60            c(i,j)=-deltaTtracer*recip_hFacC(i,j,k,bi,bj)*recip_drF(k)            c(i,j)=-deltaTX*recip_hFac(i,j,k,bi,bj)*recip_drF(k)
61       &         *KappaRT(i,j,k+1)*recip_drC(k+1)       &         *KappaRX(i,j,k+1)*recip_drC(k+1)
62            b(i,j)=1.-c(i,j)-a(i,j)            b(i,j)=1.-c(i,j)-a(i,j)
63           ENDDO           ENDDO
64          ENDDO          ENDDO
# Line 78  C--   Middle of forward sweep Line 67  C--   Middle of forward sweep
67            gam(i,j,k)=ckm1(i,j)*bet(i,j)            gam(i,j,k)=ckm1(i,j)*bet(i,j)
68            bet(i,j)=b(i,j)-a(i,j)*gam(i,j,k)            bet(i,j)=b(i,j)-a(i,j)*gam(i,j,k)
69            IF (bet(i,j).NE.0.) bet(i,j)=1. / bet(i,j)            IF (bet(i,j).NE.0.) bet(i,j)=1. / bet(i,j)
70            gTnm1(i,j,k,bi,bj)=(gTnm1(i,j,k,bi,bj)            gXnm1(i,j,k,bi,bj)=(gXnm1(i,j,k,bi,bj)
71       &       -a(i,j)*gTnm1(i,j,k-1,bi,bj))*bet(i,j)       &       -a(i,j)*gXnm1(i,j,k-1,bi,bj))*bet(i,j)
           gSnm1(i,j,k,bi,bj)=(gSnm1(i,j,k,bi,bj)  
      &       -a(i,j)*gSnm1(i,j,k-1,bi,bj))*bet(i,j)  
72           ENDDO           ENDDO
73          ENDDO          ENDDO
74         ENDDO         ENDDO
# Line 91  C--    End of forward sweep (bottom leve Line 78  C--    End of forward sweep (bottom leve
78         DO j=jMin,jMax         DO j=jMin,jMax
79          DO i=iMin,iMax          DO i=iMin,iMax
80           ckm1(i,j)=c(i,j)           ckm1(i,j)=c(i,j)
81           a(i,j)=-deltaTtracer*recip_hFacC(i,j,Nr,bi,bj)*recip_drF(Nr)           a(i,j)=-deltaTX*recip_hFac(i,j,Nr,bi,bj)*recip_drF(Nr)
82       &        *KappaRT(i,j, Nr )*recip_drC( Nr )       &        *KappaRX(i,j, Nr )*recip_drC( Nr )
83           b(i,j)=1.-a(i,j)           b(i,j)=1.-a(i,j)
84          ENDDO          ENDDO
85         ENDDO         ENDDO
# Line 101  C--    End of forward sweep (bottom leve Line 88  C--    End of forward sweep (bottom leve
88           gam(i,j,Nr)=ckm1(i,j)*bet(i,j)           gam(i,j,Nr)=ckm1(i,j)*bet(i,j)
89           bet(i,j)=b(i,j)-a(i,j)*gam(i,j,Nr)           bet(i,j)=b(i,j)-a(i,j)*gam(i,j,Nr)
90           IF (bet(i,j).NE.0.) bet(i,j)=1. / bet(i,j)           IF (bet(i,j).NE.0.) bet(i,j)=1. / bet(i,j)
91           gTnm1(i,j,Nr,bi,bj)=(gTnm1(i,j,Nr,bi,bj)           gXnm1(i,j,Nr,bi,bj)=(gXnm1(i,j,Nr,bi,bj)
92       &      -a(i,j)*gTnm1(i,j,Nr-1,bi,bj))*bet(i,j)       &      -a(i,j)*gXnm1(i,j,Nr-1,bi,bj))*bet(i,j)
          gSnm1(i,j,Nr,bi,bj)=(gSnm1(i,j,Nr,bi,bj)  
      &      -a(i,j)*gSnm1(i,j,Nr-1,bi,bj))*bet(i,j)  
93          ENDDO          ENDDO
94         ENDDO         ENDDO
95  C--    Backward sweep  C--    Backward sweep
96         DO k=Nr-1,1,-1         DO k=Nr-1,1,-1
97          DO j=jMin,jMax          DO j=jMin,jMax
98           DO i=iMin,iMax           DO i=iMin,iMax
99            gTnm1(i,j,k,bi,bj)=gTnm1(i,j,k,bi,bj)            gXnm1(i,j,k,bi,bj)=gXnm1(i,j,k,bi,bj)
100       &              -gam(i,j,k+1)*gTnm1(i,j,k+1,bi,bj)       &              -gam(i,j,k+1)*gXnm1(i,j,k+1,bi,bj)
           gSnm1(i,j,k,bi,bj)=gSnm1(i,j,k,bi,bj)  
      &              -gam(i,j,k+1)*gSnm1(i,j,k+1,bi,bj)  
101           ENDDO           ENDDO
102          ENDDO          ENDDO
103         ENDDO         ENDDO

Legend:
Removed from v.1.7  
changed lines
  Added in v.1.8

  ViewVC Help
Powered by ViewVC 1.1.22