/[MITgcm]/MITgcm/model/src/ini_cg3d.F
ViewVC logotype

Contents of /MITgcm/model/src/ini_cg3d.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.18 - (show annotations) (download)
Tue Nov 8 06:06:07 2005 UTC (18 years, 7 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57y_post, checkpoint57y_pre, checkpoint57x_post
Changes since 1.17: +58 -50 lines
- initialise all variables over the full range of indices;
- fix exchange calls for CS-grid ; fix mask for p-coordinate

1 C $Header: /u/gcmpack/MITgcm/model/src/ini_cg3d.F,v 1.17 2005/05/31 14:49:38 adcroft Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 CBOP
8 C !ROUTINE: INI_CG3D
9 C !INTERFACE:
10 SUBROUTINE INI_CG3D( myThid )
11 C !DESCRIPTION: \bv
12 C *==========================================================*
13 C | SUBROUTINE INI_CG3D
14 C | o Initialise 3d conjugate gradient solver operators.
15 C *==========================================================*
16 C | These arrays are purely a function of the basin geom.
17 C | We set then here once and them use then repeatedly.
18 C *==========================================================*
19 C \ev
20
21 C !USES:
22 IMPLICIT NONE
23 C === Global variables ===
24 #include "SIZE.h"
25 #include "EEPARAMS.h"
26 #include "PARAMS.h"
27 #include "GRID.h"
28 #include "CG3D.h"
29 #include "SOLVE_FOR_PRESSURE3D.h"
30 #ifdef ALLOW_OBCS
31 #include "OBCS.h"
32 #endif
33
34 C !INPUT/OUTPUT PARAMETERS:
35 C === Routine arguments ===
36 C myThid - Thread no. that called this routine.
37 INTEGER myThid
38
39 #ifdef ALLOW_NONHYDROSTATIC
40
41 C !LOCAL VARIABLES:
42 C === Local variables ===
43 C bi,bj - Loop counters
44 C I,J,K - Loop counters
45 C faceArea - Temporary used to hold cell face areas.
46 C myNorm - Work variable used in clculating normalisation factor
47 CHARACTER*(MAX_LEN_MBUF) msgBuf
48 INTEGER bi, bj
49 INTEGER I, J, K
50 _RL faceArea
51 _RS myNorm
52 _RL theRecip_Dr
53 _RL aU, aL, aW, aE, aN, aS, aC
54 CEOP
55
56 CcnhDebugStarts
57 c _RL phi(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
58 CcnhDebugEnds
59
60 C-- Initialise to zero over the full range of indices
61 DO bj=myByLo(myThid),myByHi(myThid)
62 DO bi=myBxLo(myThid),myBxHi(myThid)
63 DO K=1,Nr
64 C- From common bloc CG3D_R:
65 DO J=1-OLy,sNy+OLy
66 DO I=1-OLx,sNx+OLx
67 aW3d(I,J,K,bi,bj) = 0.
68 aS3d(I,J,K,bi,bj) = 0.
69 aV3d(I,J,K,bi,bj) = 0.
70 zMC (I,J,K,bi,bj) = 0.
71 zML (I,J,K,bi,bj) = 0.
72 zMU (I,J,K,bi,bj) = 0.
73 ENDDO
74 ENDDO
75 C- From common bloc CG3D_WK_R:
76 DO J=0,sNy+1
77 DO I=0,sNx+1
78 cg3d_q(I,J,K,bi,bj) = 0.
79 cg3d_r(I,J,K,bi,bj) = 0.
80 cg3d_s(I,J,K,bi,bj) = 0.
81 ENDDO
82 ENDDO
83 C- From common bloc SFP3D_COMMON_R:
84 DO J=1-OLy,sNy+OLy
85 DO I=1-OLx,sNx+OLx
86 c cg3d_x(I,J,K,bi,bj) = 0.
87 cg3d_b(I,J,K,bi,bj) = 0.
88 ENDDO
89 ENDDO
90 ENDDO
91 ENDDO
92 ENDDO
93
94 C-- Initialise laplace operator
95 C aW3d: Ax/dX
96 C aS3d: Ay/dY
97 C aV3d: Ar/dR
98 myNorm = 0. _d 0
99 DO bj=myByLo(myThid),myByHi(myThid)
100 DO bi=myBxLo(myThid),myBxHi(myThid)
101 DO K=1,Nr
102 DO J=1,sNy
103 DO I=1,sNx
104 faceArea = _dyG(I,J,bi,bj)*drF(K)
105 & *_hFacW(I,J,K,bi,bj)
106 aW3d(I,J,K,bi,bj) = faceArea*recip_dxC(I,J,bi,bj)
107 faceArea = _dxG(I,J,bi,bj)*drF(K)
108 & *_hFacS(I,J,K,bi,bj)
109 aS3d(I,J,K,bi,bj) = faceArea*recip_dyC(I,J,bi,bj)
110 myNorm = MAX(ABS(aW3d(I,J,K,bi,bj)),myNorm)
111 myNorm = MAX(ABS(aS3d(I,J,K,bi,bj)),myNorm)
112 ENDDO
113 ENDDO
114 ENDDO
115 DO K=1,1
116 DO J=1,sNy
117 DO I=1,sNx
118 aV3d(I,J,K,bi,bj) = 0.
119 myNorm = MAX(ABS(aV3d(I,J,K,bi,bj)),myNorm)
120 ENDDO
121 ENDDO
122 ENDDO
123 DO K=2,Nr
124 DO J=1,sNy
125 DO I=1,sNx
126 faceArea = _rA(I,J,bi,bj)*maskC(I,J, K ,bi,bj)
127 & *maskC(I,J,K-1,bi,bj)
128 theRecip_Dr = recip_drC(K)
129 c theRecip_Dr =
130 caja & drF(K )*_hFacC(i,j,k ,bi,bj)*0.5
131 caja & +drF(K-1)*_hFacC(i,j,k-1,bi,bj)*0.5
132 c IF ( theRecip_Dr .NE. 0. )
133 c & theRecip_Dr = 1. _d 0/theRecip_Dr
134 aV3d(I,J,K,bi,bj) = faceArea*theRecip_Dr
135 & *horiVertRatio*horiVertRatio *nh_Am2
136 myNorm = MAX(ABS(aV3d(I,J,K,bi,bj)),myNorm)
137 ENDDO
138 ENDDO
139 ENDDO
140 #ifdef ALLOW_OBCS
141 IF ( useOBCS ) THEN
142 DO K=1,Nr
143 DO I=1,sNx
144 IF (OB_Jn(I,bi,bj).NE.0) THEN
145 aS3d(I,OB_Jn(I,bi,bj),K,bi,bj)=0.
146 aS3d(I,OB_Jn(I,bi,bj)+1,K,bi,bj)=0.
147 aW3d(I,OB_Jn(I,bi,bj),K,bi,bj)=0.
148 aW3d(I+1,OB_Jn(I,bi,bj),K,bi,bj)=0.
149 aV3d(I,OB_Jn(I,bi,bj),K,bi,bj)=0.
150 ENDIF
151 IF (OB_Js(I,bi,bj).NE.0) THEN
152 aS3d(I,OB_Js(I,bi,bj)+1,K,bi,bj)=0.
153 aS3d(I,OB_Js(I,bi,bj),K,bi,bj)=0.
154 aW3d(I,OB_Js(I,bi,bj),K,bi,bj)=0.
155 aW3d(I+1,OB_Js(I,bi,bj),K,bi,bj)=0.
156 aV3d(I,OB_Js(I,bi,bj),K,bi,bj)=0.
157 ENDIF
158 ENDDO
159 DO J=1,sNy
160 IF (OB_Ie(J,bi,bj).NE.0) THEN
161 aW3d(OB_Ie(J,bi,bj),J,K,bi,bj)=0.
162 aW3d(OB_Ie(J,bi,bj)+1,J,K,bi,bj)=0.
163 aS3d(OB_Ie(J,bi,bj),J,K,bi,bj)=0.
164 aS3d(OB_Ie(J,bi,bj),J+1,K,bi,bj)=0.
165 aV3d(OB_Ie(J,bi,bj),J,K,bi,bj)=0.
166 ENDIF
167 IF (OB_Iw(J,bi,bj).NE.0) THEN
168 aW3d(OB_Iw(J,bi,bj)+1,J,K,bi,bj)=0.
169 aW3d(OB_Iw(J,bi,bj),J,K,bi,bj)=0.
170 aS3d(OB_Iw(J,bi,bj),J,K,bi,bj)=0.
171 aS3d(OB_Iw(J,bi,bj),J+1,K,bi,bj)=0.
172 aV3d(OB_Iw(J,bi,bj),J,K,bi,bj)=0.
173 ENDIF
174 ENDDO
175 ENDDO
176 ENDIF
177 #endif
178 ENDDO
179 ENDDO
180 _GLOBAL_MAX_R4( myNorm, myThid )
181 IF ( myNorm .NE. 0. _d 0 ) THEN
182 myNorm = 1. _d 0/myNorm
183 ELSE
184 myNorm = 1. _d 0
185 ENDIF
186 cg3dNorm = myNorm
187 _BEGIN_MASTER( myThid )
188 CcnhDebugStarts
189 WRITE(msgBuf,'(A,E40.25)') '// CG3D normalisation factor = '
190 & , cg3dNorm
191 CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
192 WRITE(msgBuf,*) ' '
193 CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
194 CcnhDebugEnds
195 _END_MASTER( myThid )
196 DO bj=myByLo(myThid),myByHi(myThid)
197 DO bi=myBxLo(myThid),myBxHi(myThid)
198 DO K=1,Nr
199 DO J=1,sNy
200 DO I=1,sNx
201 aW3d(I,J,K,bi,bj) = aW3d(I,J,K,bi,bj)*myNorm
202 aS3d(I,J,K,bi,bj) = aS3d(I,J,K,bi,bj)*myNorm
203 aV3d(I,J,K,bi,bj) = aV3d(I,J,K,bi,bj)*myNorm
204 ENDDO
205 ENDDO
206 ENDDO
207 ENDDO
208 ENDDO
209
210 C-- Update overlap regions
211 c _EXCH_XYZ_R4(aW3d, myThid)
212 c _EXCH_XYZ_R4(aS3d, myThid)
213 CALL EXCH_UV_XYZ_RS(aW3d,aS3d,.FALSE.,myThid)
214 _EXCH_XYZ_R4(aV3d, myThid)
215 CcnhDebugStarts
216 C CALL PLOT_FIELD_XYZRS( aW3d, 'AW3D INI_CG3D.1' , Nr, 1, myThid )
217 C CALL PLOT_FIELD_XYZRS( aS3d, 'AS3D INI_CG3D.1' , Nr, 1, myThid )
218 CcnhDebugEnds
219
220 C-- Initialise preconditioner
221 C For now PC is just the identity. Change to
222 C be LU factorization of d2/dz2 later. Note
223 C check for consistency with S/R CG3D before
224 C assuming zML is lower and zMU is upper!
225 DO bj=myByLo(myThid),myByHi(myThid)
226 DO bi=myBxLo(myThid),myBxHi(myThid)
227 DO K=1,Nr
228 DO J=1,sNy
229 DO I=1,sNx
230 aW = aW3d(I ,J,K,bi,bj)
231 aE = aW3d(I+1,J,K,bi,bj)
232 aN = aS3d(I,J+1,K,bi,bj)
233 aS = aS3d(I,J ,K,bi,bj)
234 IF ( K .NE. 1 ) THEN
235 aU = aV3d(I,J,K,bi,bj)
236 ELSE
237 aU = 0
238 ENDIF
239 IF ( K .NE. Nr ) THEN
240 aL = aV3d(I,J,K+1,bi,bj)
241 ELSE
242 aL = 0
243 ENDIF
244 aC = -aW-aE-aN-aS-aU-aL
245 IF ( K .EQ. 1 .AND. aC .NE. 0. ) THEN
246 aC = aC
247 & -freeSurfFac*myNorm*(horiVertRatio/gravity)*
248 & rA(I,J,bi,bj)/deltaTMom/deltaTMom
249 ENDIF
250 IF ( aC .NE. 0. ) THEN
251 zMC(i,j,k,bi,bj) = aC
252 zMU(i,j,k,bi,bj) = aL
253 zML(i,j,k,bi,bj) = aU
254 CcnhDebugStarts
255 C zMC(i,j,k,bi,bj) = 1.
256 C zMU(i,j,k,bi,bj) = 0.
257 C zML(i,j,k,bi,bj) = 0.
258 CcnhDebugEnds
259 ELSE
260 zMC(i,j,k,bi,bj) = 1. _d 0
261 zMU(i,j,k,bi,bj) = 0.
262 zML(i,j,k,bi,bj) = 0.
263 ENDIF
264 ENDDO
265 ENDDO
266 ENDDO
267 DO J=1,sNy
268 DO I=1,sNx
269 zMC(i,j,1,bi,bj)=
270 & 1. _d 0 / zMC(i,j,1,bi,bj)
271 zMU(i,j,1,bi,bj)=
272 & zMU(i,j,1,bi,bj)*zMC(i,j,1,bi,bj)
273 ENDDO
274 ENDDO
275 DO K=2,Nr
276 DO J=1,sNy
277 DO I=1,sNx
278 zMC(i,j,k,bi,bj) = 1. _d 0 /
279 & (zMC(i,j,k,bi,bj)-zML(i,j,k,bi,bj)*zMU(i,j,k-1,bi,bj))
280 zMU(i,j,k,bi,bj)=zMU(i,j,k,bi,bj)*zMC(i,j,k,bi,bj)
281 ENDDO
282 ENDDO
283 ENDDO
284 DO K=1,Nr
285 DO J=1,sNy
286 DO I=1,sNx
287 aW = aW3d(I ,J,K,bi,bj)
288 aE = aW3d(I+1,J,K,bi,bj)
289 aN = aS3d(I,J+1,K,bi,bj)
290 aS = aS3d(I,J ,K,bi,bj)
291 IF ( K .NE. 1 ) THEN
292 aU = aV3d(I,J,K-1,bi,bj)
293 ELSE
294 aU = 0
295 ENDIF
296 IF ( K .NE. Nr ) THEN
297 aL = aV3d(I,J,K+1,bi,bj)
298 ELSE
299 aL = 0
300 ENDIF
301 aC = -aW-aE-aN-aS-aU-aL
302 IF ( aC .EQ. 0. ) THEN
303 zMC(i,j,k,bi,bj) = 1.
304 zML(i,j,k,bi,bj) = 0.
305 zMU(i,j,k,bi,bj) = 0.
306 CcnhDebugStarts
307 C ELSE
308 C zMC(i,j,k,bi,bj) = 1.
309 C zML(i,j,k,bi,bj) = 0.
310 C zMU(i,j,k,bi,bj) = 0.
311 CcnhDEbugEnds
312 ENDIF
313 ENDDO
314 ENDDO
315 ENDDO
316 ENDDO
317 ENDDO
318 C-- Update overlap regions
319 _EXCH_XYZ_R4(zMC, myThid)
320 _EXCH_XYZ_R4(zML, myThid)
321 _EXCH_XYZ_R4(zMU, myThid)
322
323 CcnhDebugStarts
324 c DO k=1,Nr
325 c DO j=1-OLy,sNy+OLy
326 c DO i=1-OLx,sNx+OLx
327 c phi(i,j,1,1) = zMc(i,j,k,1,1)
328 c ENDDO
329 c ENDDO
330 C CALL PLOT_FIELD_XYRS( phi, 'zMC INI_CG3D.1' , 1, myThid )
331 c ENDDO
332 C CALL PLOT_FIELD_XYRS( zMU, 'zMU INI_CG3D.1' , Nr, 1, myThid )
333 C CALL PLOT_FIELD_XYRS( zML, 'zML INI_CG3D.1' , Nr, 1, myThid )
334 CcnhDebugEnds
335
336 #endif /* ALLOW_NONHYDROSTATIC */
337
338 RETURN
339 END

  ViewVC Help
Powered by ViewVC 1.1.22