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

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

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


Revision 1.33.2.2 - (show annotations) (download)
Fri Mar 23 20:46:07 2001 UTC (23 years, 2 months ago) by adcroft
Branch: pre38
Changes since 1.33.2.1: +2 -2 lines
Ooops. Including wrong OPTIONS file.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/ini_cg2d.F,v 1.33.2.1 2001/03/23 20:44:53 adcroft Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 SUBROUTINE INI_CG2D( myThid )
7 C /==========================================================\
8 C | SUBROUTINE INI_CG2D |
9 C | o Initialise 2d 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 IMPLICIT NONE
15
16 C === Global variables ===
17 #include "SIZE.h"
18 #include "EEPARAMS.h"
19 #include "PARAMS.h"
20 #include "GRID.h"
21 c#include "DYNVARS.h"
22 #include "SURFACE.h"
23 #include "CG2D.h"
24 #ifdef ALLOW_OBCS
25 #include "OBCS.h"
26 #endif
27
28 C === Routine arguments ===
29 C myThid - Thread no. that called this routine.
30 INTEGER myThid
31
32 C === Local variables ===
33 C xG, yG - Global coordinate location.
34 C zG
35 C iG, jG - Global coordinate index
36 C bi,bj - Loop counters
37 C faceArea - Temporary used to hold cell face areas.
38 C I,J,K
39 C myNorm - Work variable used in clculating normalisation factor
40 CHARACTER*(MAX_LEN_MBUF) msgBuf
41 INTEGER bi, bj
42 INTEGER I, J, K
43 _RL faceArea
44 _RS myNorm
45 _RL aC, aCw, aCs
46
47 C-- Initialise -Boyancy at surface level : Bo_surf
48 C Bo_surf is defined as d/dr(Phi_surf) and set to g/rtoz (linear free surface)
49 C with rtoz = conversion factor from r-unit to z-unit (=horiVertRatio)
50 C an acurate formulation will include P_surf and T,S_surf effects on rho_surf:
51 C z-ocean (rtoz=1) : Bo_surf = - Boyancy = gravity * rho_surf/rho_0
52 C p-atmos (rtoz=rho_c*g) : Bo_surf = (1/rho)_surf
53 C--
54 IF ( buoyancyRelation .EQ. 'ATMOSPHERIC' ) THEN
55 DO bj=myByLo(myThid),myByHi(myThid)
56 DO bi=myBxLo(myThid),myBxHi(myThid)
57 DO J=1-Oly,sNy+Oly
58 DO I=1-Olx,sNx+Olx
59 Bo_surf(I,J,bi,bj) = recip_rhoConst
60 recip_Bo(I,J,bi,bj) = rhoConst
61 ENDDO
62 ENDDO
63 ENDDO
64 ENDDO
65 ELSE
66 C- gBaro = gravity (except for External mode test with reduced gravity)
67 DO bj=myByLo(myThid),myByHi(myThid)
68 DO bi=myBxLo(myThid),myBxHi(myThid)
69 DO J=1-Oly,sNy+Oly
70 DO I=1-Olx,sNx+Olx
71 Bo_surf(I,J,bi,bj) = gBaro
72 recip_Bo(I,J,bi,bj) = 1. _d 0 / gBaro
73 ENDDO
74 ENDDO
75 ENDDO
76 ENDDO
77 ENDIF
78
79 C-- Update overlap regions
80 _EXCH_XY_R8(Bo_surf, myThid)
81 _EXCH_XY_R8(recip_Bo, myThid)
82
83 C-- Initialise laplace operator
84 C aW2d: integral in Z Ax/dX
85 C aS2d: integral in Z Ay/dY
86 myNorm = 0. _d 0
87 DO bj=myByLo(myThid),myByHi(myThid)
88 DO bi=myBxLo(myThid),myBxHi(myThid)
89 DO J=1,sNy
90 DO I=1,sNx
91 aW2d(I,J,bi,bj) = 0. _d 0
92 aS2d(I,J,bi,bj) = 0. _d 0
93 ENDDO
94 ENDDO
95 DO K=1,Nr
96 DO J=1,sNy
97 DO I=1,sNx
98 faceArea = _dyG(I,J,bi,bj)*drF(K)
99 & *_hFacW(I,J,K,bi,bj)
100 aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj) +
101 & implicSurfPress*implicDiv2DFlow
102 & *faceArea*recip_dxC(I,J,bi,bj)
103 faceArea = _dxG(I,J,bi,bj)*drF(K)
104 & *_hFacS(I,J,K,bi,bj)
105 aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj) +
106 & implicSurfPress*implicDiv2DFlow
107 & *faceArea*recip_dyC(I,J,bi,bj)
108 ENDDO
109 ENDDO
110 ENDDO
111 #ifdef ALLOW_OBCS
112 IF (useOBCS) THEN
113 DO I=1,sNx
114 IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj),bi,bj)=0.
115 IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj)+1,bi,bj)=0.
116 IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj)+1,bi,bj)=0.
117 IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj),bi,bj)=0.
118 ENDDO
119 DO J=1,sNy
120 IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj),J,bi,bj)=0.
121 IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj)+1,J,bi,bj)=0.
122 IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj)+1,J,bi,bj)=0.
123 IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj),J,bi,bj)=0.
124 ENDDO
125 ENDIF
126 #endif
127 DO J=1,sNy
128 DO I=1,sNx
129 myNorm = MAX(ABS(aW2d(I,J,bi,bj)),myNorm)
130 myNorm = MAX(ABS(aS2d(I,J,bi,bj)),myNorm)
131 ENDDO
132 ENDDO
133 ENDDO
134 ENDDO
135 _GLOBAL_MAX_R4( myNorm, myThid )
136 IF ( myNorm .NE. 0. _d 0 ) THEN
137 myNorm = 1. _d 0/myNorm
138 ELSE
139 myNorm = 1. _d 0
140 ENDIF
141 cg2dNorm = myNorm
142 _BEGIN_MASTER( myThid )
143 CcnhDebugStarts
144 WRITE(msgBuf,'(A,E40.25)') '// CG2D normalisation factor = ',
145 & cg2dNorm
146 CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
147 WRITE(msgBuf,*) ' '
148 CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
149 CcnhDebugEnds
150 _END_MASTER( myThid )
151 DO bj=myByLo(myThid),myByHi(myThid)
152 DO bi=myBxLo(myThid),myBxHi(myThid)
153 DO J=1,sNy
154 DO I=1,sNx
155 aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj)*myNorm
156 aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj)*myNorm
157 ENDDO
158 ENDDO
159 ENDDO
160 ENDDO
161
162 C-- Update overlap regions
163 CcnhDebugStarts
164 C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.1' , 1, myThid )
165 C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.1' , 1, myThid )
166 CcnhDebugEnds
167 _EXCH_XY_R4(aW2d, myThid)
168 _EXCH_XY_R4(aS2d, myThid)
169 CcnhDebugStarts
170 C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.2' , 1, myThid )
171 C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.2' , 1, myThid )
172 CcnhDebugEnds
173
174 C-- Initialise preconditioner
175 C Note. 20th May 1998
176 C I made a weird discovery! In the model paper we argue
177 C for the form of the preconditioner used here ( see
178 C A Finite-volume, Incompressible Navier-Stokes Model
179 C ...., Marshall et. al ). The algebra gives a simple
180 C 0.5 factor for the averaging of ac and aCw to get a
181 C symmettric pre-conditioner. By using a factor of 0.51
182 C i.e. scaling the off-diagonal terms in the
183 C preconditioner down slightly I managed to get the
184 C number of iterations for convergence in a test case to
185 C drop form 192 -> 134! Need to investigate this further!
186 C For now I have introduced a parameter cg2dpcOffDFac which
187 C defaults to 0.51 but can be set at runtime.
188 DO bj=myByLo(myThid),myByHi(myThid)
189 DO bi=myBxLo(myThid),myBxHi(myThid)
190 DO J=1,sNy
191 DO I=1,sNx
192 pC(I,J,bi,bj) = 1. _d 0
193 aC = -(
194 & aW2d(I,J,bi,bj) + aW2d(I+1,J ,bi,bj)
195 & +aS2d(I,J,bi,bj) + aS2d(I ,J+1,bi,bj)
196 & +freeSurfFac*myNorm*recip_Bo(I,J,bi,bj)*
197 & rA(I,J,bi,bj)/deltaTMom/deltaTMom
198 & )
199 aCs = -(
200 & aW2d(I,J-1,bi,bj) + aW2d(I+1,J-1,bi,bj)
201 & +aS2d(I,J-1,bi,bj) + aS2d(I ,J ,bi,bj)
202 & +freeSurfFac*myNorm*recip_Bo(I,J-1,bi,bj)*
203 & rA(I,J-1,bi,bj)/deltaTMom/deltaTMom
204 & )
205 aCw = -(
206 & aW2d(I-1,J,bi,bj) + aW2d(I ,J ,bi,bj)
207 & +aS2d(I-1,J,bi,bj) + aS2d(I-1,J+1,bi,bj)
208 & +freeSurfFac*myNorm*recip_Bo(I-1,J,bi,bj)*
209 & rA(I-1,J,bi,bj)/deltaTMom/deltaTMom
210 & )
211 IF ( aC .EQ. 0. ) THEN
212 pC(I,J,bi,bj) = 1. _d 0
213 ELSE
214 pC(I,J,bi,bj) = 1. _d 0 / aC
215 ENDIF
216 IF ( aC + aCw .EQ. 0. ) THEN
217 pW(I,J,bi,bj) = 0.
218 ELSE
219 pW(I,J,bi,bj) =
220 & -aW2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCw+aC))**2 )
221 ENDIF
222 IF ( aC + aCs .EQ. 0. ) THEN
223 pS(I,J,bi,bj) = 0.
224 ELSE
225 pS(I,J,bi,bj) =
226 & -aS2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCs+aC))**2 )
227 ENDIF
228 C pC(I,J,bi,bj) = 1.
229 C pW(I,J,bi,bj) = 0.
230 C pS(I,J,bi,bj) = 0.
231 ENDDO
232 ENDDO
233 ENDDO
234 ENDDO
235 C-- Update overlap regions
236 _EXCH_XY_R4(pC, myThid)
237 _EXCH_XY_R4(pW, myThid)
238 _EXCH_XY_R4(pS, myThid)
239 CcnhDebugStarts
240 C CALL PLOT_FIELD_XYRS( pC, 'pC INI_CG2D.2' , 1, myThid )
241 C CALL PLOT_FIELD_XYRS( pW, 'pW INI_CG2D.2' , 1, myThid )
242 C CALL PLOT_FIELD_XYRS( pS, 'pS INI_CG2D.2' , 1, myThid )
243 CcnhDebugEnds
244
245 RETURN
246 END

  ViewVC Help
Powered by ViewVC 1.1.22