/[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.26 - (show annotations) (download)
Sat Jan 23 00:04:03 2010 UTC (14 years, 4 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62c, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62w, checkpoint62x, checkpoint62b
Changes since 1.25: +10 -1 lines
add NH free-surface formulation (selectNHfreeSurf=1) (not fully tested)

1 C $Header: /u/gcmpack/MITgcm/model/src/ini_cg3d.F,v 1.25 2009/11/29 03:16:10 jmc 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 "SURFACE.h"
29 #include "CG3D.h"
30 #include "SOLVE_FOR_PRESSURE3D.h"
31 #ifdef ALLOW_OBCS
32 #include "OBCS.h"
33 #endif
34
35 C !INPUT/OUTPUT PARAMETERS:
36 C === Routine arguments ===
37 C myThid :: My Thread Id number
38 INTEGER myThid
39
40 #ifdef ALLOW_NONHYDROSTATIC
41
42 C !LOCAL VARIABLES:
43 C === Local variables ===
44 C bi,bj :: tile indices
45 C i,j,k :: Loop counters
46 C faceArea :: Temporary used to hold cell face areas.
47 C myNorm :: Work variable used in clculating normalisation factor
48 CHARACTER*(MAX_LEN_MBUF) msgBuf
49 INTEGER bi, bj
50 INTEGER i, j, k, ks
51 _RL faceArea
52 _RS myNorm
53 _RL theRecip_Dr
54 _RL aU, aL, aW, aE, aN, aS
55 _RL tmpFac, nh_Fac, igwFac
56 _RL tmpSurf
57 CEOP
58
59 CcnhDebugStarts
60 c _RL phi(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
61 CcnhDebugEnds
62
63 C-- Initialise to zero over the full range of indices
64 DO bj=myByLo(myThid),myByHi(myThid)
65 DO bi=myBxLo(myThid),myBxHi(myThid)
66 DO k=1,Nr
67 C- From common bloc CG3D_R:
68 DO j=1-OLy,sNy+OLy
69 DO i=1-OLx,sNx+OLx
70 aW3d(i,j,k,bi,bj) = 0.
71 aS3d(i,j,k,bi,bj) = 0.
72 aV3d(i,j,k,bi,bj) = 0.
73 aC3d(i,j,k,bi,bj) = 0.
74 zMC (i,j,k,bi,bj) = 0.
75 zML (i,j,k,bi,bj) = 0.
76 zMU (i,j,k,bi,bj) = 0.
77 ENDDO
78 ENDDO
79 C- From common bloc CG3D_WK_R:
80 DO j=0,sNy+1
81 DO i=0,sNx+1
82 cg3d_q(i,j,k,bi,bj) = 0.
83 cg3d_r(i,j,k,bi,bj) = 0.
84 cg3d_s(i,j,k,bi,bj) = 0.
85 ENDDO
86 ENDDO
87 C- From common bloc SFP3D_COMMON_R:
88 DO j=1-OLy,sNy+OLy
89 DO i=1-OLx,sNx+OLx
90 cg3d_b(i,j,k,bi,bj) = 0.
91 ENDDO
92 ENDDO
93 ENDDO
94 ENDDO
95 ENDDO
96
97 nh_Fac = 0.
98 igwFac = 0.
99 IF ( nonHydrostatic
100 & .AND. nh_Am2.NE.0. ) nh_Fac = 1. _d 0 / nh_Am2
101 IF ( implicitIntGravWave ) igwFac = 1. _d 0
102
103 IF ( use3Dsolver ) THEN
104 C-- Initialise laplace operator
105 C aW3d: Ax/dX
106 C aS3d: Ay/dY
107 C aV3d: Ar/dR
108 myNorm = 0. _d 0
109 DO bj=myByLo(myThid),myByHi(myThid)
110 DO bi=myBxLo(myThid),myBxHi(myThid)
111 DO k=1,Nr
112 DO j=1,sNy
113 DO i=1,sNx+1
114 faceArea = _dyG(i,j,bi,bj)*drF(k)
115 & *_hFacW(i,j,k,bi,bj)
116 aW3d(i,j,k,bi,bj) = faceArea*recip_dxC(i,j,bi,bj)
117 & *implicitNHPress*implicDiv2DFlow
118 myNorm = MAX(ABS(aW3d(i,j,k,bi,bj)),myNorm)
119 ENDDO
120 ENDDO
121 C deep-model: *deepFacC (faceArea), /deepFacC (recip_dx,y): => no net effect
122 DO j=1,sNy+1
123 DO i=1,sNx
124 faceArea = _dxG(i,j,bi,bj)*drF(K)
125 & *_hFacS(i,j,k,bi,bj)
126 aS3d(i,j,k,bi,bj) = faceArea*recip_dyC(i,j,bi,bj)
127 & *implicitNHPress*implicDiv2DFlow
128 myNorm = MAX(ABS(aS3d(i,j,k,bi,bj)),myNorm)
129 ENDDO
130 ENDDO
131 ENDDO
132 DO k=1,1
133 DO j=1,sNy
134 DO i=1,sNx
135 aV3d(i,j,k,bi,bj) = 0.
136 ENDDO
137 ENDDO
138 ENDDO
139 DO k=2,Nr
140 tmpFac = nh_Fac*rVel2wUnit(k)*rVel2wUnit(k)
141 & + igwFac*dBdrRef(k)*deltaTMom*dTtracerLev(k)
142 IF (tmpFac.GT.0. ) tmpFac = 1. _d 0 / tmpFac
143 DO j=1,sNy
144 DO i=1,sNx
145 faceArea = _rA(i,j,bi,bj)*maskC(i,j, k ,bi,bj)
146 & *maskC(i,j,k-1,bi,bj)
147 & *deepFac2F(k)
148 theRecip_Dr = recip_drC(k)
149 c theRecip_Dr =
150 caja & drF(k )*_hFacC(i,j,k ,bi,bj)*0.5
151 caja & +drF(k-1)*_hFacC(i,j,k-1,bi,bj)*0.5
152 c IF ( theRecip_Dr .NE. 0. )
153 c & theRecip_Dr = 1. _d 0/theRecip_Dr
154 aV3d(i,j,k,bi,bj) = faceArea*theRecip_Dr*tmpFac
155 & *implicitNHPress*implicDiv2DFlow
156 myNorm = MAX(ABS(aV3d(i,j,k,bi,bj)),myNorm)
157 ENDDO
158 ENDDO
159 ENDDO
160 #ifdef ALLOW_OBCS
161 IF ( useOBCS ) THEN
162 DO k=1,Nr
163 DO i=1,sNx
164 IF (OB_Jn(i,bi,bj).NE.0) THEN
165 aS3d( i, OB_Jn(i,bi,bj), k,bi,bj) = 0.
166 aS3d( i,OB_Jn(i,bi,bj)+1,k,bi,bj) = 0.
167 aW3d( i, OB_Jn(i,bi,bj), k,bi,bj) = 0.
168 aW3d(i+1,OB_Jn(i,bi,bj), k,bi,bj) = 0.
169 aV3d( i, OB_Jn(i,bi,bj), k,bi,bj) = 0.
170 ENDIF
171 IF (OB_Js(i,bi,bj).NE.0) THEN
172 aS3d( i,OB_Js(i,bi,bj)+1,k,bi,bj) = 0.
173 aS3d( i, OB_Js(i,bi,bj), k,bi,bj) = 0.
174 aW3d( i, OB_Js(i,bi,bj), k,bi,bj) = 0.
175 aW3d(i+1,OB_Js(i,bi,bj), k,bi,bj) = 0.
176 aV3d( i, OB_Js(i,bi,bj), k,bi,bj) = 0.
177 ENDIF
178 ENDDO
179 DO j=1,sNy
180 IF (OB_Ie(j,bi,bj).NE.0) THEN
181 aW3d( OB_Ie(j,bi,bj), j, k,bi,bj) = 0.
182 aW3d(OB_Ie(j,bi,bj)+1,j, k,bi,bj) = 0.
183 aS3d( OB_Ie(j,bi,bj), j, k,bi,bj) = 0.
184 aS3d( OB_Ie(j,bi,bj),j+1,k,bi,bj) = 0.
185 aV3d( OB_Ie(j,bi,bj), j, k,bi,bj) = 0.
186 ENDIF
187 IF (OB_Iw(j,bi,bj).NE.0) THEN
188 aW3d(OB_Iw(j,bi,bj)+1,j, k,bi,bj) = 0.
189 aW3d( OB_Iw(j,bi,bj), j, k,bi,bj) = 0.
190 aS3d( OB_Iw(j,bi,bj), j, k,bi,bj) = 0.
191 aS3d( OB_Iw(j,bi,bj),j+1,k,bi,bj) = 0.
192 aV3d( OB_Iw(j,bi,bj), j, k,bi,bj) = 0.
193 ENDIF
194 ENDDO
195 ENDDO
196 ENDIF
197 #endif
198 ENDDO
199 ENDDO
200 _GLOBAL_MAX_RS( myNorm, myThid )
201 IF ( myNorm .NE. 0. _d 0 ) THEN
202 myNorm = 1. _d 0/myNorm
203 ELSE
204 myNorm = 1. _d 0
205 ENDIF
206
207 _BEGIN_MASTER( myThid )
208 C-- set global parameter in common block:
209 cg3dNorm = myNorm
210 WRITE(msgBuf,'(2A,1PE23.16)') 'INI_CG3D: ',
211 & 'CG3D normalisation factor = ', cg3dNorm
212 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
213 & SQUEEZE_RIGHT, myThid )
214 WRITE(msgBuf,*) ' '
215 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
216 & SQUEEZE_RIGHT, myThid )
217 _END_MASTER( myThid )
218
219 DO bj=myByLo(myThid),myByHi(myThid)
220 DO bi=myBxLo(myThid),myBxHi(myThid)
221 C- Set solver main diagonal term
222 DO k=1,Nr
223 DO j=1,sNy
224 DO i=1,sNx
225 aW = aW3d( i, j, k, bi,bj)
226 aE = aW3d(i+1,j, k, bi,bj)
227 aN = aS3d( i,j+1,k, bi,bj)
228 aS = aS3d( i, j, k, bi,bj)
229 aU = aV3d( i, j, k, bi,bj)
230 IF ( k .NE. Nr ) THEN
231 aL = aV3d(i, j,k+1,bi,bj)
232 ELSE
233 aL = 0.
234 ENDIF
235 aC3d(i,j,k,bi,bj) = -aW-aE-aN-aS-aU-aL
236 ENDDO
237 ENDDO
238 ENDDO
239 C- Add free-surface source term
240 IF ( select_rStar .EQ. 0 ) THEN
241 DO j=1,sNy
242 DO i=1,sNx
243 tmpSurf = 1.
244 IF ( selectNHfreeSurf.GE.1 ) THEN
245 tmpSurf = deltaTMom*deltaTfreesurf
246 & *Bo_surf(i,j,bi,bj)*recip_drC(1)
247 & *implicitNHPress*implicDiv2DFlow
248 tmpSurf = tmpSurf / (1. _d 0 + tmpSurf )
249 ENDIF
250 ks = ksurfC(i,j,bi,bj)
251 IF ( ks.LE.Nr ) THEN
252 aC3d(i,j,ks,bi,bj) = aC3d(i,j,ks,bi,bj)
253 & - freeSurfFac*recip_Bo(i,j,bi,bj)
254 & *rA(i,j,bi,bj)*deepFac2F(ks)/deltaTMom/deltaTfreesurf
255 & *tmpSurf
256 ENDIF
257 ENDDO
258 ENDDO
259 ENDIF
260 C- Matrix solver normalisation
261 DO k=1,Nr
262 DO j=1,sNy
263 DO i=1,sNx
264 aW3d(i,j,k,bi,bj) = aW3d(i,j,k,bi,bj)*myNorm
265 aS3d(i,j,k,bi,bj) = aS3d(i,j,k,bi,bj)*myNorm
266 aV3d(i,j,k,bi,bj) = aV3d(i,j,k,bi,bj)*myNorm
267 aC3d(i,j,k,bi,bj) = aC3d(i,j,k,bi,bj)*myNorm
268 ENDDO
269 ENDDO
270 ENDDO
271 ENDDO
272 ENDDO
273
274 C-- Update overlap regions
275 CALL EXCH_UV_XYZ_RS(aW3d,aS3d,.FALSE.,myThid)
276 _EXCH_XYZ_RS(aV3d, myThid)
277 _EXCH_XYZ_RS(aC3d, myThid)
278 CcnhDebugStarts
279 C CALL PLOT_FIELD_XYZRS( aW3d, 'AW3D INI_CG3D.1' , Nr, 1, myThid )
280 C CALL PLOT_FIELD_XYZRS( aS3d, 'AS3D INI_CG3D.1' , Nr, 1, myThid )
281 CcnhDebugEnds
282
283 C-- Initialise preconditioner
284 C For now PC is just the identity. Change to
285 C be LU factorization of d2/dz2 later. Note
286 C check for consistency with S/R CG3D before
287 C assuming zML is lower and zMU is upper!
288 DO bj=myByLo(myThid),myByHi(myThid)
289 DO bi=myBxLo(myThid),myBxHi(myThid)
290 DO k=1,Nr
291 DO j=1,sNy
292 DO i=1,sNx
293 IF ( aC3d(i,j,k,bi,bj) .NE. 0. ) THEN
294 zMC(i,j,k,bi,bj) = aC3d(i,j,k,bi,bj)
295 zML(i,j,k,bi,bj) = aV3d(i,j,k,bi,bj)
296 IF ( k.NE.Nr ) THEN
297 zMU(i,j,k,bi,bj)= aV3d(i,j,k+1,bi,bj)
298 ELSE
299 zMU(i,j,k,bi,bj)= 0.
300 ENDIF
301 CcnhDebugStarts
302 C zMC(i,j,k,bi,bj) = 1.
303 C zMU(i,j,k,bi,bj) = 0.
304 C zML(i,j,k,bi,bj) = 0.
305 CcnhDebugEnds
306 ELSE
307 zMC(i,j,k,bi,bj) = 1. _d 0
308 zMU(i,j,k,bi,bj) = 0.
309 zML(i,j,k,bi,bj) = 0.
310 ENDIF
311 ENDDO
312 ENDDO
313 ENDDO
314 k = 1
315 DO j=1,sNy
316 DO i=1,sNx
317 zMC(i,j,k,bi,bj) = 1. _d 0 / zMC(i,j,k,bi,bj)
318 zMU(i,j,k,bi,bj) = zMU(i,j,k,bi,bj)*zMC(i,j,k,bi,bj)
319 ENDDO
320 ENDDO
321 DO k=2,Nr
322 DO j=1,sNy
323 DO i=1,sNx
324 zMC(i,j,k,bi,bj) = 1. _d 0 /
325 & (zMC(i,j,k,bi,bj)-zML(i,j,k,bi,bj)*zMU(i,j,k-1,bi,bj))
326 zMU(i,j,k,bi,bj) = zMU(i,j,k,bi,bj)*zMC(i,j,k,bi,bj)
327 ENDDO
328 ENDDO
329 ENDDO
330 DO k=1,Nr
331 DO j=1,sNy
332 DO i=1,sNx
333 IF ( aC3d(i,j,k,bi,bj) .EQ. 0. ) THEN
334 zMC(i,j,k,bi,bj) = 1.
335 zML(i,j,k,bi,bj) = 0.
336 zMU(i,j,k,bi,bj) = 0.
337 CcnhDebugStarts
338 C ELSE
339 C zMC(i,j,k,bi,bj) = 1.
340 C zML(i,j,k,bi,bj) = 0.
341 C zMU(i,j,k,bi,bj) = 0.
342 CcnhDEbugEnds
343 ENDIF
344 ENDDO
345 ENDDO
346 ENDDO
347 ENDDO
348 ENDDO
349 C-- Update overlap regions
350 _EXCH_XYZ_RS(zMC, myThid)
351 _EXCH_XYZ_RS(zML, myThid)
352 _EXCH_XYZ_RS(zMU, myThid)
353
354 IF ( debugLevel .GE. debLevB ) THEN
355 CALL WRITE_FLD_XYZ_RS( 'zMC',' ',zMC, 0, myThid )
356 CALL WRITE_FLD_XYZ_RS( 'zML',' ',zML, 0, myThid )
357 CALL WRITE_FLD_XYZ_RS( 'zMU',' ',zMU, 0, myThid )
358 ENDIF
359 CcnhDebugStarts
360 c DO k=1,Nr
361 c DO j=1-OLy,sNy+OLy
362 c DO i=1-OLx,sNx+OLx
363 c phi(i,j,1,1) = zMc(i,j,k,1,1)
364 c ENDDO
365 c ENDDO
366 C CALL PLOT_FIELD_XYRS( phi, 'zMC INI_CG3D.1' , 1, myThid )
367 c ENDDO
368 C CALL PLOT_FIELD_XYRS( zMU, 'zMU INI_CG3D.1' , Nr, 1, myThid )
369 C CALL PLOT_FIELD_XYRS( zML, 'zML INI_CG3D.1' , Nr, 1, myThid )
370 CcnhDebugEnds
371
372 C-- end if (use3Dsolver)
373 ENDIF
374
375 #endif /* ALLOW_NONHYDROSTATIC */
376
377 RETURN
378 END

  ViewVC Help
Powered by ViewVC 1.1.22