1 |
adcroft |
1.3 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/impldiff.F,v 1.2 1998/06/09 16:48:02 cnh Exp $ |
2 |
adcroft |
1.1 |
|
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 |
adcroft |
1.3 |
I KappaZT,KappaZS, |
11 |
adcroft |
1.1 |
I myThid ) |
12 |
|
|
implicit none |
13 |
|
|
! Common |
14 |
|
|
#include "SIZE.h" |
15 |
|
|
#include "DYNVARS.h" |
16 |
cnh |
1.2 |
#include "EEPARAMS.h" |
17 |
adcroft |
1.1 |
#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 |
adcroft |
1.3 |
_RL KappaZS(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz) |
23 |
adcroft |
1.1 |
INTEGER myThid |
24 |
|
|
C == Local variables == |
25 |
|
|
INTEGER i,j,k |
26 |
|
|
_RL a(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
27 |
|
|
_RL b(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
28 |
|
|
_RL c(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
29 |
|
|
_RL ckm1(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
30 |
|
|
_RL bet(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
31 |
|
|
_RL gam(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz) |
32 |
|
|
|
33 |
adcroft |
1.3 |
C *************************************************************** |
34 |
|
|
C ***************** **************** |
35 |
|
|
C ***************** N O T E **************** |
36 |
|
|
C ***************** **************** |
37 |
|
|
C *************************************************************** |
38 |
|
|
C |
39 |
|
|
C The implicit diffusion of SALT currently uses the diffusivity |
40 |
|
|
C of the THETA. |
41 |
|
|
C ie. KappaZS is ignored. |
42 |
|
|
C |
43 |
|
|
C This is only a temporary situation. We appreciate your |
44 |
|
|
C patience and understanding. |
45 |
|
|
C |
46 |
|
|
C Signed, The Management. |
47 |
|
|
C |
48 |
|
|
C *************************************************************** |
49 |
|
|
C *************************************************************** |
50 |
|
|
|
51 |
adcroft |
1.1 |
IF (Nz.GT.1) THEN ! Only need do anything if Nz>1 |
52 |
|
|
C Beginning of forward sweep (top level) |
53 |
|
|
DO j=jMin,jMax |
54 |
|
|
DO i=iMin,iMax |
55 |
|
|
c(i,j)=-deltaTtracer*_rhFacC(i,j,1,bi,bj)*rdzF(1) |
56 |
|
|
& *KappaZT(i,j,2)*rdzC(2) |
57 |
|
|
b(i,j)=1.-c(i,j) |
58 |
|
|
bet(i,j)=0. |
59 |
|
|
IF (b(i,j).NE.0.) bet(i,j)=1. / b(i,j) |
60 |
|
|
theta(i,j,1,bi,bj)=theta(i,j,1,bi,bj)*bet(i,j) |
61 |
adcroft |
1.3 |
salt(i,j,1,bi,bj)=salt(i,j,1,bi,bj)*bet(i,j) |
62 |
adcroft |
1.1 |
ENDDO |
63 |
|
|
ENDDO |
64 |
|
|
ENDIF |
65 |
|
|
C Middle of forward sweep |
66 |
|
|
IF (Nz.GT.2) THEN |
67 |
|
|
DO k=2,Nz-1 |
68 |
|
|
DO j=jMin,jMax |
69 |
|
|
DO i=iMin,iMax |
70 |
|
|
ckm1(i,j)=c(i,j) |
71 |
|
|
a(i,j)=-deltaTtracer*_rhFacC(i,j,k,bi,bj)*rdzF(k) |
72 |
|
|
& *KappaZT(i,j, k )*rdzC( k ) |
73 |
|
|
c(i,j)=-deltaTtracer*_rhFacC(i,j,k,bi,bj)*rdzF(k) |
74 |
|
|
& *KappaZT(i,j,k+1)*rdzC(k+1) |
75 |
|
|
b(i,j)=1.-c(i,j)-a(i,j) |
76 |
|
|
ENDDO |
77 |
|
|
ENDDO |
78 |
|
|
DO j=jMin,jMax |
79 |
|
|
DO i=iMin,iMax |
80 |
|
|
gam(i,j,k)=ckm1(i,j)*bet(i,j) |
81 |
|
|
bet(i,j)=b(i,j)-a(i,j)*gam(i,j,k) |
82 |
|
|
IF (bet(i,j).NE.0.) bet(i,j)=1. / bet(i,j) |
83 |
|
|
theta(i,j,k,bi,bj)=(theta(i,j,k,bi,bj) |
84 |
|
|
& -a(i,j)*theta(i,j,k-1,bi,bj))*bet(i,j) |
85 |
adcroft |
1.3 |
salt(i,j,k,bi,bj)=(salt(i,j,k,bi,bj) |
86 |
|
|
& -a(i,j)*salt(i,j,k-1,bi,bj))*bet(i,j) |
87 |
adcroft |
1.1 |
ENDDO |
88 |
|
|
ENDDO |
89 |
|
|
ENDDO |
90 |
|
|
ENDIF |
91 |
|
|
IF (Nz.GT.1) THEN |
92 |
|
|
C End of forward sweep (bottom level) |
93 |
|
|
DO j=jMin,jMax |
94 |
|
|
DO i=iMin,iMax |
95 |
|
|
ckm1(i,j)=c(i,j) |
96 |
|
|
a(i,j)=-deltaTtracer*_rhFacC(i,j,Nz,bi,bj)*rdzF(Nz) |
97 |
|
|
& *KappaZT(i,j, Nz )*rdzC( Nz ) |
98 |
|
|
b(i,j)=1.-a(i,j) |
99 |
|
|
ENDDO |
100 |
|
|
ENDDO |
101 |
|
|
DO j=jMin,jMax |
102 |
|
|
DO i=iMin,iMax |
103 |
|
|
gam(i,j,Nz)=ckm1(i,j)*bet(i,j) |
104 |
|
|
bet(i,j)=b(i,j)-a(i,j)*gam(i,j,Nz) |
105 |
|
|
IF (bet(i,j).NE.0.) bet(i,j)=1. / bet(i,j) |
106 |
|
|
theta(i,j,Nz,bi,bj)=(theta(i,j,Nz,bi,bj) |
107 |
|
|
& -a(i,j)*theta(i,j,Nz-1,bi,bj))*bet(i,j) |
108 |
adcroft |
1.3 |
salt(i,j,Nz,bi,bj)=(salt(i,j,Nz,bi,bj) |
109 |
|
|
& -a(i,j)*salt(i,j,Nz-1,bi,bj))*bet(i,j) |
110 |
adcroft |
1.1 |
ENDDO |
111 |
|
|
ENDDO |
112 |
|
|
C Backward sweep |
113 |
|
|
DO k=Nz-1,1,-1 |
114 |
|
|
DO j=jMin,jMax |
115 |
|
|
DO i=iMin,iMax |
116 |
|
|
theta(i,j,k,bi,bj)=theta(i,j,k,bi,bj) |
117 |
|
|
& -gam(i,j,k+1)*theta(i,j,k+1,bi,bj) |
118 |
adcroft |
1.3 |
salt(i,j,k,bi,bj)=salt(i,j,k,bi,bj) |
119 |
|
|
& -gam(i,j,k+1)*salt(i,j,k+1,bi,bj) |
120 |
adcroft |
1.1 |
ENDDO |
121 |
|
|
ENDDO |
122 |
|
|
ENDDO |
123 |
|
|
ENDIF |
124 |
|
|
|
125 |
|
|
RETURN |
126 |
|
|
END |