1 |
cnh |
1.7 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/impldiff.F,v 1.6 1998/08/22 17:51:08 cnh Exp $ |
2 |
adcroft |
1.1 |
|
3 |
cnh |
1.7 |
#include "CPP_OPTIONS.h" |
4 |
adcroft |
1.1 |
|
5 |
|
|
C /==========================================================\ |
6 |
|
|
C | S/R IMPLDIFF | |
7 |
cnh |
1.5 |
C | o Solve implicit diffusion equation for vertical | |
8 |
|
|
C | diffusivity. | |
9 |
adcroft |
1.1 |
C \==========================================================/ |
10 |
|
|
SUBROUTINE IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax, |
11 |
cnh |
1.5 |
I KappaRT,KappaRS, |
12 |
adcroft |
1.1 |
I myThid ) |
13 |
cnh |
1.5 |
IMPLICIT NONE |
14 |
|
|
C == Global data == |
15 |
adcroft |
1.1 |
#include "SIZE.h" |
16 |
|
|
#include "DYNVARS.h" |
17 |
cnh |
1.2 |
#include "EEPARAMS.h" |
18 |
adcroft |
1.1 |
#include "PARAMS.h" |
19 |
|
|
#include "GRID.h" |
20 |
cnh |
1.5 |
|
21 |
adcroft |
1.1 |
C == Routine Arguments == |
22 |
|
|
INTEGER bi,bj,iMin,iMax,jMin,jMax |
23 |
cnh |
1.6 |
_RL KappaRT(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
24 |
|
|
_RL KappaRS(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
25 |
adcroft |
1.1 |
INTEGER myThid |
26 |
cnh |
1.5 |
|
27 |
adcroft |
1.1 |
C == Local variables == |
28 |
|
|
INTEGER i,j,k |
29 |
|
|
_RL a(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
30 |
|
|
_RL b(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
31 |
|
|
_RL c(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
32 |
|
|
_RL ckm1(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
33 |
|
|
_RL bet(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
34 |
cnh |
1.6 |
_RL gam(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
35 |
adcroft |
1.1 |
|
36 |
adcroft |
1.3 |
C *************************************************************** |
37 |
|
|
C ***************** **************** |
38 |
|
|
C ***************** N O T E **************** |
39 |
|
|
C ***************** **************** |
40 |
|
|
C *************************************************************** |
41 |
|
|
C |
42 |
|
|
C The implicit diffusion of SALT currently uses the diffusivity |
43 |
|
|
C of the THETA. |
44 |
cnh |
1.5 |
C ie. KappaRS is ignored. |
45 |
adcroft |
1.3 |
C |
46 |
|
|
C *************************************************************** |
47 |
|
|
C *************************************************************** |
48 |
|
|
|
49 |
cnh |
1.6 |
IF (Nr.GT.1) THEN ! Only need do anything if Nr>1 |
50 |
cnh |
1.5 |
C-- Beginning of forward sweep (top level) |
51 |
adcroft |
1.1 |
DO j=jMin,jMax |
52 |
|
|
DO i=iMin,iMax |
53 |
cnh |
1.6 |
c(i,j)=-deltaTtracer*recip_hFacC(i,j,1,bi,bj)*recip_drF(1) |
54 |
cnh |
1.5 |
& *KappaRT(i,j,2)*recip_drC(2) |
55 |
adcroft |
1.1 |
b(i,j)=1.-c(i,j) |
56 |
|
|
bet(i,j)=0. |
57 |
|
|
IF (b(i,j).NE.0.) bet(i,j)=1. / b(i,j) |
58 |
adcroft |
1.4 |
gTnm1(i,j,1,bi,bj)=gTnm1(i,j,1,bi,bj)*bet(i,j) |
59 |
|
|
gSnm1(i,j,1,bi,bj)=gSnm1(i,j,1,bi,bj)*bet(i,j) |
60 |
adcroft |
1.1 |
ENDDO |
61 |
|
|
ENDDO |
62 |
|
|
ENDIF |
63 |
cnh |
1.5 |
C-- Middle of forward sweep |
64 |
cnh |
1.6 |
IF (Nr.GT.2) THEN |
65 |
|
|
DO k=2,Nr-1 |
66 |
adcroft |
1.1 |
DO j=jMin,jMax |
67 |
|
|
DO i=iMin,iMax |
68 |
|
|
ckm1(i,j)=c(i,j) |
69 |
cnh |
1.6 |
a(i,j)=-deltaTtracer*recip_hFacC(i,j,k,bi,bj)*recip_drF(k) |
70 |
cnh |
1.5 |
& *KappaRT(i,j, k )*recip_drC( k ) |
71 |
cnh |
1.6 |
c(i,j)=-deltaTtracer*recip_hFacC(i,j,k,bi,bj)*recip_drF(k) |
72 |
cnh |
1.5 |
& *KappaRT(i,j,k+1)*recip_drC(k+1) |
73 |
adcroft |
1.1 |
b(i,j)=1.-c(i,j)-a(i,j) |
74 |
|
|
ENDDO |
75 |
|
|
ENDDO |
76 |
|
|
DO j=jMin,jMax |
77 |
|
|
DO i=iMin,iMax |
78 |
|
|
gam(i,j,k)=ckm1(i,j)*bet(i,j) |
79 |
|
|
bet(i,j)=b(i,j)-a(i,j)*gam(i,j,k) |
80 |
|
|
IF (bet(i,j).NE.0.) bet(i,j)=1. / bet(i,j) |
81 |
adcroft |
1.4 |
gTnm1(i,j,k,bi,bj)=(gTnm1(i,j,k,bi,bj) |
82 |
|
|
& -a(i,j)*gTnm1(i,j,k-1,bi,bj))*bet(i,j) |
83 |
|
|
gSnm1(i,j,k,bi,bj)=(gSnm1(i,j,k,bi,bj) |
84 |
|
|
& -a(i,j)*gSnm1(i,j,k-1,bi,bj))*bet(i,j) |
85 |
adcroft |
1.1 |
ENDDO |
86 |
|
|
ENDDO |
87 |
|
|
ENDDO |
88 |
|
|
ENDIF |
89 |
cnh |
1.6 |
IF (Nr.GT.1) THEN |
90 |
cnh |
1.5 |
C-- End of forward sweep (bottom level) |
91 |
adcroft |
1.1 |
DO j=jMin,jMax |
92 |
|
|
DO i=iMin,iMax |
93 |
|
|
ckm1(i,j)=c(i,j) |
94 |
cnh |
1.6 |
a(i,j)=-deltaTtracer*recip_hFacC(i,j,Nr,bi,bj)*recip_drF(Nr) |
95 |
|
|
& *KappaRT(i,j, Nr )*recip_drC( Nr ) |
96 |
adcroft |
1.1 |
b(i,j)=1.-a(i,j) |
97 |
|
|
ENDDO |
98 |
|
|
ENDDO |
99 |
|
|
DO j=jMin,jMax |
100 |
|
|
DO i=iMin,iMax |
101 |
cnh |
1.6 |
gam(i,j,Nr)=ckm1(i,j)*bet(i,j) |
102 |
|
|
bet(i,j)=b(i,j)-a(i,j)*gam(i,j,Nr) |
103 |
adcroft |
1.1 |
IF (bet(i,j).NE.0.) bet(i,j)=1. / bet(i,j) |
104 |
cnh |
1.6 |
gTnm1(i,j,Nr,bi,bj)=(gTnm1(i,j,Nr,bi,bj) |
105 |
|
|
& -a(i,j)*gTnm1(i,j,Nr-1,bi,bj))*bet(i,j) |
106 |
|
|
gSnm1(i,j,Nr,bi,bj)=(gSnm1(i,j,Nr,bi,bj) |
107 |
|
|
& -a(i,j)*gSnm1(i,j,Nr-1,bi,bj))*bet(i,j) |
108 |
adcroft |
1.1 |
ENDDO |
109 |
|
|
ENDDO |
110 |
cnh |
1.5 |
C-- Backward sweep |
111 |
cnh |
1.6 |
DO k=Nr-1,1,-1 |
112 |
adcroft |
1.1 |
DO j=jMin,jMax |
113 |
|
|
DO i=iMin,iMax |
114 |
adcroft |
1.4 |
gTnm1(i,j,k,bi,bj)=gTnm1(i,j,k,bi,bj) |
115 |
|
|
& -gam(i,j,k+1)*gTnm1(i,j,k+1,bi,bj) |
116 |
|
|
gSnm1(i,j,k,bi,bj)=gSnm1(i,j,k,bi,bj) |
117 |
|
|
& -gam(i,j,k+1)*gSnm1(i,j,k+1,bi,bj) |
118 |
adcroft |
1.1 |
ENDDO |
119 |
|
|
ENDDO |
120 |
|
|
ENDDO |
121 |
|
|
ENDIF |
122 |
|
|
|
123 |
|
|
RETURN |
124 |
|
|
END |