/[MITgcm]/MITgcm/model/src/ini_cg2d.F
ViewVC logotype

Annotation of /MITgcm/model/src/ini_cg2d.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.44 - (hide annotations) (download)
Tue Oct 17 18:52:34 2006 UTC (17 years, 7 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58q_post, checkpoint58r_post
Changes since 1.43: +5 -2 lines
clean-up multi-threaded problems (reported by debugger tcheck on ACES).

1 jmc 1.44 C $Header: /u/gcmpack/MITgcm/model/src/ini_cg2d.F,v 1.43 2005/07/20 22:33:02 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     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 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 adcroft 1.26 #ifdef ALLOW_OBCS
32 adcroft 1.22 #include "OBCS.h"
33 adcroft 1.26 #endif
34 cnh 1.1
35 cnh 1.37 C !INPUT/OUTPUT PARAMETERS:
36 cnh 1.1 C === Routine arguments ===
37     C myThid - Thread no. that called this routine.
38     INTEGER myThid
39    
40 cnh 1.37 C !LOCAL VARIABLES:
41 cnh 1.1 C === Local variables ===
42     C xG, yG - Global coordinate location.
43     C zG
44     C iG, jG - Global coordinate index
45     C bi,bj - Loop counters
46     C faceArea - Temporary used to hold cell face areas.
47     C I,J,K
48 adcroft 1.34 C myNorm - Work variable used in calculating normalisation factor
49     C sumArea - Work variable used to compute the total Domain Area
50 cnh 1.1 CHARACTER*(MAX_LEN_MBUF) msgBuf
51     INTEGER bi, bj
52     INTEGER I, J, K
53 cnh 1.7 _RL faceArea
54 cnh 1.15 _RS myNorm
55 cnh 1.4 _RL aC, aCw, aCs
56 cnh 1.37 CEOP
57 jmc 1.38
58     C-- Initialize arrays in common blocs (CG2D.h) ; not really necessary
59     C but safer when EXCH do not fill all the overlap regions.
60     DO bj=myByLo(myThid),myByHi(myThid)
61     DO bi=myBxLo(myThid),myBxHi(myThid)
62     DO J=1-OLy,sNy+OLy
63     DO I=1-OLx,sNx+OLx
64     aW2d(I,J,bi,bj) = 0. _d 0
65     aS2d(I,J,bi,bj) = 0. _d 0
66     pW(I,J,bi,bj) = 0. _d 0
67     pS(I,J,bi,bj) = 0. _d 0
68     pC(I,J,bi,bj) = 0. _d 0
69     cg2d_q(I,J,bi,bj) = 0. _d 0
70     ENDDO
71     ENDDO
72     DO J=1-1,sNy+1
73     DO I=1-1,sNx+1
74     cg2d_r(I,J,bi,bj) = 0. _d 0
75     cg2d_s(I,J,bi,bj) = 0. _d 0
76     ENDDO
77     ENDDO
78     ENDDO
79     ENDDO
80 jmc 1.31
81 cnh 1.1 C-- Initialise laplace operator
82     C aW2d: integral in Z Ax/dX
83     C aS2d: integral in Z Ay/dY
84     myNorm = 0. _d 0
85     DO bj=myByLo(myThid),myByHi(myThid)
86     DO bi=myBxLo(myThid),myBxHi(myThid)
87     DO J=1,sNy
88     DO I=1,sNx
89     aW2d(I,J,bi,bj) = 0. _d 0
90     aS2d(I,J,bi,bj) = 0. _d 0
91     ENDDO
92     ENDDO
93 cnh 1.17 DO K=1,Nr
94 cnh 1.1 DO J=1,sNy
95     DO I=1,sNx
96 cnh 1.20 faceArea = _dyG(I,J,bi,bj)*drF(K)
97     & *_hFacW(I,J,K,bi,bj)
98 cnh 1.1 aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj) +
99 jmc 1.32 & implicSurfPress*implicDiv2DFlow
100     & *faceArea*recip_dxC(I,J,bi,bj)
101 cnh 1.20 faceArea = _dxG(I,J,bi,bj)*drF(K)
102     & *_hFacS(I,J,K,bi,bj)
103 cnh 1.1 aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj) +
104 jmc 1.32 & implicSurfPress*implicDiv2DFlow
105     & *faceArea*recip_dyC(I,J,bi,bj)
106 cnh 1.1 ENDDO
107     ENDDO
108     ENDDO
109 adcroft 1.26 #ifdef ALLOW_OBCS
110 adcroft 1.28 IF (useOBCS) THEN
111 adcroft 1.22 DO I=1,sNx
112     IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj),bi,bj)=0.
113     IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj)+1,bi,bj)=0.
114     IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj)+1,bi,bj)=0.
115     IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj),bi,bj)=0.
116     ENDDO
117     DO J=1,sNy
118     IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj),J,bi,bj)=0.
119     IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj)+1,J,bi,bj)=0.
120     IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj)+1,J,bi,bj)=0.
121     IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj),J,bi,bj)=0.
122     ENDDO
123     ENDIF
124 adcroft 1.26 #endif
125 cnh 1.1 DO J=1,sNy
126     DO I=1,sNx
127     myNorm = MAX(ABS(aW2d(I,J,bi,bj)),myNorm)
128     myNorm = MAX(ABS(aS2d(I,J,bi,bj)),myNorm)
129     ENDDO
130     ENDDO
131     ENDDO
132     ENDDO
133 adcroft 1.25 _GLOBAL_MAX_R4( myNorm, myThid )
134     IF ( myNorm .NE. 0. _d 0 ) THEN
135     myNorm = 1. _d 0/myNorm
136 cnh 1.1 ELSE
137     myNorm = 1. _d 0
138     ENDIF
139     DO bj=myByLo(myThid),myByHi(myThid)
140     DO bi=myBxLo(myThid),myBxHi(myThid)
141     DO J=1,sNy
142     DO I=1,sNx
143     aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj)*myNorm
144     aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj)*myNorm
145     ENDDO
146     ENDDO
147     ENDDO
148     ENDDO
149    
150     C-- Update overlap regions
151     CcnhDebugStarts
152 cnh 1.14 C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.1' , 1, myThid )
153     C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.1' , 1, myThid )
154 cnh 1.1 CcnhDebugEnds
155 adcroft 1.34 c _EXCH_XY_R4(aW2d, myThid)
156     c _EXCH_XY_R4(aS2d, myThid)
157     CALL EXCH_UV_XY_RS(aW2d,aS2d,.FALSE.,myThid)
158 cnh 1.1 CcnhDebugStarts
159 adcroft 1.24 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 cnh 1.1 CcnhDebugEnds
162    
163 jmc 1.44 _BEGIN_MASTER(myThid)
164     C-- set global parameter in common block:
165     cg2dNorm = myNorm
166 adcroft 1.34 C-- Define the solver tolerance in the appropriate Unit :
167 jmc 1.43 cg2dNormaliseRHS = cg2dTargetResWunit.LE.0.
168 adcroft 1.34 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 jmc 1.43 & * globalArea / deltaTmom
175 adcroft 1.34 ENDIF
176 jmc 1.44 _END_MASTER(myThid)
177 jmc 1.43
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 adcroft 1.34
193 cnh 1.1 C-- Initialise preconditioner
194 cnh 1.4 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 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
208     DO bi=myBxLo(myThid),myBxHi(myThid)
209     DO J=1,sNy
210     DO I=1,sNx
211     pC(I,J,bi,bj) = 1. _d 0
212 cnh 1.4 aC = -(
213     & aW2d(I,J,bi,bj) + aW2d(I+1,J ,bi,bj)
214 adcroft 1.24 & +aS2d(I,J,bi,bj) + aS2d(I ,J+1,bi,bj)
215 jmc 1.32 & +freeSurfFac*myNorm*recip_Bo(I,J,bi,bj)*
216 adcroft 1.39 & rA(I,J,bi,bj)/deltaTMom/deltaTfreesurf
217 cnh 1.4 & )
218     aCs = -(
219     & aW2d(I,J-1,bi,bj) + aW2d(I+1,J-1,bi,bj)
220     & +aS2d(I,J-1,bi,bj) + aS2d(I ,J ,bi,bj)
221 jmc 1.32 & +freeSurfFac*myNorm*recip_Bo(I,J-1,bi,bj)*
222 adcroft 1.39 & rA(I,J-1,bi,bj)/deltaTMom/deltaTfreesurf
223 cnh 1.4 & )
224     aCw = -(
225     & aW2d(I-1,J,bi,bj) + aW2d(I ,J ,bi,bj)
226     & +aS2d(I-1,J,bi,bj) + aS2d(I-1,J+1,bi,bj)
227 jmc 1.32 & +freeSurfFac*myNorm*recip_Bo(I-1,J,bi,bj)*
228 adcroft 1.39 & rA(I-1,J,bi,bj)/deltaTMom/deltaTfreesurf
229 cnh 1.4 & )
230     IF ( aC .EQ. 0. ) THEN
231 adcroft 1.27 pC(I,J,bi,bj) = 1. _d 0
232 cnh 1.4 ELSE
233     pC(I,J,bi,bj) = 1. _d 0 / aC
234     ENDIF
235     IF ( aC + aCw .EQ. 0. ) THEN
236     pW(I,J,bi,bj) = 0.
237     ELSE
238     pW(I,J,bi,bj) =
239     & -aW2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCw+aC))**2 )
240     ENDIF
241     IF ( aC + aCs .EQ. 0. ) THEN
242     pS(I,J,bi,bj) = 0.
243     ELSE
244     pS(I,J,bi,bj) =
245     & -aS2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCs+aC))**2 )
246     ENDIF
247 cnh 1.6 C pC(I,J,bi,bj) = 1.
248     C pW(I,J,bi,bj) = 0.
249     C pS(I,J,bi,bj) = 0.
250 cnh 1.1 ENDDO
251     ENDDO
252     ENDDO
253     ENDDO
254     C-- Update overlap regions
255     _EXCH_XY_R4(pC, myThid)
256 adcroft 1.34 c _EXCH_XY_R4(pW, myThid)
257     c _EXCH_XY_R4(pS, myThid)
258     CALL EXCH_UV_XY_RS(pW,pS,.FALSE.,myThid)
259 cnh 1.18 CcnhDebugStarts
260 adcroft 1.24 C CALL PLOT_FIELD_XYRS( pC, 'pC INI_CG2D.2' , 1, myThid )
261     C CALL PLOT_FIELD_XYRS( pW, 'pW INI_CG2D.2' , 1, myThid )
262     C CALL PLOT_FIELD_XYRS( pS, 'pS INI_CG2D.2' , 1, myThid )
263 cnh 1.18 CcnhDebugEnds
264 cnh 1.1
265     RETURN
266     END

  ViewVC Help
Powered by ViewVC 1.1.22