/[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.22 - (hide annotations) (download)
Tue Dec 8 19:44:28 1998 UTC (25 years, 5 months ago) by adcroft
Branch: MAIN
Changes since 1.21: +16 -1 lines
Implementation of Open Boundaries:
 o new source code: ini_obcs.F set_obcs.F apply_obcs1.F apply_obcs2.F
                    OBCS.h
 o modified code at a few points, key changes are in
    dynamcis.F the_model_main.F and ini_cg2d.F
 o documentation in OBCS.h and doc/OpenBound.*

1 adcroft 1.22 C $Header: /u/gcmpack/models/MITgcmUV/model/src/ini_cg2d.F,v 1.21 1998/11/06 22:44:46 cnh Exp $
2 cnh 1.1
3 adcroft 1.13 #include "CPP_OPTIONS.h"
4 cnh 1.1
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 cnh 1.12 #include "DYNVARS.h"
21 cnh 1.1 #include "CG2D.h"
22 adcroft 1.22 #include "OBCS.h"
23 cnh 1.1
24     C === Routine arguments ===
25     C myThid - Thread no. that called this routine.
26     INTEGER myThid
27     CEndOfInterface
28    
29     C === Local variables ===
30     C xG, yG - Global coordinate location.
31     C zG
32     C iG, jG - Global coordinate index
33     C bi,bj - Loop counters
34     C faceArea - Temporary used to hold cell face areas.
35     C I,J,K
36     C myNorm - Work variable used in clculating normalisation factor
37     CHARACTER*(MAX_LEN_MBUF) msgBuf
38     INTEGER bi, bj
39     INTEGER I, J, K
40 cnh 1.7 _RL faceArea
41 cnh 1.15 _RS myNorm
42 cnh 1.4 _RL aC, aCw, aCs
43 cnh 1.1
44     C-- Initialise laplace operator
45     C aW2d: integral in Z Ax/dX
46     C aS2d: integral in Z Ay/dY
47     myNorm = 0. _d 0
48     DO bj=myByLo(myThid),myByHi(myThid)
49     DO bi=myBxLo(myThid),myBxHi(myThid)
50     DO J=1,sNy
51     DO I=1,sNx
52     aW2d(I,J,bi,bj) = 0. _d 0
53     aS2d(I,J,bi,bj) = 0. _d 0
54     ENDDO
55     ENDDO
56 cnh 1.17 DO K=1,Nr
57 cnh 1.1 DO J=1,sNy
58     DO I=1,sNx
59 cnh 1.20 faceArea = _dyG(I,J,bi,bj)*drF(K)
60     & *_hFacW(I,J,K,bi,bj)
61 cnh 1.1 aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj) +
62 cnh 1.17 & gBaro*faceArea*recip_dxC(I,J,bi,bj)
63 cnh 1.20 faceArea = _dxG(I,J,bi,bj)*drF(K)
64     & *_hFacS(I,J,K,bi,bj)
65 cnh 1.1 aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj) +
66 cnh 1.17 & gBaro*faceArea*recip_dyC(I,J,bi,bj)
67 cnh 1.1 ENDDO
68     ENDDO
69     ENDDO
70 adcroft 1.22 IF (openBoundaries) THEN
71     DO I=1,sNx
72     IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj),bi,bj)=0.
73     IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj)+1,bi,bj)=0.
74     IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj)+1,bi,bj)=0.
75     IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj),bi,bj)=0.
76     ENDDO
77     DO J=1,sNy
78     IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj),J,bi,bj)=0.
79     IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj)+1,J,bi,bj)=0.
80     IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj)+1,J,bi,bj)=0.
81     IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj),J,bi,bj)=0.
82     ENDDO
83     ENDIF
84 cnh 1.1 DO J=1,sNy
85     DO I=1,sNx
86     myNorm = MAX(ABS(aW2d(I,J,bi,bj)),myNorm)
87     myNorm = MAX(ABS(aS2d(I,J,bi,bj)),myNorm)
88     ENDDO
89     ENDDO
90     ENDDO
91     ENDDO
92     cg2dNbuf(1,myThid) = myNorm
93 adcroft 1.16 _GLOBAL_MAX_R4( cg2dNbuf, myNorm, myThid )
94 cnh 1.1 IF ( cg2dNbuf(1,1) .NE. 0. _d 0 ) THEN
95     myNorm = 1. _d 0/cg2dNbuf(1,1)
96     ELSE
97     myNorm = 1. _d 0
98     ENDIF
99 cnh 1.5 cg2dNorm = myNorm
100 cnh 1.1 _BEGIN_MASTER( myThid )
101     CcnhDebugStarts
102 cnh 1.20 WRITE(msgBuf,'(A,E40.25)') '// CG2D normalisation factor = ',
103     & cg2dNorm
104 cnh 1.3 CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
105     WRITE(msgBuf,*) ' '
106 cnh 1.1 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 cnh 1.14 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 cnh 1.1 CcnhDebugEnds
125     _EXCH_XY_R4(aW2d, myThid)
126     _EXCH_XY_R4(aS2d, myThid)
127     CcnhDebugStarts
128 cnh 1.18 CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.2' , 1, myThid )
129     CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.2' , 1, myThid )
130 cnh 1.1 CcnhDebugEnds
131    
132     C-- Initialise preconditioner
133 cnh 1.4 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 cnh 1.1 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 cnh 1.4 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 cnh 1.19 & +freeSurfFac*myNorm* horiVertRatio*
155 cnh 1.17 & rA(I,J,bi,bj)/deltaTMom/deltaTMom
156 cnh 1.4 & )
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 cnh 1.19 & +freeSurfFac*myNorm* horiVertRatio*
161 cnh 1.17 & rA(I,J-1,bi,bj)/deltaTMom/deltaTMom
162 cnh 1.4 & )
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 cnh 1.19 & +freeSurfFac*myNorm* horiVertRatio*
167 cnh 1.17 & rA(I-1,J,bi,bj)/deltaTMom/deltaTMom
168 cnh 1.4 & )
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 cnh 1.6 C pC(I,J,bi,bj) = 1.
187     C pW(I,J,bi,bj) = 0.
188     C pS(I,J,bi,bj) = 0.
189 cnh 1.1 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 cnh 1.18 CcnhDebugStarts
198     CALL PLOT_FIELD_XYRS( pC, 'pC INI_CG2D.2' , 1, myThid )
199     CALL PLOT_FIELD_XYRS( pW, 'pW INI_CG2D.2' , 1, myThid )
200     CALL PLOT_FIELD_XYRS( pS, 'pS INI_CG2D.2' , 1, myThid )
201     CcnhDebugEnds
202 cnh 1.1
203 cnh 1.5 C-- Set default values for initial guess and RHS
204 cnh 1.4 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 cnh 1.5 cg2d_b(I,J,bi,bj) = 0. _d 0
211 cnh 1.21 #ifdef INCLUDE_CD_CODE
212 cnh 1.12 cg2d_xNM1(I,J,bi,bj) = 0. _d 0
213     #endif
214 cnh 1.4 ENDDO
215 cnh 1.1 ENDDO
216     ENDDO
217     ENDDO
218 cnh 1.4 C-- Update overlap regions
219     _EXCH_XY_R8(cg2d_x, myThid)
220 cnh 1.5 _EXCH_XY_R8(cg2d_b, myThid)
221 cnh 1.21 #ifdef INCLUDE_CD_CODE
222 cnh 1.12 _EXCH_XY_R8(cg2d_xNM1, myThid)
223     #endif
224 cnh 1.4 ENDIF
225 cnh 1.1
226     RETURN
227     END

  ViewVC Help
Powered by ViewVC 1.1.22