1 |
C $Header: /u/gcmpack/MITgcm/model/src/update_cg2d.F,v 1.6 2004/07/09 22:32:35 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "PACKAGES_CONFIG.h" |
5 |
#include "CPP_OPTIONS.h" |
6 |
|
7 |
CBOP |
8 |
C !ROUTINE: UPDATE_CG2D |
9 |
C !INTERFACE: |
10 |
SUBROUTINE UPDATE_CG2D( myTime, myIter, myThid ) |
11 |
C !DESCRIPTION: \bv |
12 |
C *==========================================================* |
13 |
C | SUBROUTINE UPDATE_CG2D |
14 |
C | o Update 2d conjugate gradient solver operators |
15 |
C | account for Free-Surf effect on total column thickness |
16 |
C *==========================================================* |
17 |
C | This routine is based on INI_CG2D, and simplified. It is |
18 |
C | only active when the non-linear free surface mode of |
19 |
C | equations is active. |
20 |
C *==========================================================* |
21 |
C \ev |
22 |
|
23 |
C !USES: |
24 |
IMPLICIT NONE |
25 |
C === Global variables === |
26 |
#include "SIZE.h" |
27 |
#include "EEPARAMS.h" |
28 |
#include "PARAMS.h" |
29 |
#include "GRID.h" |
30 |
#include "SURFACE.h" |
31 |
#include "CG2D.h" |
32 |
#ifdef ALLOW_OBCS |
33 |
#include "OBCS.h" |
34 |
#endif |
35 |
|
36 |
C !INPUT/OUTPUT PARAMETERS: |
37 |
C === Routine arguments === |
38 |
C myTime - Current time in simulation |
39 |
C myIter - Current iteration number in simulation |
40 |
C myThid - Thread number for this instance of the routine. |
41 |
_RL myTime |
42 |
INTEGER myIter |
43 |
INTEGER myThid |
44 |
|
45 |
C !LOCAL VARIABLES: |
46 |
#ifdef NONLIN_FRSURF |
47 |
C-- Note : compared to "INI_CG2D", no needs to compute again |
48 |
C the solver norn=malisation factor of the solver tolerance |
49 |
C === Local variables === |
50 |
C bi,bj :: tile indices |
51 |
C I,J,K :: Loop counters |
52 |
C faceArea :: Temporary used to hold cell face areas. |
53 |
INTEGER bi, bj |
54 |
INTEGER I, J, K, ks |
55 |
_RL faceArea |
56 |
_RL pW_tmp, pS_tmp |
57 |
LOGICAL updatePreCond |
58 |
CEOP |
59 |
|
60 |
C-- Decide when to update cg2d Preconditioner : |
61 |
IF ( cg2dPreCondFreq.EQ.0 ) THEN |
62 |
updatePreCond = .FALSE. |
63 |
ELSE |
64 |
updatePreCond = ( myIter.EQ.nIter0 ) |
65 |
IF ( MOD(myIter,cg2dPreCondFreq).EQ.0 ) updatePreCond=.TRUE. |
66 |
ENDIF |
67 |
|
68 |
C-- Initialise laplace operator |
69 |
C aW2d: integral in Z Ax/dX |
70 |
C aS2d: integral in Z Ay/dY |
71 |
DO bj=myByLo(myThid),myByHi(myThid) |
72 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
73 |
DO J=1,sNy+1 |
74 |
DO I=1,sNx+1 |
75 |
aW2d(I,J,bi,bj) = 0. _d 0 |
76 |
aS2d(I,J,bi,bj) = 0. _d 0 |
77 |
ENDDO |
78 |
ENDDO |
79 |
DO K=1,Nr |
80 |
DO J=1,sNy+1 |
81 |
DO I=1,sNx+1 |
82 |
C deep-model: *deepFacC (faceArea), /deepFacC (recip_dx,y): => no net effect |
83 |
faceArea = _dyG(I,J,bi,bj)*drF(K) |
84 |
& *_hFacW(I,J,K,bi,bj) |
85 |
aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj) |
86 |
& + faceArea*recip_dxC(I,J,bi,bj) |
87 |
faceArea = _dxG(I,J,bi,bj)*drF(K) |
88 |
& *_hFacS(I,J,K,bi,bj) |
89 |
aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj) |
90 |
& + faceArea*recip_dyC(I,J,bi,bj) |
91 |
ENDDO |
92 |
ENDDO |
93 |
ENDDO |
94 |
#ifdef ALLOW_OBCS |
95 |
IF (useOBCS) THEN |
96 |
C Note: would need loop from 1 to sNx+1 (and below, from 1 to sNy+1) |
97 |
C to get the same solver-matrix as from INI_CG2D, |
98 |
C since aS2d & aW2d are not exchanged here (but in INI_CG2D, they are). |
99 |
DO I=1,sNx |
100 |
IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj),bi,bj)=0. |
101 |
IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj)+1,bi,bj)=0. |
102 |
IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj)+1,bi,bj)=0. |
103 |
IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj),bi,bj)=0. |
104 |
ENDDO |
105 |
DO J=1,sNy |
106 |
IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj),J,bi,bj)=0. |
107 |
IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj)+1,J,bi,bj)=0. |
108 |
IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj)+1,J,bi,bj)=0. |
109 |
IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj),J,bi,bj)=0. |
110 |
ENDDO |
111 |
ENDIF |
112 |
#endif |
113 |
DO J=1,sNy+1 |
114 |
DO I=1,sNx+1 |
115 |
aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj)*cg2dNorm |
116 |
& *implicSurfPress*implicDiv2DFlow |
117 |
aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj)*cg2dNorm |
118 |
& *implicSurfPress*implicDiv2DFlow |
119 |
ENDDO |
120 |
ENDDO |
121 |
C-- compute matrix main diagonal : |
122 |
IF ( deepAtmosphere ) THEN |
123 |
DO J=1,sNy |
124 |
DO I=1,sNx |
125 |
ks = ksurfC(I,J,bi,bj) |
126 |
aC2d(I,J,bi,bj) = -( |
127 |
& aW2d(I,J,bi,bj) + aW2d(I+1,J ,bi,bj) |
128 |
& +aS2d(I,J,bi,bj) + aS2d(I ,J+1,bi,bj) |
129 |
& +freeSurfFac*cg2dNorm*recip_Bo(I,J,bi,bj)*deepFac2F(ks) |
130 |
& *rA(I,J,bi,bj)/deltaTMom/deltaTfreesurf |
131 |
& ) |
132 |
ENDDO |
133 |
ENDDO |
134 |
ELSE |
135 |
DO J=1,sNy |
136 |
DO I=1,sNx |
137 |
aC2d(I,J,bi,bj) = -( |
138 |
& aW2d(I,J,bi,bj) + aW2d(I+1,J ,bi,bj) |
139 |
& +aS2d(I,J,bi,bj) + aS2d(I ,J+1,bi,bj) |
140 |
& +freeSurfFac*cg2dNorm*recip_Bo(I,J,bi,bj) |
141 |
& *rA(I,J,bi,bj)/deltaTMom/deltaTfreesurf |
142 |
& ) |
143 |
ENDDO |
144 |
ENDDO |
145 |
ENDIF |
146 |
C- end bi,bj loops |
147 |
ENDDO |
148 |
ENDDO |
149 |
|
150 |
IF ( updatePreCond ) THEN |
151 |
C-- Update overlap regions |
152 |
CALL EXCH_XY_RS(aC2d, myThid) |
153 |
|
154 |
C-- Initialise preconditioner |
155 |
DO bj=myByLo(myThid),myByHi(myThid) |
156 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
157 |
DO J=1,sNy+1 |
158 |
DO I=1,sNx+1 |
159 |
IF ( aC2d(I,J,bi,bj) .EQ. 0. ) THEN |
160 |
pC(I,J,bi,bj) = 1. _d 0 |
161 |
ELSE |
162 |
pC(I,J,bi,bj) = 1. _d 0 / aC2d(I,J,bi,bj) |
163 |
ENDIF |
164 |
pW_tmp = aC2d(I,J,bi,bj)+aC2d(I-1,J,bi,bj) |
165 |
IF ( pW_tmp .EQ. 0. ) THEN |
166 |
pW(I,J,bi,bj) = 0. |
167 |
ELSE |
168 |
pW(I,J,bi,bj) = |
169 |
& -aW2d(I,J,bi,bj)/((cg2dpcOffDFac *pW_tmp)**2 ) |
170 |
ENDIF |
171 |
pS_tmp = aC2d(I,J,bi,bj)+aC2d(I,J-1,bi,bj) |
172 |
IF ( pS_tmp .EQ. 0. ) THEN |
173 |
pS(I,J,bi,bj) = 0. |
174 |
ELSE |
175 |
pS(I,J,bi,bj) = |
176 |
& -aS2d(I,J,bi,bj)/((cg2dpcOffDFac *pS_tmp)**2 ) |
177 |
ENDIF |
178 |
ENDDO |
179 |
ENDDO |
180 |
ENDDO |
181 |
ENDDO |
182 |
C- if update Preconditioner : end |
183 |
ENDIF |
184 |
|
185 |
#endif /* NONLIN_FRSURF */ |
186 |
|
187 |
RETURN |
188 |
END |