/[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.50 - (show annotations) (download)
Wed May 18 01:15:18 2011 UTC (13 years ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62z, checkpoint62y, checkpoint63g, checkpoint63, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63a, checkpoint63b, checkpoint63c
Changes since 1.49: +11 -25 lines
with OBCS: use maskInC to cancel cg2d-matrix coeff at OB
 and also outside OB interior region (<- was missing in previous code);
This can change cg2d Norm value and, as a consequence, can affect the results.

1 C $Header: /u/gcmpack/MITgcm/model/src/ini_cg2d.F,v 1.49 2009/12/08 00:10:26 jmc Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 CBOP
8 C !ROUTINE: INI_CG2D
9 C !INTERFACE:
10 SUBROUTINE INI_CG2D( myThid )
11
12 C !DESCRIPTION: \bv
13 C *==========================================================*
14 C | SUBROUTINE INI_CG2D
15 C | o Initialise 2d conjugate gradient solver operators.
16 C *==========================================================*
17 C | These arrays are purely a function of the basin geom.
18 C | We set then here once and them use then repeatedly.
19 C *==========================================================*
20 C \ev
21
22 C !USES:
23 IMPLICIT NONE
24 C === Global variables ===
25 #include "SIZE.h"
26 #include "EEPARAMS.h"
27 #include "PARAMS.h"
28 #include "GRID.h"
29 #include "SURFACE.h"
30 #include "CG2D.h"
31
32 C !INPUT/OUTPUT PARAMETERS:
33 C === Routine arguments ===
34 C myThid - Thread no. that called this routine.
35 INTEGER myThid
36
37 C !LOCAL VARIABLES:
38 C === Local variables ===
39 C bi,bj :: tile indices
40 C i,j,k :: Loop counters
41 C faceArea :: Temporary used to hold cell face areas.
42 C myNorm :: Work variable used in calculating normalisation factor
43 CHARACTER*(MAX_LEN_MBUF) msgBuf
44 INTEGER bi, bj
45 INTEGER i, j, k, ks
46 _RL faceArea
47 _RS myNorm
48 _RS aC, aCw, aCs
49 CEOP
50
51 C-- Initialize arrays in common blocs (CG2D.h) ; not really necessary
52 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 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 ENDDO
64 ENDDO
65 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 #ifdef ALLOW_CG2D_NSA
71 cg2d_z(i,j,bi,bj) = 0. _d 0
72 #endif /* ALLOW_CG2D_NSA */
73 #ifdef ALLOW_SRCG
74 cg2d_y(i,j,bi,bj) = 0. _d 0
75 cg2d_v(i,j,bi,bj) = 0. _d 0
76 #endif /* ALLOW_SRCG */
77 ENDDO
78 ENDDO
79 ENDDO
80 ENDDO
81
82 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 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 ENDDO
93 ENDDO
94 DO k=1,Nr
95 DO j=1,sNy
96 DO i=1,sNx
97 C deep-model: *deepFacC (faceArea), /deepFacC (recip_dx,y): => no net effect
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 DO j=1,sNy
112 DO i=1,sNx
113 #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 ENDDO
122 ENDDO
123 ENDDO
124 ENDDO
125 _GLOBAL_MAX_RS( myNorm, myThid )
126 IF ( myNorm .NE. 0. _d 0 ) THEN
127 myNorm = 1. _d 0/myNorm
128 ELSE
129 myNorm = 1. _d 0
130 ENDIF
131 DO bj=myByLo(myThid),myByHi(myThid)
132 DO bi=myBxLo(myThid),myBxHi(myThid)
133 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 ENDDO
138 ENDDO
139 ENDDO
140 ENDDO
141
142 C-- Update overlap regions
143 CcnhDebugStarts
144 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 CcnhDebugEnds
147 CALL EXCH_UV_XY_RS( aW2d, aS2d, .FALSE., myThid )
148 CcnhDebugStarts
149 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 CcnhDebugEnds
152
153 _BEGIN_MASTER(myThid)
154 C-- set global parameter in common block:
155 cg2dNorm = myNorm
156 C-- Define the solver tolerance in the appropriate Unit :
157 cg2dNormaliseRHS = cg2dTargetResWunit.LE.0.
158 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 & * globalArea / deltaTmom
165 ENDIF
166 _END_MASTER(myThid)
167
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
183 C-- Initialise preconditioner
184 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 DO bj=myByLo(myThid),myByHi(myThid)
198 DO bi=myBxLo(myThid),myBxHi(myThid)
199 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 ENDDO
237 ENDDO
238 ENDDO
239 ENDDO
240 C-- Update overlap regions
241 CALL EXCH_XY_RS( pC, myThid )
242 CALL EXCH_UV_XY_RS( pW, pS, .FALSE., myThid )
243 CcnhDebugStarts
244 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 CcnhDebugEnds
248
249 RETURN
250 END

  ViewVC Help
Powered by ViewVC 1.1.22