/[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.51 - (show annotations) (download)
Fri Jun 15 20:13:16 2012 UTC (11 years, 10 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint64, checkpoint65, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63o, checkpoint65o, checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, HEAD
Changes since 1.50: +6 -1 lines
Initialise some scalars

1 C $Header: /u/gcmpack/MITgcm/model/src/ini_cg2d.F,v 1.50 2011/05/18 01:15:18 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-- Init. scalars
83 cg2dNorm = 0. _d 0
84 cg2dNormaliseRHS = .FALSE.
85 cg2dtolerance = 0. _d 0
86
87 C-- Initialise laplace operator
88 C aW2d: integral in Z Ax/dX
89 C aS2d: integral in Z Ay/dY
90 myNorm = 0. _d 0
91 DO bj=myByLo(myThid),myByHi(myThid)
92 DO bi=myBxLo(myThid),myBxHi(myThid)
93 DO j=1,sNy
94 DO i=1,sNx
95 aW2d(i,j,bi,bj) = 0. _d 0
96 aS2d(i,j,bi,bj) = 0. _d 0
97 ENDDO
98 ENDDO
99 DO k=1,Nr
100 DO j=1,sNy
101 DO i=1,sNx
102 C deep-model: *deepFacC (faceArea), /deepFacC (recip_dx,y): => no net effect
103 faceArea = _dyG(i,j,bi,bj)*drF(k)
104 & *_hFacW(i,j,k,bi,bj)
105 aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj)
106 & + implicSurfPress*implicDiv2DFlow
107 & *faceArea*recip_dxC(i,j,bi,bj)
108 faceArea = _dxG(i,j,bi,bj)*drF(k)
109 & *_hFacS(i,j,k,bi,bj)
110 aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj)
111 & + implicSurfPress*implicDiv2DFlow
112 & *faceArea*recip_dyC(i,j,bi,bj)
113 ENDDO
114 ENDDO
115 ENDDO
116 DO j=1,sNy
117 DO i=1,sNx
118 #ifdef ALLOW_OBCS
119 aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj)
120 & *maskInC(i,j,bi,bj)*maskInC(i-1,j,bi,bj)
121 aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj)
122 & *maskInC(i,j,bi,bj)*maskInC(i,j-1,bi,bj)
123 #endif /* ALLOW_OBCS */
124 myNorm = MAX(ABS(aW2d(i,j,bi,bj)),myNorm)
125 myNorm = MAX(ABS(aS2d(i,j,bi,bj)),myNorm)
126 ENDDO
127 ENDDO
128 ENDDO
129 ENDDO
130 _GLOBAL_MAX_RS( myNorm, myThid )
131 IF ( myNorm .NE. 0. _d 0 ) THEN
132 myNorm = 1. _d 0/myNorm
133 ELSE
134 myNorm = 1. _d 0
135 ENDIF
136 DO bj=myByLo(myThid),myByHi(myThid)
137 DO bi=myBxLo(myThid),myBxHi(myThid)
138 DO j=1,sNy
139 DO i=1,sNx
140 aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj)*myNorm
141 aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj)*myNorm
142 ENDDO
143 ENDDO
144 ENDDO
145 ENDDO
146
147 C-- Update overlap regions
148 CcnhDebugStarts
149 C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.1' , 1, myThid )
150 C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.1' , 1, myThid )
151 CcnhDebugEnds
152 CALL EXCH_UV_XY_RS( aW2d, aS2d, .FALSE., myThid )
153 CcnhDebugStarts
154 C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.2' , 1, myThid )
155 C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.2' , 1, myThid )
156 CcnhDebugEnds
157
158 _BEGIN_MASTER(myThid)
159 C-- set global parameter in common block:
160 cg2dNorm = myNorm
161 C-- Define the solver tolerance in the appropriate Unit :
162 cg2dNormaliseRHS = cg2dTargetResWunit.LE.0.
163 IF (cg2dNormaliseRHS) THEN
164 C- when using a normalisation of RHS, tolerance has no unit => no conversion
165 cg2dTolerance = cg2dTargetResidual
166 ELSE
167 C- convert Target-Residual (in W unit) to cg2d-solver residual unit [m^2/s^2]
168 cg2dTolerance = cg2dNorm * cg2dTargetResWunit
169 & * globalArea / deltaTmom
170 ENDIF
171 _END_MASTER(myThid)
172
173 CcnhDebugStarts
174 _BEGIN_MASTER( myThid )
175 WRITE(msgBuf,'(2A,1PE23.16)') 'INI_CG2D: ',
176 & 'CG2D normalisation factor = ', cg2dNorm
177 CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
178 IF (.NOT.cg2dNormaliseRHS) THEN
179 WRITE(msgBuf,'(2A,1PE22.15,A,1PE16.10,A)') 'INI_CG2D: ',
180 & 'cg2dTolerance =', cg2dTolerance, ' (Area=',globalArea,')'
181 CALL PRINT_MESSAGE(msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
182 ENDIF
183 WRITE(msgBuf,*) ' '
184 CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
185 _END_MASTER( myThid )
186 CcnhDebugEnds
187
188 C-- Initialise preconditioner
189 C Note. 20th May 1998
190 C I made a weird discovery! In the model paper we argue
191 C for the form of the preconditioner used here ( see
192 C A Finite-volume, Incompressible Navier-Stokes Model
193 C ...., Marshall et. al ). The algebra gives a simple
194 C 0.5 factor for the averaging of ac and aCw to get a
195 C symmettric pre-conditioner. By using a factor of 0.51
196 C i.e. scaling the off-diagonal terms in the
197 C preconditioner down slightly I managed to get the
198 C number of iterations for convergence in a test case to
199 C drop form 192 -> 134! Need to investigate this further!
200 C For now I have introduced a parameter cg2dpcOffDFac which
201 C defaults to 0.51 but can be set at runtime.
202 DO bj=myByLo(myThid),myByHi(myThid)
203 DO bi=myBxLo(myThid),myBxHi(myThid)
204 C- calculate and store solver main diagonal :
205 DO j=0,sNy+1
206 DO i=0,sNx+1
207 ks = ksurfC(i,j,bi,bj)
208 aC2d(i,j,bi,bj) = -(
209 & aW2d(i,j,bi,bj) + aW2d(i+1,j, bi,bj)
210 & +aS2d(i,j,bi,bj) + aS2d( i,j+1,bi,bj)
211 & +freeSurfFac*myNorm*recip_Bo(i,j,bi,bj)*deepFac2F(ks)
212 & *rA(i,j,bi,bj)/deltaTMom/deltaTfreesurf
213 & )
214 ENDDO
215 ENDDO
216 DO j=1,sNy
217 DO i=1,sNx
218 aC = aC2d( i, j, bi,bj)
219 aCs = aC2d( i,j-1,bi,bj)
220 aCw = aC2d(i-1,j, bi,bj)
221 IF ( aC .EQ. 0. ) THEN
222 pC(i,j,bi,bj) = 1. _d 0
223 ELSE
224 pC(i,j,bi,bj) = 1. _d 0 / aC
225 ENDIF
226 IF ( aC + aCw .EQ. 0. ) THEN
227 pW(i,j,bi,bj) = 0.
228 ELSE
229 pW(i,j,bi,bj) =
230 & -aW2d(i,j,bi,bj)/((cg2dpcOffDFac *(aCw+aC))**2 )
231 ENDIF
232 IF ( aC + aCs .EQ. 0. ) THEN
233 pS(i,j,bi,bj) = 0.
234 ELSE
235 pS(i,j,bi,bj) =
236 & -aS2d(i,j,bi,bj)/((cg2dpcOffDFac *(aCs+aC))**2 )
237 ENDIF
238 C pC(i,j,bi,bj) = 1.
239 C pW(i,j,bi,bj) = 0.
240 C pS(i,j,bi,bj) = 0.
241 ENDDO
242 ENDDO
243 ENDDO
244 ENDDO
245 C-- Update overlap regions
246 CALL EXCH_XY_RS( pC, myThid )
247 CALL EXCH_UV_XY_RS( pW, pS, .FALSE., myThid )
248 CcnhDebugStarts
249 C CALL PLOT_FIELD_XYRS( pC, 'pC INI_CG2D.2' , 1, myThid )
250 C CALL PLOT_FIELD_XYRS( pW, 'pW INI_CG2D.2' , 1, myThid )
251 C CALL PLOT_FIELD_XYRS( pS, 'pS INI_CG2D.2' , 1, myThid )
252 CcnhDebugEnds
253
254 RETURN
255 END

  ViewVC Help
Powered by ViewVC 1.1.22