1 |
C $Header: /u/gcmpack/MITgcm/model/src/update_cg2d.F,v 1.7 2006/12/05 05:25:08 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 |
|
33 |
C !INPUT/OUTPUT PARAMETERS: |
34 |
C === Routine arguments === |
35 |
C myTime :: Current time in simulation |
36 |
C myIter :: Current iteration number in simulation |
37 |
C myThid :: Thread number for this instance of the routine. |
38 |
_RL myTime |
39 |
INTEGER myIter |
40 |
INTEGER myThid |
41 |
|
42 |
C !LOCAL VARIABLES: |
43 |
#ifdef NONLIN_FRSURF |
44 |
C-- Note : compared to "INI_CG2D", no needs to compute again |
45 |
C the solver norn=malisation factor of the solver tolerance |
46 |
C === Local variables === |
47 |
C bi,bj :: tile indices |
48 |
C i,j,k :: Loop counters |
49 |
C faceArea :: Temporary used to hold cell face areas. |
50 |
INTEGER bi, bj |
51 |
INTEGER i, j, k, ks |
52 |
_RL faceArea |
53 |
_RL pW_tmp, pS_tmp |
54 |
LOGICAL updatePreCond |
55 |
CEOP |
56 |
|
57 |
C-- Decide when to update cg2d Preconditioner : |
58 |
IF ( cg2dPreCondFreq.EQ.0 ) THEN |
59 |
updatePreCond = .FALSE. |
60 |
ELSE |
61 |
updatePreCond = ( myIter.EQ.nIter0 ) |
62 |
IF ( MOD(myIter,cg2dPreCondFreq).EQ.0 ) updatePreCond=.TRUE. |
63 |
ENDIF |
64 |
|
65 |
C-- Initialise laplace operator |
66 |
C aW2d: integral in Z Ax/dX |
67 |
C aS2d: integral in Z Ay/dY |
68 |
DO bj=myByLo(myThid),myByHi(myThid) |
69 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
70 |
DO j=1,sNy+1 |
71 |
DO i=1,sNx+1 |
72 |
aW2d(i,j,bi,bj) = 0. _d 0 |
73 |
aS2d(i,j,bi,bj) = 0. _d 0 |
74 |
ENDDO |
75 |
ENDDO |
76 |
DO k=1,Nr |
77 |
DO j=1,sNy+1 |
78 |
DO i=1,sNx+1 |
79 |
C deep-model: *deepFacC (faceArea), /deepFacC (recip_dx,y): => no net effect |
80 |
faceArea = _dyG(i,j,bi,bj)*drF(k) |
81 |
& *_hFacW(i,j,k,bi,bj) |
82 |
aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj) |
83 |
& + faceArea*recip_dxC(i,j,bi,bj) |
84 |
faceArea = _dxG(i,j,bi,bj)*drF(k) |
85 |
& *_hFacS(i,j,k,bi,bj) |
86 |
aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj) |
87 |
& + faceArea*recip_dyC(i,j,bi,bj) |
88 |
ENDDO |
89 |
ENDDO |
90 |
ENDDO |
91 |
DO j=1,sNy+1 |
92 |
DO i=1,sNx+1 |
93 |
aW2d(i,j,bi,bj) = aW2d(i,j,bi,bj)*cg2dNorm |
94 |
& *implicSurfPress*implicDiv2DFlow |
95 |
#ifdef ALLOW_OBCS |
96 |
& *maskInC(i,j,bi,bj)*maskInC(i-1,j,bi,bj) |
97 |
#endif |
98 |
aS2d(i,j,bi,bj) = aS2d(i,j,bi,bj)*cg2dNorm |
99 |
& *implicSurfPress*implicDiv2DFlow |
100 |
#ifdef ALLOW_OBCS |
101 |
& *maskInC(i,j,bi,bj)*maskInC(i,j-1,bi,bj) |
102 |
#endif |
103 |
ENDDO |
104 |
ENDDO |
105 |
C-- compute matrix main diagonal : |
106 |
IF ( deepAtmosphere ) THEN |
107 |
DO j=1,sNy |
108 |
DO i=1,sNx |
109 |
ks = ksurfC(i,j,bi,bj) |
110 |
aC2d(i,j,bi,bj) = -( |
111 |
& aW2d(i,j,bi,bj) + aW2d(i+1, j ,bi,bj) |
112 |
& +aS2d(i,j,bi,bj) + aS2d( i ,j+1,bi,bj) |
113 |
& +freeSurfFac*cg2dNorm*recip_Bo(i,j,bi,bj)*deepFac2F(ks) |
114 |
& *rA(i,j,bi,bj)/deltaTMom/deltaTfreesurf |
115 |
& ) |
116 |
ENDDO |
117 |
ENDDO |
118 |
ELSE |
119 |
DO j=1,sNy |
120 |
DO i=1,sNx |
121 |
aC2d(i,j,bi,bj) = -( |
122 |
& aW2d(i,j,bi,bj) + aW2d(i+1, j ,bi,bj) |
123 |
& +aS2d(i,j,bi,bj) + aS2d( i ,j+1,bi,bj) |
124 |
& +freeSurfFac*cg2dNorm*recip_Bo(i,j,bi,bj) |
125 |
& *rA(i,j,bi,bj)/deltaTMom/deltaTfreesurf |
126 |
& ) |
127 |
ENDDO |
128 |
ENDDO |
129 |
ENDIF |
130 |
C- end bi,bj loops |
131 |
ENDDO |
132 |
ENDDO |
133 |
|
134 |
IF ( updatePreCond ) THEN |
135 |
C-- Update overlap regions |
136 |
CALL EXCH_XY_RS(aC2d, myThid) |
137 |
|
138 |
C-- Initialise preconditioner |
139 |
DO bj=myByLo(myThid),myByHi(myThid) |
140 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
141 |
DO j=1,sNy+1 |
142 |
DO i=1,sNx+1 |
143 |
IF ( aC2d(i,j,bi,bj) .EQ. 0. ) THEN |
144 |
pC(i,j,bi,bj) = 1. _d 0 |
145 |
ELSE |
146 |
pC(i,j,bi,bj) = 1. _d 0 / aC2d(i,j,bi,bj) |
147 |
ENDIF |
148 |
pW_tmp = aC2d(i,j,bi,bj)+aC2d(i-1,j,bi,bj) |
149 |
IF ( pW_tmp .EQ. 0. ) THEN |
150 |
pW(i,j,bi,bj) = 0. |
151 |
ELSE |
152 |
pW(i,j,bi,bj) = |
153 |
& -aW2d(i,j,bi,bj)/((cg2dpcOffDFac *pW_tmp)**2 ) |
154 |
ENDIF |
155 |
pS_tmp = aC2d(i,j,bi,bj)+aC2d(i,j-1,bi,bj) |
156 |
IF ( pS_tmp .EQ. 0. ) THEN |
157 |
pS(i,j,bi,bj) = 0. |
158 |
ELSE |
159 |
pS(i,j,bi,bj) = |
160 |
& -aS2d(i,j,bi,bj)/((cg2dpcOffDFac *pS_tmp)**2 ) |
161 |
ENDIF |
162 |
ENDDO |
163 |
ENDDO |
164 |
ENDDO |
165 |
ENDDO |
166 |
C- if update Preconditioner : end |
167 |
ENDIF |
168 |
|
169 |
#endif /* NONLIN_FRSURF */ |
170 |
|
171 |
RETURN |
172 |
END |