/[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.5 by cnh, Thu Aug 20 20:21:23 1998 UTC revision 1.6 by cnh, Sat Aug 22 17:51:08 1998 UTC
# Line 20  C     == Global data == Line 20  C     == Global data ==
20    
21  C     == Routine Arguments ==  C     == Routine Arguments ==
22        INTEGER bi,bj,iMin,iMax,jMin,jMax        INTEGER bi,bj,iMin,iMax,jMin,jMax
23        _RL KappaRT(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)        _RL KappaRT(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
24        _RL KappaRS(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)        _RL KappaRS(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
25        INTEGER myThid        INTEGER myThid
26    
27  C     == Local variables ==  C     == Local variables ==
# Line 31  C     == Local variables == Line 31  C     == Local variables ==
31        _RL c(1-Olx:sNx+Olx,1-Oly:sNy+Oly)        _RL c(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
32        _RL ckm1(1-Olx:sNx+Olx,1-Oly:sNy+Oly)        _RL ckm1(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
33        _RL bet(1-Olx:sNx+Olx,1-Oly:sNy+Oly)        _RL bet(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
34        _RL gam(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)        _RL gam(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
35    
36  C ***************************************************************  C ***************************************************************
37  C *****************                              ****************  C *****************                              ****************
# Line 46  C Line 46  C
46  C ***************************************************************  C ***************************************************************
47  C ***************************************************************  C ***************************************************************
48    
49        IF (Nz.GT.1) THEN         ! Only need do anything if Nz>1        IF (Nr.GT.1) THEN         ! Only need do anything if Nr>1
50  C--    Beginning of forward sweep (top level)  C--    Beginning of forward sweep (top level)
51         DO j=jMin,jMax         DO j=jMin,jMax
52          DO i=iMin,iMax          DO i=iMin,iMax
53           c(i,j)=-deltaTtracer*_rhFacC(i,j,1,bi,bj)*recip_drF(1)           c(i,j)=-deltaTtracer*recip_hFacC(i,j,1,bi,bj)*recip_drF(1)
54       &         *KappaRT(i,j,2)*recip_drC(2)       &         *KappaRT(i,j,2)*recip_drC(2)
55           b(i,j)=1.-c(i,j)           b(i,j)=1.-c(i,j)
56           bet(i,j)=0.           bet(i,j)=0.
# Line 61  C--    Beginning of forward sweep (top l Line 61  C--    Beginning of forward sweep (top l
61         ENDDO         ENDDO
62        ENDIF        ENDIF
63  C--   Middle of forward sweep  C--   Middle of forward sweep
64        IF (Nz.GT.2) THEN        IF (Nr.GT.2) THEN
65         DO k=2,Nz-1         DO k=2,Nr-1
66          DO j=jMin,jMax          DO j=jMin,jMax
67           DO i=iMin,iMax           DO i=iMin,iMax
68            ckm1(i,j)=c(i,j)            ckm1(i,j)=c(i,j)
69            a(i,j)=-deltaTtracer*_rhFacC(i,j,k,bi,bj)*recip_drF(k)            a(i,j)=-deltaTtracer*recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
70       &         *KappaRT(i,j, k )*recip_drC( k )       &         *KappaRT(i,j, k )*recip_drC( k )
71            c(i,j)=-deltaTtracer*_rhFacC(i,j,k,bi,bj)*recip_drF(k)            c(i,j)=-deltaTtracer*recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
72       &         *KappaRT(i,j,k+1)*recip_drC(k+1)       &         *KappaRT(i,j,k+1)*recip_drC(k+1)
73            b(i,j)=1.-c(i,j)-a(i,j)            b(i,j)=1.-c(i,j)-a(i,j)
74           ENDDO           ENDDO
# Line 86  C--   Middle of forward sweep Line 86  C--   Middle of forward sweep
86          ENDDO          ENDDO
87         ENDDO         ENDDO
88        ENDIF        ENDIF
89        IF (Nz.GT.1) THEN        IF (Nr.GT.1) THEN
90  C--    End of forward sweep (bottom level)  C--    End of forward sweep (bottom level)
91         DO j=jMin,jMax         DO j=jMin,jMax
92          DO i=iMin,iMax          DO i=iMin,iMax
93           ckm1(i,j)=c(i,j)           ckm1(i,j)=c(i,j)
94           a(i,j)=-deltaTtracer*_rhFacC(i,j,Nz,bi,bj)*recip_drF(Nz)           a(i,j)=-deltaTtracer*recip_hFacC(i,j,Nr,bi,bj)*recip_drF(Nr)
95       &        *KappaRT(i,j, Nz )*recip_drC( Nz )       &        *KappaRT(i,j, Nr )*recip_drC( Nr )
96           b(i,j)=1.-a(i,j)           b(i,j)=1.-a(i,j)
97          ENDDO          ENDDO
98         ENDDO         ENDDO
99         DO j=jMin,jMax         DO j=jMin,jMax
100          DO i=iMin,iMax          DO i=iMin,iMax
101           gam(i,j,Nz)=ckm1(i,j)*bet(i,j)           gam(i,j,Nr)=ckm1(i,j)*bet(i,j)
102           bet(i,j)=b(i,j)-a(i,j)*gam(i,j,Nz)           bet(i,j)=b(i,j)-a(i,j)*gam(i,j,Nr)
103           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)
104           gTnm1(i,j,Nz,bi,bj)=(gTnm1(i,j,Nz,bi,bj)           gTnm1(i,j,Nr,bi,bj)=(gTnm1(i,j,Nr,bi,bj)
105       &      -a(i,j)*gTnm1(i,j,Nz-1,bi,bj))*bet(i,j)       &      -a(i,j)*gTnm1(i,j,Nr-1,bi,bj))*bet(i,j)
106           gSnm1(i,j,Nz,bi,bj)=(gSnm1(i,j,Nz,bi,bj)           gSnm1(i,j,Nr,bi,bj)=(gSnm1(i,j,Nr,bi,bj)
107       &      -a(i,j)*gSnm1(i,j,Nz-1,bi,bj))*bet(i,j)       &      -a(i,j)*gSnm1(i,j,Nr-1,bi,bj))*bet(i,j)
108          ENDDO          ENDDO
109         ENDDO         ENDDO
110  C--    Backward sweep  C--    Backward sweep
111         DO k=Nz-1,1,-1         DO k=Nr-1,1,-1
112          DO j=jMin,jMax          DO j=jMin,jMax
113           DO i=iMin,iMax           DO i=iMin,iMax
114            gTnm1(i,j,k,bi,bj)=gTnm1(i,j,k,bi,bj)            gTnm1(i,j,k,bi,bj)=gTnm1(i,j,k,bi,bj)

Legend:
Removed from v.1.5  
changed lines
  Added in v.1.6

  ViewVC Help
Powered by ViewVC 1.1.22