/[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.1 - (show annotations) (download)
Mon Jun 1 22:27:14 1998 UTC (26 years ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint5, checkpoint6
Implemented implicit vertical diffusion (tracers only).
Involved introducing a "total" diffusivity array (local 3D)
calculated by calc_diffusivity().
Made some small changes to time-stepping algorithm.
Switched on by setting implicitZdiffusion.
(note: *Not* fully tested with topography. But when switched off
this does produce identical results)

1 C $Header$
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 "PARAMS.h"
17 #include "GRID.h"
18 C == Routine Arguments ==
19 INTEGER bi,bj,iMin,iMax,jMin,jMax
20 _RL KappaZT(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)
21 INTEGER myThid
22 C == Local variables ==
23 INTEGER i,j,k
24 _RL a(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
25 _RL b(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
26 _RL c(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
27 _RL ckm1(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
28 _RL bet(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
29 _RL gam(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz)
30
31 IF (Nz.GT.1) THEN ! Only need do anything if Nz>1
32 C Beginning of forward sweep (top level)
33 DO j=jMin,jMax
34 DO i=iMin,iMax
35 c(i,j)=-deltaTtracer*_rhFacC(i,j,1,bi,bj)*rdzF(1)
36 & *KappaZT(i,j,2)*rdzC(2)
37 b(i,j)=1.-c(i,j)
38 bet(i,j)=0.
39 IF (b(i,j).NE.0.) bet(i,j)=1. / b(i,j)
40 theta(i,j,1,bi,bj)=theta(i,j,1,bi,bj)*bet(i,j)
41 ENDDO
42 ENDDO
43 ENDIF
44 C Middle of forward sweep
45 IF (Nz.GT.2) THEN
46 DO k=2,Nz-1
47 DO j=jMin,jMax
48 DO i=iMin,iMax
49 ckm1(i,j)=c(i,j)
50 a(i,j)=-deltaTtracer*_rhFacC(i,j,k,bi,bj)*rdzF(k)
51 & *KappaZT(i,j, k )*rdzC( k )
52 c(i,j)=-deltaTtracer*_rhFacC(i,j,k,bi,bj)*rdzF(k)
53 & *KappaZT(i,j,k+1)*rdzC(k+1)
54 b(i,j)=1.-c(i,j)-a(i,j)
55 ENDDO
56 ENDDO
57 DO j=jMin,jMax
58 DO i=iMin,iMax
59 gam(i,j,k)=ckm1(i,j)*bet(i,j)
60 bet(i,j)=b(i,j)-a(i,j)*gam(i,j,k)
61 IF (bet(i,j).NE.0.) bet(i,j)=1. / bet(i,j)
62 theta(i,j,k,bi,bj)=(theta(i,j,k,bi,bj)
63 & -a(i,j)*theta(i,j,k-1,bi,bj))*bet(i,j)
64 ENDDO
65 ENDDO
66 ENDDO
67 ENDIF
68 IF (Nz.GT.1) THEN
69 C End of forward sweep (bottom level)
70 DO j=jMin,jMax
71 DO i=iMin,iMax
72 ckm1(i,j)=c(i,j)
73 a(i,j)=-deltaTtracer*_rhFacC(i,j,Nz,bi,bj)*rdzF(Nz)
74 & *KappaZT(i,j, Nz )*rdzC( Nz )
75 b(i,j)=1.-a(i,j)
76 ENDDO
77 ENDDO
78 DO j=jMin,jMax
79 DO i=iMin,iMax
80 gam(i,j,Nz)=ckm1(i,j)*bet(i,j)
81 bet(i,j)=b(i,j)-a(i,j)*gam(i,j,Nz)
82 IF (bet(i,j).NE.0.) bet(i,j)=1. / bet(i,j)
83 theta(i,j,Nz,bi,bj)=(theta(i,j,Nz,bi,bj)
84 & -a(i,j)*theta(i,j,Nz-1,bi,bj))*bet(i,j)
85 ENDDO
86 ENDDO
87 C Backward sweep
88 DO k=Nz-1,1,-1
89 DO j=jMin,jMax
90 DO i=iMin,iMax
91 theta(i,j,k,bi,bj)=theta(i,j,k,bi,bj)
92 & -gam(i,j,k+1)*theta(i,j,k+1,bi,bj)
93 ENDDO
94 ENDDO
95 ENDDO
96 ENDIF
97
98 RETURN
99 END

  ViewVC Help
Powered by ViewVC 1.1.22