/[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.25 - (show annotations) (download)
Tue May 18 17:36:55 1999 UTC (25 years ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint22
Changes since 1.24: +4 -5 lines
Changed number of arguments to GLOBAL_SUM and GLOBAL_MAX to two.
Instigated by Ralf.

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

  ViewVC Help
Powered by ViewVC 1.1.22