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

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

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


Revision 1.2 - (show annotations) (download)
Tue Jun 9 16:48:02 1998 UTC (25 years, 11 months ago) by cnh
Branch: MAIN
Changes since 1.1: +2 -1 lines
Changes to support topography, hydrography and
forcing from files

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/impldiff.F,v 1.1 1998/06/01 22:27:14 adcroft Exp $
2
3 #include "CPP_EEOPTIONS.h"
4
5 C /==========================================================\
6 C | S/R IMPLDIFF |
7 C | o Step model fields forward in time |
8 C \==========================================================/
9 SUBROUTINE IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax,
10 I KappaZT,
11 I myThid )
12 implicit none
13 ! Common
14 #include "SIZE.h"
15 #include "DYNVARS.h"
16 #include "EEPARAMS.h"
17 #include "PARAMS.h"
18 #include "GRID.h"
19 C == Routine Arguments ==
20 INTEGER bi,bj,iMin,iMax,jMin,jMax
21 _RL KappaZT(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)
22 INTEGER myThid
23 C == Local variables ==
24 INTEGER i,j,k
25 _RL a(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
26 _RL b(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
27 _RL c(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
28 _RL ckm1(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
29 _RL bet(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
30 _RL gam(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)
31
32 IF (Nz.GT.1) THEN ! Only need do anything if Nz>1
33 C Beginning of forward sweep (top level)
34 DO j=jMin,jMax
35 DO i=iMin,iMax
36 c(i,j)=-deltaTtracer*_rhFacC(i,j,1,bi,bj)*rdzF(1)
37 & *KappaZT(i,j,2)*rdzC(2)
38 b(i,j)=1.-c(i,j)
39 bet(i,j)=0.
40 IF (b(i,j).NE.0.) bet(i,j)=1. / b(i,j)
41 theta(i,j,1,bi,bj)=theta(i,j,1,bi,bj)*bet(i,j)
42 ENDDO
43 ENDDO
44 ENDIF
45 C Middle of forward sweep
46 IF (Nz.GT.2) THEN
47 DO k=2,Nz-1
48 DO j=jMin,jMax
49 DO i=iMin,iMax
50 ckm1(i,j)=c(i,j)
51 a(i,j)=-deltaTtracer*_rhFacC(i,j,k,bi,bj)*rdzF(k)
52 & *KappaZT(i,j, k )*rdzC( k )
53 c(i,j)=-deltaTtracer*_rhFacC(i,j,k,bi,bj)*rdzF(k)
54 & *KappaZT(i,j,k+1)*rdzC(k+1)
55 b(i,j)=1.-c(i,j)-a(i,j)
56 ENDDO
57 ENDDO
58 DO j=jMin,jMax
59 DO i=iMin,iMax
60 gam(i,j,k)=ckm1(i,j)*bet(i,j)
61 bet(i,j)=b(i,j)-a(i,j)*gam(i,j,k)
62 IF (bet(i,j).NE.0.) bet(i,j)=1. / bet(i,j)
63 theta(i,j,k,bi,bj)=(theta(i,j,k,bi,bj)
64 & -a(i,j)*theta(i,j,k-1,bi,bj))*bet(i,j)
65 ENDDO
66 ENDDO
67 ENDDO
68 ENDIF
69 IF (Nz.GT.1) THEN
70 C End of forward sweep (bottom level)
71 DO j=jMin,jMax
72 DO i=iMin,iMax
73 ckm1(i,j)=c(i,j)
74 a(i,j)=-deltaTtracer*_rhFacC(i,j,Nz,bi,bj)*rdzF(Nz)
75 & *KappaZT(i,j, Nz )*rdzC( Nz )
76 b(i,j)=1.-a(i,j)
77 ENDDO
78 ENDDO
79 DO j=jMin,jMax
80 DO i=iMin,iMax
81 gam(i,j,Nz)=ckm1(i,j)*bet(i,j)
82 bet(i,j)=b(i,j)-a(i,j)*gam(i,j,Nz)
83 IF (bet(i,j).NE.0.) bet(i,j)=1. / bet(i,j)
84 theta(i,j,Nz,bi,bj)=(theta(i,j,Nz,bi,bj)
85 & -a(i,j)*theta(i,j,Nz-1,bi,bj))*bet(i,j)
86 ENDDO
87 ENDDO
88 C Backward sweep
89 DO k=Nz-1,1,-1
90 DO j=jMin,jMax
91 DO i=iMin,iMax
92 theta(i,j,k,bi,bj)=theta(i,j,k,bi,bj)
93 & -gam(i,j,k+1)*theta(i,j,k+1,bi,bj)
94 ENDDO
95 ENDDO
96 ENDDO
97 ENDIF
98
99 RETURN
100 END

  ViewVC Help
Powered by ViewVC 1.1.22