/[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.22 - (show annotations) (download)
Fri Oct 22 20:14:46 2004 UTC (19 years, 7 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint56b_post, checkpoint57, checkpoint56, checkpoint55i_post, checkpoint57a_post, checkpoint55j_post, checkpoint56a_post, checkpoint56c_post, checkpoint57a_pre
Changes since 1.21: +2 -2 lines
fix the Nr=1 Pb (2nd time) in a different way, so that it will not be removed

1 C $Header: /u/gcmpack/MITgcm/model/src/impldiff.F,v 1.21 2004/02/18 22:23:18 heimbach 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 cph(
58 cph Not good for TAF: may create irreducible control flow graph
59 cph IF (Nr.LE.1) RETURN
60 cph)
61
62 C-- Initialise
63 DO k=1,Nr
64 DO j=jMin,jMax
65 DO i=iMin,iMax
66 gYNm1(i,j,k,bi,bj) = 0. _d 0
67 ENDDO
68 ENDDO
69 ENDDO
70
71 C-- Old aLower
72 DO j=jMin,jMax
73 DO i=iMin,iMax
74 a(i,j,1) = 0. _d 0
75 ENDDO
76 ENDDO
77 DO k=2,Nr
78 DO j=jMin,jMax
79 DO i=iMin,iMax
80 a(i,j,k) = -deltaTX*recip_hFac(i,j,k,bi,bj)*recip_drF(k)
81 & *KappaRX(i,j, k )*recip_drC( k )
82 IF (recip_hFac(i,j,k-1,bi,bj).EQ.0.) a(i,j,k)=0.
83 ENDDO
84 ENDDO
85 ENDDO
86
87 C-- Old aUpper
88 DO k=1,Nr-1
89 DO j=jMin,jMax
90 DO i=iMin,iMax
91 c(i,j,k) = -deltaTX*recip_hFac(i,j,k,bi,bj)*recip_drF(k)
92 & *KappaRX(i,j,k+1)*recip_drC(k+1)
93 IF (recip_hFac(i,j,k+1,bi,bj).EQ.0.) c(i,j,k)=0.
94 ENDDO
95 ENDDO
96 ENDDO
97 DO j=jMin,jMax
98 DO i=iMin,iMax
99 c(i,j,Nr) = 0. _d 0
100 ENDDO
101 ENDDO
102
103 C-- Old aCenter
104 DO k=1,Nr
105 DO j=jMin,jMax
106 DO i=iMin,iMax
107 b(i,j,k) = 1. _d 0 - c(i,j,k) - a(i,j,k)
108 ENDDO
109 ENDDO
110 ENDDO
111
112 C-- Old and new gam, bet are the same
113 DO k=1,Nr
114 DO j=jMin,jMax
115 DO i=iMin,iMax
116 bet(i,j,k) = 1. _d 0
117 gam(i,j,k) = 0. _d 0
118 ENDDO
119 ENDDO
120 ENDDO
121
122 C-- Only need do anything if Nr>1
123 IF (Nr.GT.1) THEN
124
125 k = 1
126 C-- Beginning of forward sweep (top level)
127 DO j=jMin,jMax
128 DO i=iMin,iMax
129 IF (b(i,j,1).NE.0.) bet(i,j,1) = 1. _d 0 / b(i,j,1)
130 ENDDO
131 ENDDO
132
133 ENDIF
134
135 C-- Middle of forward sweep
136 IF (Nr.GE.2) THEN
137
138 CADJ loop = sequential
139 DO k=2,Nr
140
141 DO j=jMin,jMax
142 DO i=iMin,iMax
143 gam(i,j,k) = c(i,j,k-1)*bet(i,j,k-1)
144 IF ( ( b(i,j,k) - a(i,j,k)*gam(i,j,k) ) .NE. 0.)
145 & bet(i,j,k) = 1. _d 0 / ( b(i,j,k) - a(i,j,k)*gam(i,j,k) )
146 ENDDO
147 ENDDO
148
149 ENDDO
150
151 ENDIF
152
153
154 DO j=jMin,jMax
155 DO i=iMin,iMax
156 gYNm1(i,j,1,bi,bj) = gXNm1(i,j,1,bi,bj)*bet(i,j,1)
157 ENDDO
158 ENDDO
159 DO k=2,Nr
160 DO j=jMin,jMax
161 DO i=iMin,iMax
162 gYnm1(i,j,k,bi,bj) = bet(i,j,k)*
163 & (gXnm1(i,j,k,bi,bj) - a(i,j,k)*gYnm1(i,j,k-1,bi,bj))
164 ENDDO
165 ENDDO
166 ENDDO
167
168
169 C-- Backward sweep
170 CADJ loop = sequential
171 DO k=Nr-1,1,-1
172 DO j=jMin,jMax
173 DO i=iMin,iMax
174 gYnm1(i,j,k,bi,bj)=gYnm1(i,j,k,bi,bj)
175 & -gam(i,j,k+1)*gYnm1(i,j,k+1,bi,bj)
176 ENDDO
177 ENDDO
178 ENDDO
179
180 DO k=1,Nr
181 DO j=jMin,jMax
182 DO i=iMin,iMax
183 gXnm1(i,j,k,bi,bj)=gYnm1(i,j,k,bi,bj)
184 ENDDO
185 ENDDO
186 ENDDO
187
188 RETURN
189 END

  ViewVC Help
Powered by ViewVC 1.1.22