1 |
C $Header: /u/gcmpack/MITgcm/model/src/impldiff.F,v 1.17 2001/09/26 18:09:15 cnh Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CPP_OPTIONS.h" |
5 |
|
6 |
CBOP |
7 |
C !ROUTINE: IMPLDIFF |
8 |
C !INTERFACE: |
9 |
SUBROUTINE IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax, |
10 |
I deltaTX,KappaRX,recip_hFac, |
11 |
U gXNm1, |
12 |
I myThid ) |
13 |
C !DESCRIPTION: \bv |
14 |
C *==========================================================* |
15 |
C | S/R IMPLDIFF |
16 |
C | o Solve implicit diffusion equation for vertical |
17 |
C | diffusivity. |
18 |
C *==========================================================* |
19 |
C | o Recoded from 2d intermediate fields to 3d to reduce |
20 |
C | TAMC storage |
21 |
C | o Fixed missing masks for fields a(), c() |
22 |
C *==========================================================* |
23 |
C \ev |
24 |
|
25 |
C !USES: |
26 |
IMPLICIT NONE |
27 |
C == Global data == |
28 |
#include "SIZE.h" |
29 |
#include "DYNVARS.h" |
30 |
#include "EEPARAMS.h" |
31 |
#include "PARAMS.h" |
32 |
#include "GRID.h" |
33 |
#ifdef ALLOW_AUTODIFF_TAMC |
34 |
#include "tamc_keys.h" |
35 |
#endif |
36 |
|
37 |
C !INPUT/OUTPUT PARAMETERS: |
38 |
C == Routine Arguments == |
39 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
40 |
_RL deltaTX |
41 |
_RL KappaRX(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
42 |
_RS recip_hFac(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr,nSx,nSy) |
43 |
_RL gXnm1(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr,nSx,nSy) |
44 |
INTEGER myThid |
45 |
|
46 |
C !LOCAL VARIABLES: |
47 |
C == Local variables == |
48 |
INTEGER i,j,k |
49 |
_RL gYnm1(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr,nSx,nSy) |
50 |
_RL a(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
51 |
_RL b(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
52 |
_RL c(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
53 |
_RL bet(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
54 |
_RL gam(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
55 |
CEOP |
56 |
|
57 |
C-- Initialise |
58 |
DO k=1,Nr |
59 |
DO j=jMin,jMax |
60 |
DO i=iMin,iMax |
61 |
gYNm1(i,j,k,bi,bj) = 0. _d 0 |
62 |
ENDDO |
63 |
ENDDO |
64 |
ENDDO |
65 |
|
66 |
C-- Old aLower |
67 |
DO j=jMin,jMax |
68 |
DO i=iMin,iMax |
69 |
a(i,j,1) = 0. _d 0 |
70 |
ENDDO |
71 |
ENDDO |
72 |
DO k=2,Nr |
73 |
DO j=jMin,jMax |
74 |
DO i=iMin,iMax |
75 |
a(i,j,k) = -deltaTX*recip_hFac(i,j,k,bi,bj)*recip_drF(k) |
76 |
& *KappaRX(i,j, k )*recip_drC( k ) |
77 |
IF (recip_hFac(i,j,k-1,bi,bj).EQ.0.) a(i,j,k)=0. |
78 |
ENDDO |
79 |
ENDDO |
80 |
ENDDO |
81 |
|
82 |
C-- Old aUpper |
83 |
DO k=1,Nr-1 |
84 |
DO j=jMin,jMax |
85 |
DO i=iMin,iMax |
86 |
c(i,j,k) = -deltaTX*recip_hFac(i,j,k,bi,bj)*recip_drF(k) |
87 |
& *KappaRX(i,j,k+1)*recip_drC(k+1) |
88 |
IF (recip_hFac(i,j,k+1,bi,bj).EQ.0.) c(i,j,k)=0. |
89 |
ENDDO |
90 |
ENDDO |
91 |
ENDDO |
92 |
DO j=jMin,jMax |
93 |
DO i=iMin,iMax |
94 |
c(i,j,Nr) = 0. _d 0 |
95 |
ENDDO |
96 |
ENDDO |
97 |
|
98 |
C-- Old aCenter |
99 |
DO k=1,Nr |
100 |
DO j=jMin,jMax |
101 |
DO i=iMin,iMax |
102 |
b(i,j,k) = 1. _d 0 - c(i,j,k) - a(i,j,k) |
103 |
ENDDO |
104 |
ENDDO |
105 |
ENDDO |
106 |
|
107 |
C-- Old and new gam, bet are the same |
108 |
DO k=1,Nr |
109 |
DO j=jMin,jMax |
110 |
DO i=iMin,iMax |
111 |
bet(i,j,k) = 0. _d 0 |
112 |
gam(i,j,k) = 0. _d 0 |
113 |
ENDDO |
114 |
ENDDO |
115 |
ENDDO |
116 |
|
117 |
C-- Only need do anything if Nr>1 |
118 |
IF (Nr.GT.1) THEN |
119 |
|
120 |
k = 1 |
121 |
C-- Beginning of forward sweep (top level) |
122 |
DO j=jMin,jMax |
123 |
DO i=iMin,iMax |
124 |
IF (b(i,j,1).NE.0.) bet(i,j,1) = 1. _d 0 / b(i,j,1) |
125 |
ENDDO |
126 |
ENDDO |
127 |
|
128 |
ENDIF |
129 |
|
130 |
C-- Middle of forward sweep |
131 |
IF (Nr.GT.2) THEN |
132 |
|
133 |
CADJ loop = sequential |
134 |
DO k=2,Nr |
135 |
|
136 |
DO j=jMin,jMax |
137 |
DO i=iMin,iMax |
138 |
gam(i,j,k) = c(i,j,k-1)*bet(i,j,k-1) |
139 |
IF ( ( b(i,j,k) - a(i,j,k)*gam(i,j,k) ) .NE. 0.) |
140 |
& bet(i,j,k) = 1. _d 0 / ( b(i,j,k) - a(i,j,k)*gam(i,j,k) ) |
141 |
ENDDO |
142 |
ENDDO |
143 |
|
144 |
ENDDO |
145 |
|
146 |
ENDIF |
147 |
|
148 |
|
149 |
DO j=jMin,jMax |
150 |
DO i=iMin,iMax |
151 |
gYNm1(i,j,1,bi,bj) = gXNm1(i,j,1,bi,bj)*bet(i,j,1) |
152 |
ENDDO |
153 |
ENDDO |
154 |
DO k=2,Nr |
155 |
DO j=jMin,jMax |
156 |
DO i=iMin,iMax |
157 |
gYnm1(i,j,k,bi,bj) = bet(i,j,k)* |
158 |
& (gXnm1(i,j,k,bi,bj) - a(i,j,k)*gYnm1(i,j,k-1,bi,bj)) |
159 |
ENDDO |
160 |
ENDDO |
161 |
ENDDO |
162 |
|
163 |
|
164 |
C-- Backward sweep |
165 |
CADJ loop = sequential |
166 |
DO k=Nr-1,1,-1 |
167 |
DO j=jMin,jMax |
168 |
DO i=iMin,iMax |
169 |
gYnm1(i,j,k,bi,bj)=gYnm1(i,j,k,bi,bj) |
170 |
& -gam(i,j,k+1)*gYnm1(i,j,k+1,bi,bj) |
171 |
ENDDO |
172 |
ENDDO |
173 |
ENDDO |
174 |
|
175 |
DO k=1,Nr |
176 |
DO j=jMin,jMax |
177 |
DO i=iMin,iMax |
178 |
gXnm1(i,j,k,bi,bj)=gYnm1(i,j,k,bi,bj) |
179 |
ENDDO |
180 |
ENDDO |
181 |
ENDDO |
182 |
|
183 |
RETURN |
184 |
END |