/[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.47 - (show annotations) (download)
Sat Jul 11 22:00:40 2009 UTC (14 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61x, checkpoint61y
Changes since 1.46: +12 -13 lines
use simple EXCH (overlap 1 and no Corner Exch): this reduces number of
 EXCH calls by 2 if using exch2).

1 C $Header: /u/gcmpack/MITgcm/model/src/ini_cg2d.F,v 1.46 2009/04/28 18:01:14 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 #ifdef ALLOW_OBCS
32 #include "OBCS.h"
33 #endif
34
35 C !INPUT/OUTPUT PARAMETERS:
36 C === Routine arguments ===
37 C myThid - Thread no. that called this routine.
38 INTEGER myThid
39
40 C !LOCAL VARIABLES:
41 C === Local variables ===
42 C bi,bj :: tile indices
43 C I,J,K :: Loop counters
44 C faceArea - Temporary used to hold cell face areas.
45 C myNorm - Work variable used in calculating normalisation factor
46 C sumArea - Work variable used to compute the total Domain Area
47 CHARACTER*(MAX_LEN_MBUF) msgBuf
48 INTEGER bi, bj
49 INTEGER I, J, K, ks
50 _RL faceArea
51 _RS myNorm
52 _RS aC, aCw, aCs
53 CEOP
54
55 C-- Initialize arrays in common blocs (CG2D.h) ; not really necessary
56 C but safer when EXCH do not fill all the overlap regions.
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 aW2d(I,J,bi,bj) = 0. _d 0
62 aS2d(I,J,bi,bj) = 0. _d 0
63 aC2d(I,J,bi,bj) = 0. _d 0
64 pW(I,J,bi,bj) = 0. _d 0
65 pS(I,J,bi,bj) = 0. _d 0
66 pC(I,J,bi,bj) = 0. _d 0
67 ENDDO
68 ENDDO
69 DO J=1-1,sNy+1
70 DO I=1-1,sNx+1
71 cg2d_q(I,J,bi,bj) = 0. _d 0
72 cg2d_r(I,J,bi,bj) = 0. _d 0
73 cg2d_s(I,J,bi,bj) = 0. _d 0
74 #ifdef ALLOW_CG2D_NSA
75 cg2d_z(I,J,bi,bj) = 0. _d 0
76 #endif /* ALLOW_CG2D_NSA */
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 #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_RS( 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 DO bj=myByLo(myThid),myByHi(myThid)
142 DO bi=myBxLo(myThid),myBxHi(myThid)
143 DO J=1,sNy
144 DO I=1,sNx
145 aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj)*myNorm
146 aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj)*myNorm
147 ENDDO
148 ENDDO
149 ENDDO
150 ENDDO
151
152 C-- Update overlap regions
153 CcnhDebugStarts
154 C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.1' , 1, myThid )
155 C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.1' , 1, myThid )
156 CcnhDebugEnds
157 CALL EXCH_UV_XY_RS( aW2d, aS2d, .FALSE., myThid )
158 CcnhDebugStarts
159 C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.2' , 1, myThid )
160 C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.2' , 1, myThid )
161 CcnhDebugEnds
162
163 _BEGIN_MASTER(myThid)
164 C-- set global parameter in common block:
165 cg2dNorm = myNorm
166 C-- Define the solver tolerance in the appropriate Unit :
167 cg2dNormaliseRHS = cg2dTargetResWunit.LE.0.
168 IF (cg2dNormaliseRHS) THEN
169 C- when using a normalisation of RHS, tolerance has no unit => no conversion
170 cg2dTolerance = cg2dTargetResidual
171 ELSE
172 C- convert Target-Residual (in W unit) to cg2d-solver residual unit [m^2/s^2]
173 cg2dTolerance = cg2dNorm * cg2dTargetResWunit
174 & * globalArea / deltaTmom
175 ENDIF
176 _END_MASTER(myThid)
177
178 CcnhDebugStarts
179 _BEGIN_MASTER( myThid )
180 WRITE(msgBuf,'(2A,1PE23.16)') 'INI_CG2D: ',
181 & 'CG2D normalisation factor = ', cg2dNorm
182 CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
183 IF (.NOT.cg2dNormaliseRHS) THEN
184 WRITE(msgBuf,'(2A,1PE22.15,A,1PE16.10,A)') 'INI_CG2D: ',
185 & 'cg2dTolerance =', cg2dTolerance, ' (Area=',globalArea,')'
186 CALL PRINT_MESSAGE(msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
187 ENDIF
188 WRITE(msgBuf,*) ' '
189 CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
190 _END_MASTER( myThid )
191 CcnhDebugEnds
192
193 C-- Initialise preconditioner
194 C Note. 20th May 1998
195 C I made a weird discovery! In the model paper we argue
196 C for the form of the preconditioner used here ( see
197 C A Finite-volume, Incompressible Navier-Stokes Model
198 C ...., Marshall et. al ). The algebra gives a simple
199 C 0.5 factor for the averaging of ac and aCw to get a
200 C symmettric pre-conditioner. By using a factor of 0.51
201 C i.e. scaling the off-diagonal terms in the
202 C preconditioner down slightly I managed to get the
203 C number of iterations for convergence in a test case to
204 C drop form 192 -> 134! Need to investigate this further!
205 C For now I have introduced a parameter cg2dpcOffDFac which
206 C defaults to 0.51 but can be set at runtime.
207 DO bj=myByLo(myThid),myByHi(myThid)
208 DO bi=myBxLo(myThid),myBxHi(myThid)
209 DO J=1,sNy
210 DO I=1,sNx
211 ks = ksurfC(I,J,bi,bj)
212 pC(I,J,bi,bj) = 1. _d 0
213 aC = -(
214 & aW2d(I,J,bi,bj) + aW2d(I+1,J ,bi,bj)
215 & +aS2d(I,J,bi,bj) + aS2d(I ,J+1,bi,bj)
216 & +freeSurfFac*myNorm*recip_Bo(I,J,bi,bj)*deepFac2F(ks)
217 & *rA(I,J,bi,bj)/deltaTMom/deltaTfreesurf
218 & )
219 aCs = -(
220 & aW2d(I,J-1,bi,bj) + aW2d(I+1,J-1,bi,bj)
221 & +aS2d(I,J-1,bi,bj) + aS2d(I ,J ,bi,bj)
222 & +freeSurfFac*myNorm*recip_Bo(I,J-1,bi,bj)*deepFac2F(ks)
223 & *rA(I,J-1,bi,bj)/deltaTMom/deltaTfreesurf
224 & )
225 aCw = -(
226 & aW2d(I-1,J,bi,bj) + aW2d(I ,J ,bi,bj)
227 & +aS2d(I-1,J,bi,bj) + aS2d(I-1,J+1,bi,bj)
228 & +freeSurfFac*myNorm*recip_Bo(I-1,J,bi,bj)*deepFac2F(ks)
229 & *rA(I-1,J,bi,bj)/deltaTMom/deltaTfreesurf
230 & )
231 IF ( aC .EQ. 0. ) THEN
232 pC(I,J,bi,bj) = 1. _d 0
233 ELSE
234 pC(I,J,bi,bj) = 1. _d 0 / aC
235 ENDIF
236 IF ( aC + aCw .EQ. 0. ) THEN
237 pW(I,J,bi,bj) = 0.
238 ELSE
239 pW(I,J,bi,bj) =
240 & -aW2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCw+aC))**2 )
241 ENDIF
242 IF ( aC + aCs .EQ. 0. ) THEN
243 pS(I,J,bi,bj) = 0.
244 ELSE
245 pS(I,J,bi,bj) =
246 & -aS2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCs+aC))**2 )
247 ENDIF
248 C- store solver main diagonal :
249 aC2d(I,J,bi,bj) = aC
250 C pC(I,J,bi,bj) = 1.
251 C pW(I,J,bi,bj) = 0.
252 C pS(I,J,bi,bj) = 0.
253 ENDDO
254 ENDDO
255 ENDDO
256 ENDDO
257 C-- Update overlap regions
258 CALL EXCH_XY_RS( pC, myThid )
259 CALL EXCH_UV_XY_RS( pW, pS, .FALSE., myThid )
260 CcnhDebugStarts
261 C CALL PLOT_FIELD_XYRS( pC, 'pC INI_CG2D.2' , 1, myThid )
262 C CALL PLOT_FIELD_XYRS( pW, 'pW INI_CG2D.2' , 1, myThid )
263 C CALL PLOT_FIELD_XYRS( pS, 'pS INI_CG2D.2' , 1, myThid )
264 CcnhDebugEnds
265
266 RETURN
267 END

  ViewVC Help
Powered by ViewVC 1.1.22