/[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.20 - (hide annotations) (download)
Wed Oct 28 03:11:37 1998 UTC (25 years, 7 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint17, checkpoint16
Changes since 1.19: +7 -4 lines
Changes to support
 - g77 compilation under Linux
 - LR(1) form of 64-bit is D or E for constants
 - Modified adjoint of exch with adjoint variables
   acuumulated.

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

  ViewVC Help
Powered by ViewVC 1.1.22