1 |
jmc |
1.50 |
C $Header: /u/gcmpack/MITgcm/model/src/ini_cg2d.F,v 1.49 2009/12/08 00:10:26 jmc Exp $ |
2 |
jmc |
1.30 |
C $Name: $ |
3 |
cnh |
1.1 |
|
4 |
edhill |
1.42 |
#include "PACKAGES_CONFIG.h" |
5 |
adcroft |
1.13 |
#include "CPP_OPTIONS.h" |
6 |
cnh |
1.1 |
|
7 |
cnh |
1.37 |
CBOP |
8 |
|
|
C !ROUTINE: INI_CG2D |
9 |
|
|
C !INTERFACE: |
10 |
cnh |
1.1 |
SUBROUTINE INI_CG2D( myThid ) |
11 |
cnh |
1.37 |
|
12 |
|
|
C !DESCRIPTION: \bv |
13 |
|
|
C *==========================================================* |
14 |
jmc |
1.47 |
C | SUBROUTINE INI_CG2D |
15 |
|
|
C | o Initialise 2d conjugate gradient solver operators. |
16 |
cnh |
1.37 |
C *==========================================================* |
17 |
jmc |
1.47 |
C | These arrays are purely a function of the basin geom. |
18 |
|
|
C | We set then here once and them use then repeatedly. |
19 |
cnh |
1.37 |
C *==========================================================* |
20 |
|
|
C \ev |
21 |
|
|
|
22 |
|
|
C !USES: |
23 |
adcroft |
1.23 |
IMPLICIT NONE |
24 |
cnh |
1.1 |
C === Global variables === |
25 |
|
|
#include "SIZE.h" |
26 |
|
|
#include "EEPARAMS.h" |
27 |
|
|
#include "PARAMS.h" |
28 |
|
|
#include "GRID.h" |
29 |
jmc |
1.31 |
#include "SURFACE.h" |
30 |
adcroft |
1.34 |
#include "CG2D.h" |
31 |
cnh |
1.1 |
|
32 |
cnh |
1.37 |
C !INPUT/OUTPUT PARAMETERS: |
33 |
cnh |
1.1 |
C === Routine arguments === |
34 |
|
|
C myThid - Thread no. that called this routine. |
35 |
|
|
INTEGER myThid |
36 |
|
|
|
37 |
cnh |
1.37 |
C !LOCAL VARIABLES: |
38 |
cnh |
1.1 |
C === Local variables === |
39 |
jmc |
1.45 |
C bi,bj :: tile indices |
40 |
jmc |
1.49 |
C i,j,k :: Loop counters |
41 |
jmc |
1.50 |
C faceArea :: Temporary used to hold cell face areas. |
42 |
|
|
C myNorm :: Work variable used in calculating normalisation factor |
43 |
cnh |
1.1 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
44 |
|
|
INTEGER bi, bj |
45 |
jmc |
1.49 |
INTEGER i, j, k, ks |
46 |
cnh |
1.7 |
_RL faceArea |
47 |
cnh |
1.15 |
_RS myNorm |
48 |
jmc |
1.45 |
_RS aC, aCw, aCs |
49 |
cnh |
1.37 |
CEOP |
50 |
jmc |
1.38 |
|
51 |
jmc |
1.45 |
C-- Initialize arrays in common blocs (CG2D.h) ; not really necessary |
52 |
jmc |
1.38 |
C but safer when EXCH do not fill all the overlap regions. |
53 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
54 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
55 |
jmc |
1.49 |
DO j=1-OLy,sNy+OLy |
56 |
|
|
DO i=1-OLx,sNx+OLx |
57 |
|
|
aW2d(i,j,bi,bj) = 0. _d 0 |
58 |
|
|
aS2d(i,j,bi,bj) = 0. _d 0 |
59 |
|
|
aC2d(i,j,bi,bj) = 0. _d 0 |
60 |
|
|
pW(i,j,bi,bj) = 0. _d 0 |
61 |
|
|
pS(i,j,bi,bj) = 0. _d 0 |
62 |
|
|
pC(i,j,bi,bj) = 0. _d 0 |
63 |
jmc |
1.38 |
ENDDO |
64 |
|
|
ENDDO |
65 |
jmc |
1.49 |
DO j=1-1,sNy+1 |
66 |
|
|
DO i=1-1,sNx+1 |
67 |
|
|
cg2d_q(i,j,bi,bj) = 0. _d 0 |
68 |
|
|
cg2d_r(i,j,bi,bj) = 0. _d 0 |
69 |
|
|
cg2d_s(i,j,bi,bj) = 0. _d 0 |
70 |
jmc |
1.47 |
#ifdef ALLOW_CG2D_NSA |
71 |
jmc |
1.49 |
cg2d_z(i,j,bi,bj) = 0. _d 0 |
72 |
jmc |
1.47 |
#endif /* ALLOW_CG2D_NSA */ |
73 |
mlosch |
1.48 |
#ifdef ALLOW_SRCG |
74 |
jmc |
1.49 |
cg2d_y(i,j,bi,bj) = 0. _d 0 |
75 |
|
|
cg2d_v(i,j,bi,bj) = 0. _d 0 |
76 |
mlosch |
1.48 |
#endif /* ALLOW_SRCG */ |
77 |
jmc |
1.38 |
ENDDO |
78 |
|
|
ENDDO |
79 |
|
|
ENDDO |
80 |
|
|
ENDDO |
81 |
jmc |
1.31 |
|
82 |
cnh |
1.1 |
C-- Initialise laplace operator |
83 |
|
|
C aW2d: integral in Z Ax/dX |
84 |
|
|
C aS2d: integral in Z Ay/dY |
85 |
|
|
myNorm = 0. _d 0 |
86 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
87 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
88 |
jmc |
1.49 |
DO j=1,sNy |
89 |
|
|
DO i=1,sNx |
90 |
|
|
aW2d(i,j,bi,bj) = 0. _d 0 |
91 |
|
|
aS2d(i,j,bi,bj) = 0. _d 0 |
92 |
cnh |
1.1 |
ENDDO |
93 |
|
|
ENDDO |
94 |
jmc |
1.49 |
DO k=1,Nr |
95 |
|
|
DO j=1,sNy |
96 |
|
|
DO i=1,sNx |
97 |
jmc |
1.45 |
C deep-model: *deepFacC (faceArea), /deepFacC (recip_dx,y): => no net effect |
98 |
jmc |
1.49 |
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 |
jmc |
1.45 |
& + implicSurfPress*implicDiv2DFlow |
102 |
jmc |
1.49 |
& *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 |
jmc |
1.45 |
& + implicSurfPress*implicDiv2DFlow |
107 |
jmc |
1.49 |
& *faceArea*recip_dyC(i,j,bi,bj) |
108 |
cnh |
1.1 |
ENDDO |
109 |
|
|
ENDDO |
110 |
|
|
ENDDO |
111 |
jmc |
1.49 |
DO j=1,sNy |
112 |
|
|
DO i=1,sNx |
113 |
jmc |
1.50 |
#ifdef ALLOW_OBCS |
114 |
|
|
aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj) |
115 |
|
|
& *maskInC(i,j,bi,bj)*maskInC(i-1,j,bi,bj) |
116 |
|
|
aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj) |
117 |
|
|
& *maskInC(i,j,bi,bj)*maskInC(i,j-1,bi,bj) |
118 |
|
|
#endif /* ALLOW_OBCS */ |
119 |
|
|
myNorm = MAX(ABS(aW2d(i,j,bi,bj)),myNorm) |
120 |
|
|
myNorm = MAX(ABS(aS2d(i,j,bi,bj)),myNorm) |
121 |
cnh |
1.1 |
ENDDO |
122 |
|
|
ENDDO |
123 |
|
|
ENDDO |
124 |
|
|
ENDDO |
125 |
jmc |
1.46 |
_GLOBAL_MAX_RS( myNorm, myThid ) |
126 |
adcroft |
1.25 |
IF ( myNorm .NE. 0. _d 0 ) THEN |
127 |
|
|
myNorm = 1. _d 0/myNorm |
128 |
cnh |
1.1 |
ELSE |
129 |
|
|
myNorm = 1. _d 0 |
130 |
|
|
ENDIF |
131 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
132 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
133 |
jmc |
1.49 |
DO j=1,sNy |
134 |
|
|
DO i=1,sNx |
135 |
|
|
aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj)*myNorm |
136 |
|
|
aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj)*myNorm |
137 |
cnh |
1.1 |
ENDDO |
138 |
|
|
ENDDO |
139 |
|
|
ENDDO |
140 |
|
|
ENDDO |
141 |
jmc |
1.45 |
|
142 |
cnh |
1.1 |
C-- Update overlap regions |
143 |
|
|
CcnhDebugStarts |
144 |
cnh |
1.14 |
C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.1' , 1, myThid ) |
145 |
|
|
C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.1' , 1, myThid ) |
146 |
cnh |
1.1 |
CcnhDebugEnds |
147 |
jmc |
1.47 |
CALL EXCH_UV_XY_RS( aW2d, aS2d, .FALSE., myThid ) |
148 |
cnh |
1.1 |
CcnhDebugStarts |
149 |
adcroft |
1.24 |
C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.2' , 1, myThid ) |
150 |
|
|
C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.2' , 1, myThid ) |
151 |
cnh |
1.1 |
CcnhDebugEnds |
152 |
|
|
|
153 |
jmc |
1.44 |
_BEGIN_MASTER(myThid) |
154 |
|
|
C-- set global parameter in common block: |
155 |
|
|
cg2dNorm = myNorm |
156 |
adcroft |
1.34 |
C-- Define the solver tolerance in the appropriate Unit : |
157 |
jmc |
1.43 |
cg2dNormaliseRHS = cg2dTargetResWunit.LE.0. |
158 |
adcroft |
1.34 |
IF (cg2dNormaliseRHS) THEN |
159 |
|
|
C- when using a normalisation of RHS, tolerance has no unit => no conversion |
160 |
|
|
cg2dTolerance = cg2dTargetResidual |
161 |
|
|
ELSE |
162 |
|
|
C- convert Target-Residual (in W unit) to cg2d-solver residual unit [m^2/s^2] |
163 |
|
|
cg2dTolerance = cg2dNorm * cg2dTargetResWunit |
164 |
jmc |
1.43 |
& * globalArea / deltaTmom |
165 |
adcroft |
1.34 |
ENDIF |
166 |
jmc |
1.44 |
_END_MASTER(myThid) |
167 |
jmc |
1.43 |
|
168 |
|
|
CcnhDebugStarts |
169 |
|
|
_BEGIN_MASTER( myThid ) |
170 |
|
|
WRITE(msgBuf,'(2A,1PE23.16)') 'INI_CG2D: ', |
171 |
|
|
& 'CG2D normalisation factor = ', cg2dNorm |
172 |
|
|
CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1) |
173 |
|
|
IF (.NOT.cg2dNormaliseRHS) THEN |
174 |
|
|
WRITE(msgBuf,'(2A,1PE22.15,A,1PE16.10,A)') 'INI_CG2D: ', |
175 |
|
|
& 'cg2dTolerance =', cg2dTolerance, ' (Area=',globalArea,')' |
176 |
|
|
CALL PRINT_MESSAGE(msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1) |
177 |
|
|
ENDIF |
178 |
|
|
WRITE(msgBuf,*) ' ' |
179 |
|
|
CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1) |
180 |
|
|
_END_MASTER( myThid ) |
181 |
|
|
CcnhDebugEnds |
182 |
jmc |
1.45 |
|
183 |
cnh |
1.1 |
C-- Initialise preconditioner |
184 |
cnh |
1.4 |
C Note. 20th May 1998 |
185 |
|
|
C I made a weird discovery! In the model paper we argue |
186 |
|
|
C for the form of the preconditioner used here ( see |
187 |
|
|
C A Finite-volume, Incompressible Navier-Stokes Model |
188 |
|
|
C ...., Marshall et. al ). The algebra gives a simple |
189 |
|
|
C 0.5 factor for the averaging of ac and aCw to get a |
190 |
|
|
C symmettric pre-conditioner. By using a factor of 0.51 |
191 |
|
|
C i.e. scaling the off-diagonal terms in the |
192 |
|
|
C preconditioner down slightly I managed to get the |
193 |
|
|
C number of iterations for convergence in a test case to |
194 |
|
|
C drop form 192 -> 134! Need to investigate this further! |
195 |
|
|
C For now I have introduced a parameter cg2dpcOffDFac which |
196 |
|
|
C defaults to 0.51 but can be set at runtime. |
197 |
cnh |
1.1 |
DO bj=myByLo(myThid),myByHi(myThid) |
198 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
199 |
jmc |
1.49 |
C- calculate and store solver main diagonal : |
200 |
|
|
DO j=0,sNy+1 |
201 |
|
|
DO i=0,sNx+1 |
202 |
|
|
ks = ksurfC(i,j,bi,bj) |
203 |
|
|
aC2d(i,j,bi,bj) = -( |
204 |
|
|
& aW2d(i,j,bi,bj) + aW2d(i+1,j, bi,bj) |
205 |
|
|
& +aS2d(i,j,bi,bj) + aS2d( i,j+1,bi,bj) |
206 |
|
|
& +freeSurfFac*myNorm*recip_Bo(i,j,bi,bj)*deepFac2F(ks) |
207 |
|
|
& *rA(i,j,bi,bj)/deltaTMom/deltaTfreesurf |
208 |
|
|
& ) |
209 |
|
|
ENDDO |
210 |
|
|
ENDDO |
211 |
|
|
DO j=1,sNy |
212 |
|
|
DO i=1,sNx |
213 |
|
|
aC = aC2d( i, j, bi,bj) |
214 |
|
|
aCs = aC2d( i,j-1,bi,bj) |
215 |
|
|
aCw = aC2d(i-1,j, bi,bj) |
216 |
|
|
IF ( aC .EQ. 0. ) THEN |
217 |
|
|
pC(i,j,bi,bj) = 1. _d 0 |
218 |
|
|
ELSE |
219 |
|
|
pC(i,j,bi,bj) = 1. _d 0 / aC |
220 |
|
|
ENDIF |
221 |
|
|
IF ( aC + aCw .EQ. 0. ) THEN |
222 |
|
|
pW(i,j,bi,bj) = 0. |
223 |
|
|
ELSE |
224 |
|
|
pW(i,j,bi,bj) = |
225 |
|
|
& -aW2d(i,j,bi,bj)/((cg2dpcOffDFac *(aCw+aC))**2 ) |
226 |
|
|
ENDIF |
227 |
|
|
IF ( aC + aCs .EQ. 0. ) THEN |
228 |
|
|
pS(i,j,bi,bj) = 0. |
229 |
|
|
ELSE |
230 |
|
|
pS(i,j,bi,bj) = |
231 |
|
|
& -aS2d(i,j,bi,bj)/((cg2dpcOffDFac *(aCs+aC))**2 ) |
232 |
|
|
ENDIF |
233 |
|
|
C pC(i,j,bi,bj) = 1. |
234 |
|
|
C pW(i,j,bi,bj) = 0. |
235 |
|
|
C pS(i,j,bi,bj) = 0. |
236 |
cnh |
1.1 |
ENDDO |
237 |
|
|
ENDDO |
238 |
|
|
ENDDO |
239 |
|
|
ENDDO |
240 |
|
|
C-- Update overlap regions |
241 |
jmc |
1.47 |
CALL EXCH_XY_RS( pC, myThid ) |
242 |
|
|
CALL EXCH_UV_XY_RS( pW, pS, .FALSE., myThid ) |
243 |
cnh |
1.18 |
CcnhDebugStarts |
244 |
adcroft |
1.24 |
C CALL PLOT_FIELD_XYRS( pC, 'pC INI_CG2D.2' , 1, myThid ) |
245 |
|
|
C CALL PLOT_FIELD_XYRS( pW, 'pW INI_CG2D.2' , 1, myThid ) |
246 |
|
|
C CALL PLOT_FIELD_XYRS( pS, 'pS INI_CG2D.2' , 1, myThid ) |
247 |
cnh |
1.18 |
CcnhDebugEnds |
248 |
cnh |
1.1 |
|
249 |
|
|
RETURN |
250 |
|
|
END |