/[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.13 - (show annotations) (download)
Wed Nov 29 22:38:31 2000 UTC (23 years, 6 months ago) by adcroft
Branch: MAIN
CVS Tags: branch-atmos-merge-freeze, branch-atmos-merge-start, branch-atmos-merge-shapiro, checkpoint33, branch-atmos-merge-zonalfilt, branch-atmos-merge-phase5, branch-atmos-merge-phase4, branch-atmos-merge-phase7, branch-atmos-merge-phase6, branch-atmos-merge-phase1, branch-atmos-merge-phase3, branch-atmos-merge-phase2
Branch point for: branch-atmos-merge
Changes since 1.12: +5 -1 lines
Bug fix for lower boundaries over topography of non-zero
height. Thanks to CAE for spotting it.

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

  ViewVC Help
Powered by ViewVC 1.1.22