/[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.37 - (hide annotations) (download)
Wed Sep 26 18:09:15 2001 UTC (22 years, 8 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint43a-release1mods, release1_b1, checkpoint43, icebear5, icebear4, icebear3, icebear2, release1-branch_tutorials, chkpt44a_post, chkpt44c_pre, ecco_c44_e19, ecco_c44_e18, ecco_c44_e17, ecco_c44_e16, release1-branch-end, checkpoint44b_post, ecco_ice2, ecco_ice1, ecco_c44_e22, ecco_c44_e25, chkpt44a_pre, ecco_c44_e23, ecco_c44_e20, ecco_c44_e21, ecco_c44_e26, ecco_c44_e27, ecco_c44_e24, ecco-branch-mod1, ecco-branch-mod2, ecco-branch-mod3, ecco-branch-mod4, ecco-branch-mod5, release1_beta1, checkpoint44b_pre, checkpoint42, checkpoint41, checkpoint44, chkpt44c_post, release1-branch_branchpoint
Branch point for: c24_e25_ice, release1-branch, release1, ecco-branch, icebear, release1_coupled
Changes since 1.36: +19 -9 lines
Bringing comments up to data and formatting for document extraction.

1 cnh 1.37 C $Header: /u/gcmpack/models/MITgcmUV/model/src/ini_cg2d.F,v 1.36 2001/07/06 21:39:37 jmc Exp $
2 jmc 1.30 C $Name: $
3 cnh 1.1
4 adcroft 1.13 #include "CPP_OPTIONS.h"
5 cnh 1.1
6 cnh 1.37 CBOP
7     C !ROUTINE: INI_CG2D
8     C !INTERFACE:
9 cnh 1.1 SUBROUTINE INI_CG2D( myThid )
10 cnh 1.37
11     C !DESCRIPTION: \bv
12     C *==========================================================*
13     C | SUBROUTINE INI_CG2D
14     C | o Initialise 2d conjugate gradient solver operators.
15     C *==========================================================*
16     C | These arrays are purely a function of the basin geom.
17     C | We set then here once and them use then repeatedly.
18     C *==========================================================*
19     C \ev
20    
21     C !USES:
22 adcroft 1.23 IMPLICIT NONE
23 cnh 1.1 C === Global variables ===
24     #include "SIZE.h"
25     #include "EEPARAMS.h"
26     #include "PARAMS.h"
27     #include "GRID.h"
28 adcroft 1.34 c#include "DYNVARS.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     C xG, yG - Global coordinate location.
43     C zG
44     C iG, jG - Global coordinate index
45     C bi,bj - Loop counters
46     C faceArea - Temporary used to hold cell face areas.
47     C I,J,K
48 adcroft 1.34 C myNorm - Work variable used in calculating normalisation factor
49     C sumArea - Work variable used to compute the total Domain Area
50 cnh 1.1 CHARACTER*(MAX_LEN_MBUF) msgBuf
51     INTEGER bi, bj
52     INTEGER I, J, K
53 cnh 1.7 _RL faceArea
54 cnh 1.15 _RS myNorm
55 adcroft 1.34 _RL sumArea
56 cnh 1.4 _RL aC, aCw, aCs
57 cnh 1.37 CEOP
58 jmc 1.31
59 cnh 1.1 C-- Initialise laplace operator
60     C aW2d: integral in Z Ax/dX
61     C aS2d: integral in Z Ay/dY
62     myNorm = 0. _d 0
63     DO bj=myByLo(myThid),myByHi(myThid)
64     DO bi=myBxLo(myThid),myBxHi(myThid)
65     DO J=1,sNy
66     DO I=1,sNx
67     aW2d(I,J,bi,bj) = 0. _d 0
68     aS2d(I,J,bi,bj) = 0. _d 0
69     ENDDO
70     ENDDO
71 cnh 1.17 DO K=1,Nr
72 cnh 1.1 DO J=1,sNy
73     DO I=1,sNx
74 cnh 1.20 faceArea = _dyG(I,J,bi,bj)*drF(K)
75     & *_hFacW(I,J,K,bi,bj)
76 cnh 1.1 aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj) +
77 jmc 1.32 & implicSurfPress*implicDiv2DFlow
78     & *faceArea*recip_dxC(I,J,bi,bj)
79 cnh 1.20 faceArea = _dxG(I,J,bi,bj)*drF(K)
80     & *_hFacS(I,J,K,bi,bj)
81 cnh 1.1 aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj) +
82 jmc 1.32 & implicSurfPress*implicDiv2DFlow
83     & *faceArea*recip_dyC(I,J,bi,bj)
84 cnh 1.1 ENDDO
85     ENDDO
86     ENDDO
87 adcroft 1.26 #ifdef ALLOW_OBCS
88 adcroft 1.28 IF (useOBCS) THEN
89 adcroft 1.22 DO I=1,sNx
90     IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj),bi,bj)=0.
91     IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj)+1,bi,bj)=0.
92     IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj)+1,bi,bj)=0.
93     IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj),bi,bj)=0.
94     ENDDO
95     DO J=1,sNy
96     IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj),J,bi,bj)=0.
97     IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj)+1,J,bi,bj)=0.
98     IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj)+1,J,bi,bj)=0.
99     IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj),J,bi,bj)=0.
100     ENDDO
101     ENDIF
102 adcroft 1.26 #endif
103 cnh 1.1 DO J=1,sNy
104     DO I=1,sNx
105     myNorm = MAX(ABS(aW2d(I,J,bi,bj)),myNorm)
106     myNorm = MAX(ABS(aS2d(I,J,bi,bj)),myNorm)
107     ENDDO
108     ENDDO
109     ENDDO
110     ENDDO
111 adcroft 1.25 _GLOBAL_MAX_R4( myNorm, myThid )
112     IF ( myNorm .NE. 0. _d 0 ) THEN
113     myNorm = 1. _d 0/myNorm
114 cnh 1.1 ELSE
115     myNorm = 1. _d 0
116     ENDIF
117 cnh 1.5 cg2dNorm = myNorm
118 cnh 1.1 _BEGIN_MASTER( myThid )
119     CcnhDebugStarts
120 cnh 1.20 WRITE(msgBuf,'(A,E40.25)') '// CG2D normalisation factor = ',
121     & cg2dNorm
122 cnh 1.3 CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
123     WRITE(msgBuf,*) ' '
124 cnh 1.1 CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1)
125     CcnhDebugEnds
126     _END_MASTER( myThid )
127     DO bj=myByLo(myThid),myByHi(myThid)
128     DO bi=myBxLo(myThid),myBxHi(myThid)
129     DO J=1,sNy
130     DO I=1,sNx
131     aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj)*myNorm
132     aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj)*myNorm
133     ENDDO
134     ENDDO
135     ENDDO
136     ENDDO
137    
138     C-- Update overlap regions
139     CcnhDebugStarts
140 cnh 1.14 C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.1' , 1, myThid )
141     C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.1' , 1, myThid )
142 cnh 1.1 CcnhDebugEnds
143 adcroft 1.34 c _EXCH_XY_R4(aW2d, myThid)
144     c _EXCH_XY_R4(aS2d, myThid)
145     CALL EXCH_UV_XY_RS(aW2d,aS2d,.FALSE.,myThid)
146 cnh 1.1 CcnhDebugStarts
147 adcroft 1.24 C CALL PLOT_FIELD_XYRS( aW2d, 'AW2D INI_CG2D.2' , 1, myThid )
148     C CALL PLOT_FIELD_XYRS( aS2d, 'AS2D INI_CG2D.2' , 1, myThid )
149 cnh 1.1 CcnhDebugEnds
150    
151 adcroft 1.34 C-- Define the solver tolerance in the appropriate Unit :
152     cg2dNormaliseRHS = cg2dTargetResWunit.LE.0
153     IF (cg2dNormaliseRHS) THEN
154     C- when using a normalisation of RHS, tolerance has no unit => no conversion
155     cg2dTolerance = cg2dTargetResidual
156     ELSE
157     C- compute the total Area of the domain :
158     sumArea = 0.
159     DO bj=myByLo(myThid),myByHi(myThid)
160     DO bi=myBxLo(myThid),myBxHi(myThid)
161     DO j=1,sNy
162     DO i=1,sNx
163     IF (Ro_surf(i,j,bi,bj).GT.R_low(i,j,bi,bj)) THEN
164     sumArea = sumArea + rA(i,j,bi,bj)
165     ENDIF
166     ENDDO
167     ENDDO
168     ENDDO
169     ENDDO
170 adcroft 1.35 c WRITE(*,*) ' mythid, sumArea = ', mythid, sumArea
171 adcroft 1.34 _GLOBAL_SUM_R8( sumArea, myThid )
172     C- convert Target-Residual (in W unit) to cg2d-solver residual unit [m^2/s^2]
173     cg2dTolerance = cg2dNorm * cg2dTargetResWunit
174     & * sumArea / deltaTMom
175 adcroft 1.35 WRITE(*,'(2A,1P2E22.14)') ' ini_cg2d: ',
176 adcroft 1.34 & 'sumArea,cg2dTolerance =', sumArea,cg2dTolerance
177     ENDIF
178    
179 cnh 1.1 C-- Initialise preconditioner
180 cnh 1.4 C Note. 20th May 1998
181     C I made a weird discovery! In the model paper we argue
182     C for the form of the preconditioner used here ( see
183     C A Finite-volume, Incompressible Navier-Stokes Model
184     C ...., Marshall et. al ). The algebra gives a simple
185     C 0.5 factor for the averaging of ac and aCw to get a
186     C symmettric pre-conditioner. By using a factor of 0.51
187     C i.e. scaling the off-diagonal terms in the
188     C preconditioner down slightly I managed to get the
189     C number of iterations for convergence in a test case to
190     C drop form 192 -> 134! Need to investigate this further!
191     C For now I have introduced a parameter cg2dpcOffDFac which
192     C defaults to 0.51 but can be set at runtime.
193 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
194     DO bi=myBxLo(myThid),myBxHi(myThid)
195     DO J=1,sNy
196     DO I=1,sNx
197     pC(I,J,bi,bj) = 1. _d 0
198 cnh 1.4 aC = -(
199     & aW2d(I,J,bi,bj) + aW2d(I+1,J ,bi,bj)
200 adcroft 1.24 & +aS2d(I,J,bi,bj) + aS2d(I ,J+1,bi,bj)
201 jmc 1.32 & +freeSurfFac*myNorm*recip_Bo(I,J,bi,bj)*
202 cnh 1.17 & rA(I,J,bi,bj)/deltaTMom/deltaTMom
203 cnh 1.4 & )
204     aCs = -(
205     & aW2d(I,J-1,bi,bj) + aW2d(I+1,J-1,bi,bj)
206     & +aS2d(I,J-1,bi,bj) + aS2d(I ,J ,bi,bj)
207 jmc 1.32 & +freeSurfFac*myNorm*recip_Bo(I,J-1,bi,bj)*
208 cnh 1.17 & rA(I,J-1,bi,bj)/deltaTMom/deltaTMom
209 cnh 1.4 & )
210     aCw = -(
211     & aW2d(I-1,J,bi,bj) + aW2d(I ,J ,bi,bj)
212     & +aS2d(I-1,J,bi,bj) + aS2d(I-1,J+1,bi,bj)
213 jmc 1.32 & +freeSurfFac*myNorm*recip_Bo(I-1,J,bi,bj)*
214 cnh 1.17 & rA(I-1,J,bi,bj)/deltaTMom/deltaTMom
215 cnh 1.4 & )
216     IF ( aC .EQ. 0. ) THEN
217 adcroft 1.27 pC(I,J,bi,bj) = 1. _d 0
218 cnh 1.4 ELSE
219     pC(I,J,bi,bj) = 1. _d 0 / aC
220     ENDIF
221     IF ( aC + aCw .EQ. 0. ) THEN
222     pW(I,J,bi,bj) = 0.
223     ELSE
224     pW(I,J,bi,bj) =
225     & -aW2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCw+aC))**2 )
226     ENDIF
227     IF ( aC + aCs .EQ. 0. ) THEN
228     pS(I,J,bi,bj) = 0.
229     ELSE
230     pS(I,J,bi,bj) =
231     & -aS2d(I ,J ,bi,bj)/((cg2dpcOffDFac *(aCs+aC))**2 )
232     ENDIF
233 cnh 1.6 C pC(I,J,bi,bj) = 1.
234     C pW(I,J,bi,bj) = 0.
235     C pS(I,J,bi,bj) = 0.
236 cnh 1.1 ENDDO
237     ENDDO
238     ENDDO
239     ENDDO
240     C-- Update overlap regions
241     _EXCH_XY_R4(pC, myThid)
242 adcroft 1.34 c _EXCH_XY_R4(pW, myThid)
243     c _EXCH_XY_R4(pS, myThid)
244     CALL EXCH_UV_XY_RS(pW,pS,.FALSE.,myThid)
245 cnh 1.18 CcnhDebugStarts
246 adcroft 1.24 C CALL PLOT_FIELD_XYRS( pC, 'pC INI_CG2D.2' , 1, myThid )
247     C CALL PLOT_FIELD_XYRS( pW, 'pW INI_CG2D.2' , 1, myThid )
248     C CALL PLOT_FIELD_XYRS( pS, 'pS INI_CG2D.2' , 1, myThid )
249 cnh 1.18 CcnhDebugEnds
250 cnh 1.1
251     RETURN
252     END

  ViewVC Help
Powered by ViewVC 1.1.22