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