1 |
adcroft |
1.34 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/ini_cg2d.F,v 1.33.2.4 2001/04/06 13:10:54 jmc Exp $ |
2 |
jmc |
1.30 |
C $Name: $ |
3 |
cnh |
1.1 |
|
4 |
adcroft |
1.13 |
#include "CPP_OPTIONS.h" |
5 |
cnh |
1.1 |
|
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 |
adcroft |
1.23 |
IMPLICIT NONE |
15 |
cnh |
1.1 |
|
16 |
|
|
C === Global variables === |
17 |
|
|
#include "SIZE.h" |
18 |
|
|
#include "EEPARAMS.h" |
19 |
|
|
#include "PARAMS.h" |
20 |
|
|
#include "GRID.h" |
21 |
adcroft |
1.34 |
c#include "DYNVARS.h" |
22 |
jmc |
1.31 |
#include "SURFACE.h" |
23 |
adcroft |
1.34 |
#include "CG2D.h" |
24 |
adcroft |
1.26 |
#ifdef ALLOW_OBCS |
25 |
adcroft |
1.22 |
#include "OBCS.h" |
26 |
adcroft |
1.26 |
#endif |
27 |
cnh |
1.1 |
|
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 |
adcroft |
1.34 |
C myNorm - Work variable used in calculating normalisation factor |
40 |
|
|
C sumArea - Work variable used to compute the total Domain Area |
41 |
cnh |
1.1 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
42 |
|
|
INTEGER bi, bj |
43 |
|
|
INTEGER I, J, K |
44 |
cnh |
1.7 |
_RL faceArea |
45 |
cnh |
1.15 |
_RS myNorm |
46 |
adcroft |
1.34 |
_RL sumArea |
47 |
cnh |
1.4 |
_RL aC, aCw, aCs |
48 |
cnh |
1.1 |
|
49 |
jmc |
1.31 |
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 |
jmc |
1.33 |
C- gBaro = gravity (except for External mode test with reduced gravity) |
69 |
jmc |
1.31 |
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 |
jmc |
1.33 |
Bo_surf(I,J,bi,bj) = gBaro |
74 |
|
|
recip_Bo(I,J,bi,bj) = 1. _d 0 / gBaro |
75 |
jmc |
1.31 |
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 |
cnh |
1.1 |
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 |
cnh |
1.17 |
DO K=1,Nr |
98 |
cnh |
1.1 |
DO J=1,sNy |
99 |
|
|
DO I=1,sNx |
100 |
cnh |
1.20 |
faceArea = _dyG(I,J,bi,bj)*drF(K) |
101 |
|
|
& *_hFacW(I,J,K,bi,bj) |
102 |
cnh |
1.1 |
aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj) + |
103 |
jmc |
1.32 |
& implicSurfPress*implicDiv2DFlow |
104 |
|
|
& *faceArea*recip_dxC(I,J,bi,bj) |
105 |
cnh |
1.20 |
faceArea = _dxG(I,J,bi,bj)*drF(K) |
106 |
|
|
& *_hFacS(I,J,K,bi,bj) |
107 |
cnh |
1.1 |
aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj) + |
108 |
jmc |
1.32 |
& implicSurfPress*implicDiv2DFlow |
109 |
|
|
& *faceArea*recip_dyC(I,J,bi,bj) |
110 |
cnh |
1.1 |
ENDDO |
111 |
|
|
ENDDO |
112 |
|
|
ENDDO |
113 |
adcroft |
1.26 |
#ifdef ALLOW_OBCS |
114 |
adcroft |
1.28 |
IF (useOBCS) THEN |
115 |
adcroft |
1.22 |
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 |
adcroft |
1.26 |
#endif |
129 |
cnh |
1.1 |
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 |
adcroft |
1.25 |
_GLOBAL_MAX_R4( myNorm, myThid ) |
138 |
|
|
IF ( myNorm .NE. 0. _d 0 ) THEN |
139 |
|
|
myNorm = 1. _d 0/myNorm |
140 |
cnh |
1.1 |
ELSE |
141 |
|
|
myNorm = 1. _d 0 |
142 |
|
|
ENDIF |
143 |
cnh |
1.5 |
cg2dNorm = myNorm |
144 |
cnh |
1.1 |
_BEGIN_MASTER( myThid ) |
145 |
|
|
CcnhDebugStarts |
146 |
cnh |
1.20 |
WRITE(msgBuf,'(A,E40.25)') '// CG2D normalisation factor = ', |
147 |
|
|
& cg2dNorm |
148 |
cnh |
1.3 |
CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1) |
149 |
|
|
WRITE(msgBuf,*) ' ' |
150 |
cnh |
1.1 |
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 |
cnh |
1.14 |
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 |
cnh |
1.1 |
CcnhDebugEnds |
169 |
adcroft |
1.34 |
c _EXCH_XY_R4(aW2d, myThid) |
170 |
|
|
c _EXCH_XY_R4(aS2d, myThid) |
171 |
|
|
CALL EXCH_UV_XY_RS(aW2d,aS2d,.FALSE.,myThid) |
172 |
cnh |
1.1 |
CcnhDebugStarts |
173 |
adcroft |
1.24 |
C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.2' , 1, myThid ) |
174 |
|
|
C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.2' , 1, myThid ) |
175 |
cnh |
1.1 |
CcnhDebugEnds |
176 |
|
|
|
177 |
adcroft |
1.34 |
C-- Define the solver tolerance in the appropriate Unit : |
178 |
|
|
cg2dNormaliseRHS = cg2dTargetResWunit.LE.0 |
179 |
|
|
IF (cg2dNormaliseRHS) THEN |
180 |
|
|
C- when using a normalisation of RHS, tolerance has no unit => no conversion |
181 |
|
|
cg2dTolerance = cg2dTargetResidual |
182 |
|
|
ELSE |
183 |
|
|
C- compute the total Area of the domain : |
184 |
|
|
sumArea = 0. |
185 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
186 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
187 |
|
|
DO j=1,sNy |
188 |
|
|
DO i=1,sNx |
189 |
|
|
IF (Ro_surf(i,j,bi,bj).GT.R_low(i,j,bi,bj)) THEN |
190 |
|
|
sumArea = sumArea + rA(i,j,bi,bj) |
191 |
|
|
ENDIF |
192 |
|
|
ENDDO |
193 |
|
|
ENDDO |
194 |
|
|
ENDDO |
195 |
|
|
ENDDO |
196 |
|
|
c WRITE(0,*) ' mythid, sumArea = ', mythid, sumArea |
197 |
|
|
_GLOBAL_SUM_R8( sumArea, myThid ) |
198 |
|
|
C- convert Target-Residual (in W unit) to cg2d-solver residual unit [m^2/s^2] |
199 |
|
|
cg2dTolerance = cg2dNorm * cg2dTargetResWunit |
200 |
|
|
& * sumArea / deltaTMom |
201 |
|
|
WRITE(0,'(2A,1P2E22.14)') ' ini_cg2d: ', |
202 |
|
|
& 'sumArea,cg2dTolerance =', sumArea,cg2dTolerance |
203 |
|
|
ENDIF |
204 |
|
|
|
205 |
cnh |
1.1 |
C-- Initialise preconditioner |
206 |
cnh |
1.4 |
C Note. 20th May 1998 |
207 |
|
|
C I made a weird discovery! In the model paper we argue |
208 |
|
|
C for the form of the preconditioner used here ( see |
209 |
|
|
C A Finite-volume, Incompressible Navier-Stokes Model |
210 |
|
|
C ...., Marshall et. al ). The algebra gives a simple |
211 |
|
|
C 0.5 factor for the averaging of ac and aCw to get a |
212 |
|
|
C symmettric pre-conditioner. By using a factor of 0.51 |
213 |
|
|
C i.e. scaling the off-diagonal terms in the |
214 |
|
|
C preconditioner down slightly I managed to get the |
215 |
|
|
C number of iterations for convergence in a test case to |
216 |
|
|
C drop form 192 -> 134! Need to investigate this further! |
217 |
|
|
C For now I have introduced a parameter cg2dpcOffDFac which |
218 |
|
|
C defaults to 0.51 but can be set at runtime. |
219 |
cnh |
1.1 |
DO bj=myByLo(myThid),myByHi(myThid) |
220 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
221 |
|
|
DO J=1,sNy |
222 |
|
|
DO I=1,sNx |
223 |
|
|
pC(I,J,bi,bj) = 1. _d 0 |
224 |
cnh |
1.4 |
aC = -( |
225 |
|
|
& aW2d(I,J,bi,bj) + aW2d(I+1,J ,bi,bj) |
226 |
adcroft |
1.24 |
& +aS2d(I,J,bi,bj) + aS2d(I ,J+1,bi,bj) |
227 |
jmc |
1.32 |
& +freeSurfFac*myNorm*recip_Bo(I,J,bi,bj)* |
228 |
cnh |
1.17 |
& rA(I,J,bi,bj)/deltaTMom/deltaTMom |
229 |
cnh |
1.4 |
& ) |
230 |
|
|
aCs = -( |
231 |
|
|
& aW2d(I,J-1,bi,bj) + aW2d(I+1,J-1,bi,bj) |
232 |
|
|
& +aS2d(I,J-1,bi,bj) + aS2d(I ,J ,bi,bj) |
233 |
jmc |
1.32 |
& +freeSurfFac*myNorm*recip_Bo(I,J-1,bi,bj)* |
234 |
cnh |
1.17 |
& rA(I,J-1,bi,bj)/deltaTMom/deltaTMom |
235 |
cnh |
1.4 |
& ) |
236 |
|
|
aCw = -( |
237 |
|
|
& aW2d(I-1,J,bi,bj) + aW2d(I ,J ,bi,bj) |
238 |
|
|
& +aS2d(I-1,J,bi,bj) + aS2d(I-1,J+1,bi,bj) |
239 |
jmc |
1.32 |
& +freeSurfFac*myNorm*recip_Bo(I-1,J,bi,bj)* |
240 |
cnh |
1.17 |
& rA(I-1,J,bi,bj)/deltaTMom/deltaTMom |
241 |
cnh |
1.4 |
& ) |
242 |
|
|
IF ( aC .EQ. 0. ) THEN |
243 |
adcroft |
1.27 |
pC(I,J,bi,bj) = 1. _d 0 |
244 |
cnh |
1.4 |
ELSE |
245 |
|
|
pC(I,J,bi,bj) = 1. _d 0 / aC |
246 |
|
|
ENDIF |
247 |
|
|
IF ( aC + aCw .EQ. 0. ) THEN |
248 |
|
|
pW(I,J,bi,bj) = 0. |
249 |
|
|
ELSE |
250 |
|
|
pW(I,J,bi,bj) = |
251 |
|
|
& -aW2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCw+aC))**2 ) |
252 |
|
|
ENDIF |
253 |
|
|
IF ( aC + aCs .EQ. 0. ) THEN |
254 |
|
|
pS(I,J,bi,bj) = 0. |
255 |
|
|
ELSE |
256 |
|
|
pS(I,J,bi,bj) = |
257 |
|
|
& -aS2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCs+aC))**2 ) |
258 |
|
|
ENDIF |
259 |
cnh |
1.6 |
C pC(I,J,bi,bj) = 1. |
260 |
|
|
C pW(I,J,bi,bj) = 0. |
261 |
|
|
C pS(I,J,bi,bj) = 0. |
262 |
cnh |
1.1 |
ENDDO |
263 |
|
|
ENDDO |
264 |
|
|
ENDDO |
265 |
|
|
ENDDO |
266 |
|
|
C-- Update overlap regions |
267 |
|
|
_EXCH_XY_R4(pC, myThid) |
268 |
adcroft |
1.34 |
c _EXCH_XY_R4(pW, myThid) |
269 |
|
|
c _EXCH_XY_R4(pS, myThid) |
270 |
|
|
CALL EXCH_UV_XY_RS(pW,pS,.FALSE.,myThid) |
271 |
cnh |
1.18 |
CcnhDebugStarts |
272 |
adcroft |
1.24 |
C CALL PLOT_FIELD_XYRS( pC, 'pC INI_CG2D.2' , 1, myThid ) |
273 |
|
|
C CALL PLOT_FIELD_XYRS( pW, 'pW INI_CG2D.2' , 1, myThid ) |
274 |
|
|
C CALL PLOT_FIELD_XYRS( pS, 'pS INI_CG2D.2' , 1, myThid ) |
275 |
cnh |
1.18 |
CcnhDebugEnds |
276 |
cnh |
1.1 |
|
277 |
|
|
RETURN |
278 |
|
|
END |