/[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.4 by adcroft, Tue Jun 16 15:20:08 1998 UTC revision 1.5 by cnh, Thu Aug 20 20:21:23 1998 UTC
# Line 4  C $Header$ Line 4  C $Header$
4    
5  C     /==========================================================\  C     /==========================================================\
6  C     | S/R IMPLDIFF                                             |  C     | S/R IMPLDIFF                                             |
7  C     | o Step model fields forward in time                      |  C     | o Solve implicit diffusion equation for vertical         |
8    C     |   diffusivity.                                           |
9  C     \==========================================================/  C     \==========================================================/
10        SUBROUTINE IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,        SUBROUTINE IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,
11       I                     KappaZT,KappaZS,       I                     KappaRT,KappaRS,
12       I                     myThid )       I                     myThid )
13        implicit none        IMPLICIT NONE
14  ! Common  C     == Global data ==
15  #include "SIZE.h"  #include "SIZE.h"
16  #include "DYNVARS.h"  #include "DYNVARS.h"
17  #include "EEPARAMS.h"  #include "EEPARAMS.h"
18  #include "PARAMS.h"  #include "PARAMS.h"
19  #include "GRID.h"  #include "GRID.h"
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 KappaZT(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)        _RL KappaRT(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)
24        _RL KappaZS(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)        _RL KappaRS(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)
25        INTEGER myThid        INTEGER myThid
26    
27  C     == Local variables ==  C     == Local variables ==
28        INTEGER i,j,k        INTEGER i,j,k
29        _RL a(1-Olx:sNx+Olx,1-Oly:sNy+Oly)        _RL a(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
# Line 38  C ************************************** Line 41  C **************************************
41  C  C
42  C The implicit diffusion of SALT currently uses the diffusivity  C The implicit diffusion of SALT currently uses the diffusivity
43  C of the THETA.  C of the THETA.
44  C                  ie. KappaZS is ignored.  C                  ie. KappaRS is ignored.
 C  
 C This is only a temporary situation. We appreciate your  
 C patience and understanding.  
 C  
 C                              Signed,  The Management.  
45  C  C
46  C ***************************************************************  C ***************************************************************
47  C ***************************************************************  C ***************************************************************
48    
49        IF (Nz.GT.1) THEN         ! Only need do anything if Nz>1        IF (Nz.GT.1) THEN         ! Only need do anything if Nz>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)*rdzF(1)           c(i,j)=-deltaTtracer*_rhFacC(i,j,1,bi,bj)*recip_drF(1)
54       &         *KappaZT(i,j,2)*rdzC(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.
57           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)
# Line 62  C Beginning of forward sweep (top level) Line 60  C Beginning of forward sweep (top level)
60          ENDDO          ENDDO
61         ENDDO         ENDDO
62        ENDIF        ENDIF
63  C Middle of forward sweep  C--   Middle of forward sweep
64        IF (Nz.GT.2) THEN        IF (Nz.GT.2) THEN
65         DO k=2,Nz-1         DO k=2,Nz-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)*rdzF(k)            a(i,j)=-deltaTtracer*_rhFacC(i,j,k,bi,bj)*recip_drF(k)
70       &         *KappaZT(i,j, k )*rdzC( k )       &         *KappaRT(i,j, k )*recip_drC( k )
71            c(i,j)=-deltaTtracer*_rhFacC(i,j,k,bi,bj)*rdzF(k)            c(i,j)=-deltaTtracer*_rhFacC(i,j,k,bi,bj)*recip_drF(k)
72       &         *KappaZT(i,j,k+1)*rdzC(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
75          ENDDO          ENDDO
# Line 89  C Middle of forward sweep Line 87  C Middle of forward sweep
87         ENDDO         ENDDO
88        ENDIF        ENDIF
89        IF (Nz.GT.1) THEN        IF (Nz.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)*rdzF(Nz)           a(i,j)=-deltaTtracer*_rhFacC(i,j,Nz,bi,bj)*recip_drF(Nz)
95       &        *KappaZT(i,j, Nz )*rdzC( Nz )       &        *KappaRT(i,j, Nz )*recip_drC( Nz )
96           b(i,j)=1.-a(i,j)           b(i,j)=1.-a(i,j)
97          ENDDO          ENDDO
98         ENDDO         ENDDO
# Line 109  C End of forward sweep (bottom level) Line 107  C End of forward sweep (bottom level)
107       &      -a(i,j)*gSnm1(i,j,Nz-1,bi,bj))*bet(i,j)       &      -a(i,j)*gSnm1(i,j,Nz-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=Nz-1,1,-1
112          DO j=jMin,jMax          DO j=jMin,jMax
113           DO i=iMin,iMax           DO i=iMin,iMax

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

  ViewVC Help
Powered by ViewVC 1.1.22