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

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

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


Revision 1.1 - (hide annotations) (download)
Mon Mar 22 15:54:04 1999 UTC (25 years, 2 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint20, checkpoint21
Modifications for non-hydrostatic ability + updates for open-boundaries.

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

  ViewVC Help
Powered by ViewVC 1.1.22