/[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.29 - (hide annotations) (download)
Sun Feb 4 14:38:47 2001 UTC (23 years, 4 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint35
Changes since 1.28: +2 -1 lines
Made sure each .F and .h file had
the CVS keywords Header and Name at its start.
Most had header but very few currently have Name, so
lots of changes!

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

  ViewVC Help
Powered by ViewVC 1.1.22