/[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.36 - (show annotations) (download)
Fri Jul 6 21:39:37 2001 UTC (22 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint40pre3, checkpoint40pre7, checkpoint40pre6, checkpoint40pre9, checkpoint40pre8, checkpoint40pre2, checkpoint40pre4, checkpoint40pre5, checkpoint40
Changes since 1.35: +1 -37 lines
 compute Bo_surf(Po_ground,Tref) inside routine INI_LINEAR_PHISURF

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/ini_cg2d.F,v 1.35 2001/06/07 15:53:21 adcroft Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 SUBROUTINE INI_CG2D( myThid )
7 C /==========================================================\
8 C | SUBROUTINE INI_CG2D |
9 C | o Initialise 2d conjugate gradient solver operators. |
10 C |==========================================================|
11 C | These arrays are purely a function of the basin geom. |
12 C | We set then here once and them use then repeatedly. |
13 C \==========================================================/
14 IMPLICIT NONE
15
16 C === Global variables ===
17 #include "SIZE.h"
18 #include "EEPARAMS.h"
19 #include "PARAMS.h"
20 #include "GRID.h"
21 c#include "DYNVARS.h"
22 #include "SURFACE.h"
23 #include "CG2D.h"
24 #ifdef ALLOW_OBCS
25 #include "OBCS.h"
26 #endif
27
28 C === Routine arguments ===
29 C myThid - Thread no. that called this routine.
30 INTEGER myThid
31
32 C === Local variables ===
33 C xG, yG - Global coordinate location.
34 C zG
35 C iG, jG - Global coordinate index
36 C bi,bj - Loop counters
37 C faceArea - Temporary used to hold cell face areas.
38 C I,J,K
39 C myNorm - Work variable used in calculating normalisation factor
40 C sumArea - Work variable used to compute the total Domain Area
41 CHARACTER*(MAX_LEN_MBUF) msgBuf
42 INTEGER bi, bj
43 INTEGER I, J, K
44 _RL faceArea
45 _RS myNorm
46 _RL sumArea
47 _RL aC, aCw, aCs
48
49 C-- Initialise laplace operator
50 C aW2d: integral in Z Ax/dX
51 C aS2d: integral in Z Ay/dY
52 myNorm = 0. _d 0
53 DO bj=myByLo(myThid),myByHi(myThid)
54 DO bi=myBxLo(myThid),myBxHi(myThid)
55 DO J=1,sNy
56 DO I=1,sNx
57 aW2d(I,J,bi,bj) = 0. _d 0
58 aS2d(I,J,bi,bj) = 0. _d 0
59 ENDDO
60 ENDDO
61 DO K=1,Nr
62 DO J=1,sNy
63 DO I=1,sNx
64 faceArea = _dyG(I,J,bi,bj)*drF(K)
65 & *_hFacW(I,J,K,bi,bj)
66 aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj) +
67 & implicSurfPress*implicDiv2DFlow
68 & *faceArea*recip_dxC(I,J,bi,bj)
69 faceArea = _dxG(I,J,bi,bj)*drF(K)
70 & *_hFacS(I,J,K,bi,bj)
71 aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj) +
72 & implicSurfPress*implicDiv2DFlow
73 & *faceArea*recip_dyC(I,J,bi,bj)
74 ENDDO
75 ENDDO
76 ENDDO
77 #ifdef ALLOW_OBCS
78 IF (useOBCS) THEN
79 DO I=1,sNx
80 IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj),bi,bj)=0.
81 IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj)+1,bi,bj)=0.
82 IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj)+1,bi,bj)=0.
83 IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj),bi,bj)=0.
84 ENDDO
85 DO J=1,sNy
86 IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj),J,bi,bj)=0.
87 IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj)+1,J,bi,bj)=0.
88 IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj)+1,J,bi,bj)=0.
89 IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj),J,bi,bj)=0.
90 ENDDO
91 ENDIF
92 #endif
93 DO J=1,sNy
94 DO I=1,sNx
95 myNorm = MAX(ABS(aW2d(I,J,bi,bj)),myNorm)
96 myNorm = MAX(ABS(aS2d(I,J,bi,bj)),myNorm)
97 ENDDO
98 ENDDO
99 ENDDO
100 ENDDO
101 _GLOBAL_MAX_R4( myNorm, myThid )
102 IF ( myNorm .NE. 0. _d 0 ) THEN
103 myNorm = 1. _d 0/myNorm
104 ELSE
105 myNorm = 1. _d 0
106 ENDIF
107 cg2dNorm = myNorm
108 _BEGIN_MASTER( myThid )
109 CcnhDebugStarts
110 WRITE(msgBuf,'(A,E40.25)') '// CG2D normalisation factor = ',
111 & cg2dNorm
112 CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
113 WRITE(msgBuf,*) ' '
114 CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
115 CcnhDebugEnds
116 _END_MASTER( myThid )
117 DO bj=myByLo(myThid),myByHi(myThid)
118 DO bi=myBxLo(myThid),myBxHi(myThid)
119 DO J=1,sNy
120 DO I=1,sNx
121 aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj)*myNorm
122 aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj)*myNorm
123 ENDDO
124 ENDDO
125 ENDDO
126 ENDDO
127
128 C-- Update overlap regions
129 CcnhDebugStarts
130 C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.1' , 1, myThid )
131 C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.1' , 1, myThid )
132 CcnhDebugEnds
133 c _EXCH_XY_R4(aW2d, myThid)
134 c _EXCH_XY_R4(aS2d, myThid)
135 CALL EXCH_UV_XY_RS(aW2d,aS2d,.FALSE.,myThid)
136 CcnhDebugStarts
137 C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.2' , 1, myThid )
138 C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.2' , 1, myThid )
139 CcnhDebugEnds
140
141 C-- Define the solver tolerance in the appropriate Unit :
142 cg2dNormaliseRHS = cg2dTargetResWunit.LE.0
143 IF (cg2dNormaliseRHS) THEN
144 C- when using a normalisation of RHS, tolerance has no unit => no conversion
145 cg2dTolerance = cg2dTargetResidual
146 ELSE
147 C- compute the total Area of the domain :
148 sumArea = 0.
149 DO bj=myByLo(myThid),myByHi(myThid)
150 DO bi=myBxLo(myThid),myBxHi(myThid)
151 DO j=1,sNy
152 DO i=1,sNx
153 IF (Ro_surf(i,j,bi,bj).GT.R_low(i,j,bi,bj)) THEN
154 sumArea = sumArea + rA(i,j,bi,bj)
155 ENDIF
156 ENDDO
157 ENDDO
158 ENDDO
159 ENDDO
160 c WRITE(*,*) ' mythid, sumArea = ', mythid, sumArea
161 _GLOBAL_SUM_R8( sumArea, myThid )
162 C- convert Target-Residual (in W unit) to cg2d-solver residual unit [m^2/s^2]
163 cg2dTolerance = cg2dNorm * cg2dTargetResWunit
164 & * sumArea / deltaTMom
165 WRITE(*,'(2A,1P2E22.14)') ' ini_cg2d: ',
166 & 'sumArea,cg2dTolerance =', sumArea,cg2dTolerance
167 ENDIF
168
169 C-- Initialise preconditioner
170 C Note. 20th May 1998
171 C I made a weird discovery! In the model paper we argue
172 C for the form of the preconditioner used here ( see
173 C A Finite-volume, Incompressible Navier-Stokes Model
174 C ...., Marshall et. al ). The algebra gives a simple
175 C 0.5 factor for the averaging of ac and aCw to get a
176 C symmettric pre-conditioner. By using a factor of 0.51
177 C i.e. scaling the off-diagonal terms in the
178 C preconditioner down slightly I managed to get the
179 C number of iterations for convergence in a test case to
180 C drop form 192 -> 134! Need to investigate this further!
181 C For now I have introduced a parameter cg2dpcOffDFac which
182 C defaults to 0.51 but can be set at runtime.
183 DO bj=myByLo(myThid),myByHi(myThid)
184 DO bi=myBxLo(myThid),myBxHi(myThid)
185 DO J=1,sNy
186 DO I=1,sNx
187 pC(I,J,bi,bj) = 1. _d 0
188 aC = -(
189 & aW2d(I,J,bi,bj) + aW2d(I+1,J ,bi,bj)
190 & +aS2d(I,J,bi,bj) + aS2d(I ,J+1,bi,bj)
191 & +freeSurfFac*myNorm*recip_Bo(I,J,bi,bj)*
192 & rA(I,J,bi,bj)/deltaTMom/deltaTMom
193 & )
194 aCs = -(
195 & aW2d(I,J-1,bi,bj) + aW2d(I+1,J-1,bi,bj)
196 & +aS2d(I,J-1,bi,bj) + aS2d(I ,J ,bi,bj)
197 & +freeSurfFac*myNorm*recip_Bo(I,J-1,bi,bj)*
198 & rA(I,J-1,bi,bj)/deltaTMom/deltaTMom
199 & )
200 aCw = -(
201 & aW2d(I-1,J,bi,bj) + aW2d(I ,J ,bi,bj)
202 & +aS2d(I-1,J,bi,bj) + aS2d(I-1,J+1,bi,bj)
203 & +freeSurfFac*myNorm*recip_Bo(I-1,J,bi,bj)*
204 & rA(I-1,J,bi,bj)/deltaTMom/deltaTMom
205 & )
206 IF ( aC .EQ. 0. ) THEN
207 pC(I,J,bi,bj) = 1. _d 0
208 ELSE
209 pC(I,J,bi,bj) = 1. _d 0 / aC
210 ENDIF
211 IF ( aC + aCw .EQ. 0. ) THEN
212 pW(I,J,bi,bj) = 0.
213 ELSE
214 pW(I,J,bi,bj) =
215 & -aW2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCw+aC))**2 )
216 ENDIF
217 IF ( aC + aCs .EQ. 0. ) THEN
218 pS(I,J,bi,bj) = 0.
219 ELSE
220 pS(I,J,bi,bj) =
221 & -aS2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCs+aC))**2 )
222 ENDIF
223 C pC(I,J,bi,bj) = 1.
224 C pW(I,J,bi,bj) = 0.
225 C pS(I,J,bi,bj) = 0.
226 ENDDO
227 ENDDO
228 ENDDO
229 ENDDO
230 C-- Update overlap regions
231 _EXCH_XY_R4(pC, myThid)
232 c _EXCH_XY_R4(pW, myThid)
233 c _EXCH_XY_R4(pS, myThid)
234 CALL EXCH_UV_XY_RS(pW,pS,.FALSE.,myThid)
235 CcnhDebugStarts
236 C CALL PLOT_FIELD_XYRS( pC, 'pC INI_CG2D.2' , 1, myThid )
237 C CALL PLOT_FIELD_XYRS( pW, 'pW INI_CG2D.2' , 1, myThid )
238 C CALL PLOT_FIELD_XYRS( pS, 'pS INI_CG2D.2' , 1, myThid )
239 CcnhDebugEnds
240
241 RETURN
242 END

  ViewVC Help
Powered by ViewVC 1.1.22