/[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.49 - (hide annotations) (download)
Tue Dec 8 00:10:26 2009 UTC (14 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62c, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62w, checkpoint62x, checkpoint62, checkpoint62b
Changes since 1.48: +89 -96 lines
fix deep-atmosphere with topography (+ re-arrange do-loops)

1 jmc 1.49 C $Header: /u/gcmpack/MITgcm/model/src/ini_cg2d.F,v 1.48 2009/11/23 16:15:54 mlosch Exp $
2 jmc 1.30 C $Name: $
3 cnh 1.1
4 edhill 1.42 #include "PACKAGES_CONFIG.h"
5 adcroft 1.13 #include "CPP_OPTIONS.h"
6 cnh 1.1
7 cnh 1.37 CBOP
8     C !ROUTINE: INI_CG2D
9     C !INTERFACE:
10 cnh 1.1 SUBROUTINE INI_CG2D( myThid )
11 cnh 1.37
12     C !DESCRIPTION: \bv
13     C *==========================================================*
14 jmc 1.47 C | SUBROUTINE INI_CG2D
15     C | o Initialise 2d conjugate gradient solver operators.
16 cnh 1.37 C *==========================================================*
17 jmc 1.47 C | These arrays are purely a function of the basin geom.
18     C | We set then here once and them use then repeatedly.
19 cnh 1.37 C *==========================================================*
20     C \ev
21    
22     C !USES:
23 adcroft 1.23 IMPLICIT NONE
24 cnh 1.1 C === Global variables ===
25     #include "SIZE.h"
26     #include "EEPARAMS.h"
27     #include "PARAMS.h"
28     #include "GRID.h"
29 jmc 1.31 #include "SURFACE.h"
30 adcroft 1.34 #include "CG2D.h"
31 adcroft 1.26 #ifdef ALLOW_OBCS
32 adcroft 1.22 #include "OBCS.h"
33 adcroft 1.26 #endif
34 cnh 1.1
35 cnh 1.37 C !INPUT/OUTPUT PARAMETERS:
36 cnh 1.1 C === Routine arguments ===
37     C myThid - Thread no. that called this routine.
38     INTEGER myThid
39    
40 cnh 1.37 C !LOCAL VARIABLES:
41 cnh 1.1 C === Local variables ===
42 jmc 1.45 C bi,bj :: tile indices
43 jmc 1.49 C i,j,k :: Loop counters
44 cnh 1.1 C faceArea - Temporary used to hold cell face areas.
45 adcroft 1.34 C myNorm - Work variable used in calculating normalisation factor
46 jmc 1.45 C sumArea - Work variable used to compute the total Domain Area
47 cnh 1.1 CHARACTER*(MAX_LEN_MBUF) msgBuf
48     INTEGER bi, bj
49 jmc 1.49 INTEGER i, j, k, ks
50 cnh 1.7 _RL faceArea
51 cnh 1.15 _RS myNorm
52 jmc 1.45 _RS aC, aCw, aCs
53 cnh 1.37 CEOP
54 jmc 1.38
55 jmc 1.45 C-- Initialize arrays in common blocs (CG2D.h) ; not really necessary
56 jmc 1.38 C but safer when EXCH do not fill all the overlap regions.
57     DO bj=myByLo(myThid),myByHi(myThid)
58     DO bi=myBxLo(myThid),myBxHi(myThid)
59 jmc 1.49 DO j=1-OLy,sNy+OLy
60     DO i=1-OLx,sNx+OLx
61     aW2d(i,j,bi,bj) = 0. _d 0
62     aS2d(i,j,bi,bj) = 0. _d 0
63     aC2d(i,j,bi,bj) = 0. _d 0
64     pW(i,j,bi,bj) = 0. _d 0
65     pS(i,j,bi,bj) = 0. _d 0
66     pC(i,j,bi,bj) = 0. _d 0
67 jmc 1.38 ENDDO
68     ENDDO
69 jmc 1.49 DO j=1-1,sNy+1
70     DO i=1-1,sNx+1
71     cg2d_q(i,j,bi,bj) = 0. _d 0
72     cg2d_r(i,j,bi,bj) = 0. _d 0
73     cg2d_s(i,j,bi,bj) = 0. _d 0
74 jmc 1.47 #ifdef ALLOW_CG2D_NSA
75 jmc 1.49 cg2d_z(i,j,bi,bj) = 0. _d 0
76 jmc 1.47 #endif /* ALLOW_CG2D_NSA */
77 mlosch 1.48 #ifdef ALLOW_SRCG
78 jmc 1.49 cg2d_y(i,j,bi,bj) = 0. _d 0
79     cg2d_v(i,j,bi,bj) = 0. _d 0
80 mlosch 1.48 #endif /* ALLOW_SRCG */
81 jmc 1.38 ENDDO
82     ENDDO
83     ENDDO
84     ENDDO
85 jmc 1.31
86 cnh 1.1 C-- Initialise laplace operator
87     C aW2d: integral in Z Ax/dX
88     C aS2d: integral in Z Ay/dY
89     myNorm = 0. _d 0
90     DO bj=myByLo(myThid),myByHi(myThid)
91     DO bi=myBxLo(myThid),myBxHi(myThid)
92 jmc 1.49 DO j=1,sNy
93     DO i=1,sNx
94     aW2d(i,j,bi,bj) = 0. _d 0
95     aS2d(i,j,bi,bj) = 0. _d 0
96 cnh 1.1 ENDDO
97     ENDDO
98 jmc 1.49 DO k=1,Nr
99     DO j=1,sNy
100     DO i=1,sNx
101 jmc 1.45 C deep-model: *deepFacC (faceArea), /deepFacC (recip_dx,y): => no net effect
102 jmc 1.49 faceArea = _dyG(i,j,bi,bj)*drF(k)
103     & *_hFacW(i,j,k,bi,bj)
104     aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj)
105 jmc 1.45 & + implicSurfPress*implicDiv2DFlow
106 jmc 1.49 & *faceArea*recip_dxC(i,j,bi,bj)
107     faceArea = _dxG(i,j,bi,bj)*drF(k)
108     & *_hFacS(i,j,k,bi,bj)
109     aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj)
110 jmc 1.45 & + implicSurfPress*implicDiv2DFlow
111 jmc 1.49 & *faceArea*recip_dyC(i,j,bi,bj)
112 cnh 1.1 ENDDO
113     ENDDO
114     ENDDO
115 adcroft 1.26 #ifdef ALLOW_OBCS
116 adcroft 1.28 IF (useOBCS) THEN
117 jmc 1.49 DO i=1,sNx
118     IF (OB_Jn(i,bi,bj).NE.0) aS2d(i, OB_Jn(i,bi,bj), bi,bj)=0.
119     IF (OB_Jn(i,bi,bj).NE.0) aS2d(i,OB_Jn(i,bi,bj)+1,bi,bj)=0.
120     IF (OB_Js(i,bi,bj).NE.0) aS2d(i,OB_Js(i,bi,bj)+1,bi,bj)=0.
121     IF (OB_Js(i,bi,bj).NE.0) aS2d(i, OB_Js(i,bi,bj), bi,bj)=0.
122 adcroft 1.22 ENDDO
123 jmc 1.49 DO j=1,sNy
124     IF (OB_Ie(j,bi,bj).NE.0) aW2d(OB_Ie(j,bi,bj), j, bi,bj)=0.
125     IF (OB_Ie(j,bi,bj).NE.0) aW2d(OB_Ie(j,bi,bj)+1,j,bi,bj)=0.
126     IF (OB_Iw(j,bi,bj).NE.0) aW2d(OB_Iw(j,bi,bj)+1,j,bi,bj)=0.
127     IF (OB_Iw(j,bi,bj).NE.0) aW2d(OB_Iw(j,bi,bj), j, bi,bj)=0.
128 adcroft 1.22 ENDDO
129     ENDIF
130 adcroft 1.26 #endif
131 jmc 1.49 DO j=1,sNy
132     DO i=1,sNx
133     myNorm = MAX(ABS(aW2d(i,j,bi,bj)),myNorm)
134     myNorm = MAX(ABS(aS2d(i,j,bi,bj)),myNorm)
135 cnh 1.1 ENDDO
136     ENDDO
137     ENDDO
138     ENDDO
139 jmc 1.46 _GLOBAL_MAX_RS( myNorm, myThid )
140 adcroft 1.25 IF ( myNorm .NE. 0. _d 0 ) THEN
141     myNorm = 1. _d 0/myNorm
142 cnh 1.1 ELSE
143     myNorm = 1. _d 0
144     ENDIF
145     DO bj=myByLo(myThid),myByHi(myThid)
146     DO bi=myBxLo(myThid),myBxHi(myThid)
147 jmc 1.49 DO j=1,sNy
148     DO i=1,sNx
149     aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj)*myNorm
150     aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj)*myNorm
151 cnh 1.1 ENDDO
152     ENDDO
153     ENDDO
154     ENDDO
155 jmc 1.45
156 cnh 1.1 C-- Update overlap regions
157     CcnhDebugStarts
158 cnh 1.14 C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.1' , 1, myThid )
159     C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.1' , 1, myThid )
160 cnh 1.1 CcnhDebugEnds
161 jmc 1.47 CALL EXCH_UV_XY_RS( aW2d, aS2d, .FALSE., myThid )
162 cnh 1.1 CcnhDebugStarts
163 adcroft 1.24 C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.2' , 1, myThid )
164     C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.2' , 1, myThid )
165 cnh 1.1 CcnhDebugEnds
166    
167 jmc 1.44 _BEGIN_MASTER(myThid)
168     C-- set global parameter in common block:
169     cg2dNorm = myNorm
170 adcroft 1.34 C-- Define the solver tolerance in the appropriate Unit :
171 jmc 1.43 cg2dNormaliseRHS = cg2dTargetResWunit.LE.0.
172 adcroft 1.34 IF (cg2dNormaliseRHS) THEN
173     C- when using a normalisation of RHS, tolerance has no unit => no conversion
174     cg2dTolerance = cg2dTargetResidual
175     ELSE
176     C- convert Target-Residual (in W unit) to cg2d-solver residual unit [m^2/s^2]
177     cg2dTolerance = cg2dNorm * cg2dTargetResWunit
178 jmc 1.43 & * globalArea / deltaTmom
179 adcroft 1.34 ENDIF
180 jmc 1.44 _END_MASTER(myThid)
181 jmc 1.43
182     CcnhDebugStarts
183     _BEGIN_MASTER( myThid )
184     WRITE(msgBuf,'(2A,1PE23.16)') 'INI_CG2D: ',
185     & 'CG2D normalisation factor = ', cg2dNorm
186     CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
187     IF (.NOT.cg2dNormaliseRHS) THEN
188     WRITE(msgBuf,'(2A,1PE22.15,A,1PE16.10,A)') 'INI_CG2D: ',
189     & 'cg2dTolerance =', cg2dTolerance, ' (Area=',globalArea,')'
190     CALL PRINT_MESSAGE(msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
191     ENDIF
192     WRITE(msgBuf,*) ' '
193     CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
194     _END_MASTER( myThid )
195     CcnhDebugEnds
196 jmc 1.45
197 cnh 1.1 C-- Initialise preconditioner
198 cnh 1.4 C Note. 20th May 1998
199     C I made a weird discovery! In the model paper we argue
200     C for the form of the preconditioner used here ( see
201     C A Finite-volume, Incompressible Navier-Stokes Model
202     C ...., Marshall et. al ). The algebra gives a simple
203     C 0.5 factor for the averaging of ac and aCw to get a
204     C symmettric pre-conditioner. By using a factor of 0.51
205     C i.e. scaling the off-diagonal terms in the
206     C preconditioner down slightly I managed to get the
207     C number of iterations for convergence in a test case to
208     C drop form 192 -> 134! Need to investigate this further!
209     C For now I have introduced a parameter cg2dpcOffDFac which
210     C defaults to 0.51 but can be set at runtime.
211 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
212     DO bi=myBxLo(myThid),myBxHi(myThid)
213 jmc 1.49 C- calculate and store solver main diagonal :
214     DO j=0,sNy+1
215     DO i=0,sNx+1
216     ks = ksurfC(i,j,bi,bj)
217     aC2d(i,j,bi,bj) = -(
218     & aW2d(i,j,bi,bj) + aW2d(i+1,j, bi,bj)
219     & +aS2d(i,j,bi,bj) + aS2d( i,j+1,bi,bj)
220     & +freeSurfFac*myNorm*recip_Bo(i,j,bi,bj)*deepFac2F(ks)
221     & *rA(i,j,bi,bj)/deltaTMom/deltaTfreesurf
222     & )
223     ENDDO
224     ENDDO
225     DO j=1,sNy
226     DO i=1,sNx
227     aC = aC2d( i, j, bi,bj)
228     aCs = aC2d( i,j-1,bi,bj)
229     aCw = aC2d(i-1,j, bi,bj)
230     IF ( aC .EQ. 0. ) THEN
231     pC(i,j,bi,bj) = 1. _d 0
232     ELSE
233     pC(i,j,bi,bj) = 1. _d 0 / aC
234     ENDIF
235     IF ( aC + aCw .EQ. 0. ) THEN
236     pW(i,j,bi,bj) = 0.
237     ELSE
238     pW(i,j,bi,bj) =
239     & -aW2d(i,j,bi,bj)/((cg2dpcOffDFac *(aCw+aC))**2 )
240     ENDIF
241     IF ( aC + aCs .EQ. 0. ) THEN
242     pS(i,j,bi,bj) = 0.
243     ELSE
244     pS(i,j,bi,bj) =
245     & -aS2d(i,j,bi,bj)/((cg2dpcOffDFac *(aCs+aC))**2 )
246     ENDIF
247     C pC(i,j,bi,bj) = 1.
248     C pW(i,j,bi,bj) = 0.
249     C pS(i,j,bi,bj) = 0.
250 cnh 1.1 ENDDO
251     ENDDO
252     ENDDO
253     ENDDO
254     C-- Update overlap regions
255 jmc 1.47 CALL EXCH_XY_RS( pC, myThid )
256     CALL EXCH_UV_XY_RS( pW, pS, .FALSE., myThid )
257 cnh 1.18 CcnhDebugStarts
258 adcroft 1.24 C CALL PLOT_FIELD_XYRS( pC, 'pC INI_CG2D.2' , 1, myThid )
259     C CALL PLOT_FIELD_XYRS( pW, 'pW INI_CG2D.2' , 1, myThid )
260     C CALL PLOT_FIELD_XYRS( pS, 'pS INI_CG2D.2' , 1, myThid )
261 cnh 1.18 CcnhDebugEnds
262 cnh 1.1
263     RETURN
264     END

  ViewVC Help
Powered by ViewVC 1.1.22