1 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/ini_cg2d.F,v 1.27.2.1 2001/01/30 21:02:59 adcroft Exp $ |
2 |
|
3 |
#include "CPP_OPTIONS.h" |
4 |
|
5 |
CStartOfInterface |
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 |
#include "DYNVARS.h" |
22 |
#include "CG2D.h" |
23 |
#ifdef ALLOW_OBCS |
24 |
#include "OBCS.h" |
25 |
#endif |
26 |
|
27 |
C === Routine arguments === |
28 |
C myThid - Thread no. that called this routine. |
29 |
INTEGER myThid |
30 |
CEndOfInterface |
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 clculating normalisation factor |
40 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
41 |
INTEGER bi, bj |
42 |
INTEGER I, J, K |
43 |
_RL faceArea |
44 |
_RS myNorm |
45 |
_RL aC, aCw, aCs |
46 |
|
47 |
C-- Initialise laplace operator |
48 |
C aW2d: integral in Z Ax/dX |
49 |
C aS2d: integral in Z Ay/dY |
50 |
myNorm = 0. _d 0 |
51 |
DO bj=myByLo(myThid),myByHi(myThid) |
52 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
53 |
DO J=1,sNy |
54 |
DO I=1,sNx |
55 |
aW2d(I,J,bi,bj) = 0. _d 0 |
56 |
aS2d(I,J,bi,bj) = 0. _d 0 |
57 |
ENDDO |
58 |
ENDDO |
59 |
DO K=1,Nr |
60 |
DO J=1,sNy |
61 |
DO I=1,sNx |
62 |
faceArea = _dyG(I,J,bi,bj)*drF(K) |
63 |
& *_hFacW(I,J,K,bi,bj) |
64 |
aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj) + |
65 |
& gBaro*faceArea*recip_dxC(I,J,bi,bj) |
66 |
faceArea = _dxG(I,J,bi,bj)*drF(K) |
67 |
& *_hFacS(I,J,K,bi,bj) |
68 |
aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj) + |
69 |
& gBaro*faceArea*recip_dyC(I,J,bi,bj) |
70 |
ENDDO |
71 |
ENDDO |
72 |
ENDDO |
73 |
#ifdef ALLOW_OBCS |
74 |
IF (useOBCS) THEN |
75 |
DO I=1,sNx |
76 |
IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj),bi,bj)=0. |
77 |
IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj)+1,bi,bj)=0. |
78 |
IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj)+1,bi,bj)=0. |
79 |
IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj),bi,bj)=0. |
80 |
ENDDO |
81 |
DO J=1,sNy |
82 |
IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj),J,bi,bj)=0. |
83 |
IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj)+1,J,bi,bj)=0. |
84 |
IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj)+1,J,bi,bj)=0. |
85 |
IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj),J,bi,bj)=0. |
86 |
ENDDO |
87 |
ENDIF |
88 |
#endif |
89 |
DO J=1,sNy |
90 |
DO I=1,sNx |
91 |
myNorm = MAX(ABS(aW2d(I,J,bi,bj)),myNorm) |
92 |
myNorm = MAX(ABS(aS2d(I,J,bi,bj)),myNorm) |
93 |
ENDDO |
94 |
ENDDO |
95 |
ENDDO |
96 |
ENDDO |
97 |
_GLOBAL_MAX_R4( myNorm, myThid ) |
98 |
IF ( myNorm .NE. 0. _d 0 ) THEN |
99 |
myNorm = 1. _d 0/myNorm |
100 |
ELSE |
101 |
myNorm = 1. _d 0 |
102 |
ENDIF |
103 |
cg2dNorm = myNorm |
104 |
_BEGIN_MASTER( myThid ) |
105 |
CcnhDebugStarts |
106 |
WRITE(msgBuf,'(A,E40.25)') '// CG2D normalisation factor = ', |
107 |
& cg2dNorm |
108 |
CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1) |
109 |
WRITE(msgBuf,*) ' ' |
110 |
CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1) |
111 |
CcnhDebugEnds |
112 |
_END_MASTER( myThid ) |
113 |
DO bj=myByLo(myThid),myByHi(myThid) |
114 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
115 |
DO J=1,sNy |
116 |
DO I=1,sNx |
117 |
aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj)*myNorm |
118 |
aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj)*myNorm |
119 |
ENDDO |
120 |
ENDDO |
121 |
ENDDO |
122 |
ENDDO |
123 |
|
124 |
C-- Update overlap regions |
125 |
CcnhDebugStarts |
126 |
C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.1' , 1, myThid ) |
127 |
C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.1' , 1, myThid ) |
128 |
CcnhDebugEnds |
129 |
_EXCH_XY_R4(aW2d, myThid) |
130 |
_EXCH_XY_R4(aS2d, myThid) |
131 |
CcnhDebugStarts |
132 |
C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.2' , 1, myThid ) |
133 |
C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.2' , 1, myThid ) |
134 |
CcnhDebugEnds |
135 |
|
136 |
C-- Initialise preconditioner |
137 |
C Note. 20th May 1998 |
138 |
C I made a weird discovery! In the model paper we argue |
139 |
C for the form of the preconditioner used here ( see |
140 |
C A Finite-volume, Incompressible Navier-Stokes Model |
141 |
C ...., Marshall et. al ). The algebra gives a simple |
142 |
C 0.5 factor for the averaging of ac and aCw to get a |
143 |
C symmettric pre-conditioner. By using a factor of 0.51 |
144 |
C i.e. scaling the off-diagonal terms in the |
145 |
C preconditioner down slightly I managed to get the |
146 |
C number of iterations for convergence in a test case to |
147 |
C drop form 192 -> 134! Need to investigate this further! |
148 |
C For now I have introduced a parameter cg2dpcOffDFac which |
149 |
C defaults to 0.51 but can be set at runtime. |
150 |
DO bj=myByLo(myThid),myByHi(myThid) |
151 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
152 |
DO J=1,sNy |
153 |
DO I=1,sNx |
154 |
pC(I,J,bi,bj) = 1. _d 0 |
155 |
aC = -( |
156 |
& aW2d(I,J,bi,bj) + aW2d(I+1,J ,bi,bj) |
157 |
& +aS2d(I,J,bi,bj) + aS2d(I ,J+1,bi,bj) |
158 |
& +freeSurfFac*myNorm* horiVertRatio* |
159 |
& rA(I,J,bi,bj)/deltaTMom/deltaTMom |
160 |
& ) |
161 |
aCs = -( |
162 |
& aW2d(I,J-1,bi,bj) + aW2d(I+1,J-1,bi,bj) |
163 |
& +aS2d(I,J-1,bi,bj) + aS2d(I ,J ,bi,bj) |
164 |
& +freeSurfFac*myNorm* horiVertRatio* |
165 |
& rA(I,J-1,bi,bj)/deltaTMom/deltaTMom |
166 |
& ) |
167 |
aCw = -( |
168 |
& aW2d(I-1,J,bi,bj) + aW2d(I ,J ,bi,bj) |
169 |
& +aS2d(I-1,J,bi,bj) + aS2d(I-1,J+1,bi,bj) |
170 |
& +freeSurfFac*myNorm* horiVertRatio* |
171 |
& rA(I-1,J,bi,bj)/deltaTMom/deltaTMom |
172 |
& ) |
173 |
IF ( aC .EQ. 0. ) THEN |
174 |
pC(I,J,bi,bj) = 1. _d 0 |
175 |
ELSE |
176 |
pC(I,J,bi,bj) = 1. _d 0 / aC |
177 |
ENDIF |
178 |
IF ( aC + aCw .EQ. 0. ) THEN |
179 |
pW(I,J,bi,bj) = 0. |
180 |
ELSE |
181 |
pW(I,J,bi,bj) = |
182 |
& -aW2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCw+aC))**2 ) |
183 |
ENDIF |
184 |
IF ( aC + aCs .EQ. 0. ) THEN |
185 |
pS(I,J,bi,bj) = 0. |
186 |
ELSE |
187 |
pS(I,J,bi,bj) = |
188 |
& -aS2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCs+aC))**2 ) |
189 |
ENDIF |
190 |
C pC(I,J,bi,bj) = 1. |
191 |
C pW(I,J,bi,bj) = 0. |
192 |
C pS(I,J,bi,bj) = 0. |
193 |
ENDDO |
194 |
ENDDO |
195 |
ENDDO |
196 |
ENDDO |
197 |
C-- Update overlap regions |
198 |
_EXCH_XY_R4(pC, myThid) |
199 |
_EXCH_XY_R4(pW, myThid) |
200 |
_EXCH_XY_R4(pS, myThid) |
201 |
CcnhDebugStarts |
202 |
C CALL PLOT_FIELD_XYRS( pC, 'pC INI_CG2D.2' , 1, myThid ) |
203 |
C CALL PLOT_FIELD_XYRS( pW, 'pW INI_CG2D.2' , 1, myThid ) |
204 |
C CALL PLOT_FIELD_XYRS( pS, 'pS INI_CG2D.2' , 1, myThid ) |
205 |
CcnhDebugEnds |
206 |
|
207 |
C-- Set default values for initial guess and RHS |
208 |
IF ( startTime .EQ. 0 ) THEN |
209 |
DO bj=myByLo(myThid),myByHi(myThid) |
210 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
211 |
DO J=1,sNy |
212 |
DO I=1,sNx |
213 |
cg2d_x(I,J,bi,bj) = 0. _d 0 |
214 |
cg2d_b(I,J,bi,bj) = 0. _d 0 |
215 |
#ifdef INCLUDE_CD_CODE |
216 |
cg2d_xNM1(I,J,bi,bj) = 0. _d 0 |
217 |
#endif |
218 |
ENDDO |
219 |
ENDDO |
220 |
ENDDO |
221 |
ENDDO |
222 |
C-- Update overlap regions |
223 |
_EXCH_XY_R8(cg2d_x, myThid) |
224 |
_EXCH_XY_R8(cg2d_b, myThid) |
225 |
#ifdef INCLUDE_CD_CODE |
226 |
_EXCH_XY_R8(cg2d_xNM1, myThid) |
227 |
#endif |
228 |
ENDIF |
229 |
|
230 |
RETURN |
231 |
END |