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