/[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.15 - (show annotations) (download)
Sun Feb 4 14:38:47 2001 UTC (23 years, 4 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint38, pre38tag1, c37_adj, pre38-close, checkpoint37, checkpoint36, checkpoint35
Branch point for: pre38
Changes since 1.14: +2 -1 lines
Made sure each .F and .h file had
the CVS keywords Header and Name at its start.
Most had header but very few currently have Name, so
lots of changes!

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

  ViewVC Help
Powered by ViewVC 1.1.22