1 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/ini_cg2d.F,v 1.8 1998/05/27 05:18:39 cnh Exp $ |
2 |
|
3 |
#include "CPP_EEOPTIONS.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 |
|
15 |
C === Global variables === |
16 |
#include "SIZE.h" |
17 |
#include "EEPARAMS.h" |
18 |
#include "PARAMS.h" |
19 |
#include "GRID.h" |
20 |
#include "CG2D.h" |
21 |
|
22 |
C === Routine arguments === |
23 |
C myThid - Thread no. that called this routine. |
24 |
INTEGER myThid |
25 |
CEndOfInterface |
26 |
|
27 |
C === Local variables === |
28 |
C xG, yG - Global coordinate location. |
29 |
C zG |
30 |
C iG, jG - Global coordinate index |
31 |
C bi,bj - Loop counters |
32 |
C faceArea - Temporary used to hold cell face areas. |
33 |
C I,J,K |
34 |
C myNorm - Work variable used in clculating normalisation factor |
35 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
36 |
INTEGER bi, bj |
37 |
INTEGER I, J, K |
38 |
_RL faceArea |
39 |
_RL myNorm |
40 |
_RL aC, aCw, aCs |
41 |
|
42 |
C-- Initialise laplace operator |
43 |
C aW2d: integral in Z Ax/dX |
44 |
C aS2d: integral in Z Ay/dY |
45 |
myNorm = 0. _d 0 |
46 |
DO bj=myByLo(myThid),myByHi(myThid) |
47 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
48 |
DO J=1,sNy |
49 |
DO I=1,sNx |
50 |
aW2d(I,J,bi,bj) = 0. _d 0 |
51 |
aS2d(I,J,bi,bj) = 0. _d 0 |
52 |
ENDDO |
53 |
ENDDO |
54 |
DO K=1,Nz |
55 |
DO J=1,sNy |
56 |
DO I=1,sNx |
57 |
faceArea = _dyG(I,J,bi,bj)*dzF(K)*HFacW(I,J,K,bi,bj) |
58 |
aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj) + |
59 |
& gBaro*faceArea*rDxC(I,J,bi,bj) |
60 |
faceArea = _dxG(I,J,bi,bj)*dzF(K)*HFacS(I,J,K,bi,bj) |
61 |
aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj) + |
62 |
& gBaro*faceArea*rDyC(I,J,bi,bj) |
63 |
ENDDO |
64 |
ENDDO |
65 |
ENDDO |
66 |
DO J=1,sNy |
67 |
DO I=1,sNx |
68 |
myNorm = MAX(ABS(aW2d(I,J,bi,bj)),myNorm) |
69 |
myNorm = MAX(ABS(aS2d(I,J,bi,bj)),myNorm) |
70 |
ENDDO |
71 |
ENDDO |
72 |
ENDDO |
73 |
ENDDO |
74 |
cg2dNbuf(1,myThid) = myNorm |
75 |
_GLOBAL_MAX_R8( cg2dNbuf, myNorm, myThid ) |
76 |
IF ( cg2dNbuf(1,1) .NE. 0. _d 0 ) THEN |
77 |
myNorm = 1. _d 0/cg2dNbuf(1,1) |
78 |
ELSE |
79 |
myNorm = 1. _d 0 |
80 |
ENDIF |
81 |
cg2dNorm = myNorm |
82 |
_BEGIN_MASTER( myThid ) |
83 |
CcnhDebugStarts |
84 |
WRITE(msgBuf,'(A,F)') '// CG2D normalisation factor = ', cg2dNorm |
85 |
CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1) |
86 |
WRITE(msgBuf,*) ' ' |
87 |
CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1) |
88 |
CcnhDebugEnds |
89 |
_END_MASTER( myThid ) |
90 |
DO bj=myByLo(myThid),myByHi(myThid) |
91 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
92 |
DO J=1,sNy |
93 |
DO I=1,sNx |
94 |
aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj)*myNorm |
95 |
aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj)*myNorm |
96 |
ENDDO |
97 |
ENDDO |
98 |
ENDDO |
99 |
ENDDO |
100 |
|
101 |
C-- Update overlap regions |
102 |
CcnhDebugStarts |
103 |
C CALL PLOT_FIELD_XYR8( aW2d, 'AW2D INI_CG2D.1' , 1, myThid ) |
104 |
C CALL PLOT_FIELD_XYR8( aS2d, 'AS2D INI_CG2D.1' , 1, myThid ) |
105 |
CcnhDebugEnds |
106 |
_EXCH_XY_R4(aW2d, myThid) |
107 |
_EXCH_XY_R4(aS2d, myThid) |
108 |
CcnhDebugStarts |
109 |
C CALL PLOT_FIELD_XYR8( aW2d, 'AW2D INI_CG2D.2' , 1, myThid ) |
110 |
C CALL PLOT_FIELD_XYR8( aS2d, 'AS2D INI_CG2D.2' , 1, myThid ) |
111 |
CcnhDebugEnds |
112 |
|
113 |
C-- Initialise preconditioner |
114 |
C Note. 20th May 1998 |
115 |
C I made a weird discovery! In the model paper we argue |
116 |
C for the form of the preconditioner used here ( see |
117 |
C A Finite-volume, Incompressible Navier-Stokes Model |
118 |
C ...., Marshall et. al ). The algebra gives a simple |
119 |
C 0.5 factor for the averaging of ac and aCw to get a |
120 |
C symmettric pre-conditioner. By using a factor of 0.51 |
121 |
C i.e. scaling the off-diagonal terms in the |
122 |
C preconditioner down slightly I managed to get the |
123 |
C number of iterations for convergence in a test case to |
124 |
C drop form 192 -> 134! Need to investigate this further! |
125 |
C For now I have introduced a parameter cg2dpcOffDFac which |
126 |
C defaults to 0.51 but can be set at runtime. |
127 |
DO bj=myByLo(myThid),myByHi(myThid) |
128 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
129 |
DO J=1,sNy |
130 |
DO I=1,sNx |
131 |
pC(I,J,bi,bj) = 1. _d 0 |
132 |
aC = -( |
133 |
& aW2d(I,J,bi,bj) + aW2d(I+1,J ,bi,bj) |
134 |
& +aS2d(I,J,bi,bj) + aS2D(I ,J+1,bi,bj) |
135 |
& +freeSurfFac*myNorm* |
136 |
& zA(I,J,bi,bj)/deltaTMom/deltaTMom |
137 |
& ) |
138 |
aCs = -( |
139 |
& aW2d(I,J-1,bi,bj) + aW2d(I+1,J-1,bi,bj) |
140 |
& +aS2d(I,J-1,bi,bj) + aS2d(I ,J ,bi,bj) |
141 |
& +freeSurfFac*myNorm* |
142 |
& zA(I,J-1,bi,bj)/deltaTMom/deltaTMom |
143 |
& ) |
144 |
aCw = -( |
145 |
& aW2d(I-1,J,bi,bj) + aW2d(I ,J ,bi,bj) |
146 |
& +aS2d(I-1,J,bi,bj) + aS2d(I-1,J+1,bi,bj) |
147 |
& +freeSurfFac*myNorm* |
148 |
& zA(I-1,J,bi,bj)/deltaTMom/deltaTMom |
149 |
& ) |
150 |
IF ( aC .EQ. 0. ) THEN |
151 |
pC(I,J,bi,bj) = 0. _d 0 |
152 |
ELSE |
153 |
pC(I,J,bi,bj) = 1. _d 0 / aC |
154 |
ENDIF |
155 |
IF ( aC + aCw .EQ. 0. ) THEN |
156 |
pW(I,J,bi,bj) = 0. |
157 |
ELSE |
158 |
pW(I,J,bi,bj) = |
159 |
& -aW2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCw+aC))**2 ) |
160 |
ENDIF |
161 |
IF ( aC + aCs .EQ. 0. ) THEN |
162 |
pS(I,J,bi,bj) = 0. |
163 |
ELSE |
164 |
pS(I,J,bi,bj) = |
165 |
& -aS2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCs+aC))**2 ) |
166 |
ENDIF |
167 |
C pC(I,J,bi,bj) = 1. |
168 |
C pW(I,J,bi,bj) = 0. |
169 |
C pS(I,J,bi,bj) = 0. |
170 |
ENDDO |
171 |
ENDDO |
172 |
ENDDO |
173 |
ENDDO |
174 |
C-- Update overlap regions |
175 |
_EXCH_XY_R4(pC, myThid) |
176 |
_EXCH_XY_R4(pW, myThid) |
177 |
_EXCH_XY_R4(pS, myThid) |
178 |
|
179 |
C-- Set default values for initial guess and RHS |
180 |
IF ( startTime .EQ. 0 ) THEN |
181 |
DO bj=myByLo(myThid),myByHi(myThid) |
182 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
183 |
DO J=1,sNy |
184 |
DO I=1,sNx |
185 |
cg2d_x(I,J,bi,bj) = 0. _d 0 |
186 |
cg2d_b(I,J,bi,bj) = 0. _d 0 |
187 |
ENDDO |
188 |
ENDDO |
189 |
ENDDO |
190 |
ENDDO |
191 |
C-- Update overlap regions |
192 |
_EXCH_XY_R8(cg2d_x, myThid) |
193 |
_EXCH_XY_R8(cg2d_b, myThid) |
194 |
ENDIF |
195 |
|
196 |
RETURN |
197 |
END |