1 |
jmc |
1.3 |
C $Header: /home/ubuntu/mnt/e9_copy/MITgcm/model/src/update_cg2d.F,v 1.3 2002/02/10 20:04:51 jmc Exp $ |
2 |
jmc |
1.1 |
C $Name: $ |
3 |
|
|
|
4 |
|
|
#include "CPP_OPTIONS.h" |
5 |
cnh |
1.2 |
CBOP |
6 |
|
|
C !ROUTINE: UPDATE_CG2D |
7 |
|
|
C !INTERFACE: |
8 |
|
|
SUBROUTINE UPDATE_CG2D( myTime, myIter, myThid ) |
9 |
|
|
C !DESCRIPTION: \bv |
10 |
|
|
C *==========================================================* |
11 |
|
|
C | SUBROUTINE UPDATE_CG2D |
12 |
|
|
C | o Update 2d conjugate gradient solver operators |
13 |
|
|
C | account for Free-Surf effect on total column thickness |
14 |
|
|
C *==========================================================* |
15 |
|
|
C | This routine is based on INI_CG2D, and simplified. It is |
16 |
|
|
C | only active when the non-linear free surface mode of |
17 |
|
|
C | equations is active. |
18 |
|
|
C *==========================================================* |
19 |
|
|
C \ev |
20 |
jmc |
1.1 |
|
21 |
cnh |
1.2 |
C !USES: |
22 |
jmc |
1.1 |
IMPLICIT NONE |
23 |
|
|
C === Global variables === |
24 |
|
|
#include "SIZE.h" |
25 |
|
|
#include "EEPARAMS.h" |
26 |
|
|
#include "PARAMS.h" |
27 |
|
|
#include "GRID.h" |
28 |
|
|
#include "SURFACE.h" |
29 |
|
|
#include "CG2D.h" |
30 |
|
|
#ifdef ALLOW_OBCS |
31 |
|
|
#include "OBCS.h" |
32 |
|
|
#endif |
33 |
|
|
|
34 |
cnh |
1.2 |
C !INPUT/OUTPUT PARAMETERS: |
35 |
jmc |
1.1 |
C === Routine arguments === |
36 |
|
|
C myTime - Current time in simulation |
37 |
|
|
C myIter - Current iteration number in simulation |
38 |
|
|
C myThid - Thread number for this instance of the routine. |
39 |
|
|
_RL myTime |
40 |
|
|
INTEGER myIter |
41 |
|
|
INTEGER myThid |
42 |
cnh |
1.2 |
|
43 |
|
|
C !LOCAL VARIABLES: |
44 |
jmc |
1.1 |
#ifdef NONLIN_FRSURF |
45 |
|
|
C-- Note : compared to "INI_CG2D", no needs to compute again |
46 |
|
|
C the solver norn=malisation factor of the solver tolerance |
47 |
|
|
C === Local variables === |
48 |
|
|
C bi,bj - Loop counters |
49 |
|
|
C I,J,K |
50 |
|
|
C faceArea - Temporary used to hold cell face areas. |
51 |
|
|
INTEGER bi, bj |
52 |
|
|
INTEGER I, J, K |
53 |
|
|
_RL faceArea |
54 |
jmc |
1.3 |
_RL pW_tmp, pS_tmp |
55 |
cnh |
1.2 |
CEOP |
56 |
jmc |
1.1 |
|
57 |
|
|
C-- Initialise laplace operator |
58 |
|
|
C aW2d: integral in Z Ax/dX |
59 |
|
|
C aS2d: integral in Z Ay/dY |
60 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
61 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
62 |
jmc |
1.3 |
DO J=1,sNy+1 |
63 |
|
|
DO I=1,sNx+1 |
64 |
jmc |
1.1 |
aW2d(I,J,bi,bj) = 0. _d 0 |
65 |
|
|
aS2d(I,J,bi,bj) = 0. _d 0 |
66 |
|
|
ENDDO |
67 |
|
|
ENDDO |
68 |
|
|
DO K=1,Nr |
69 |
jmc |
1.3 |
DO J=1,sNy+1 |
70 |
|
|
DO I=1,sNx+1 |
71 |
jmc |
1.1 |
faceArea = _dyG(I,J,bi,bj)*drF(K) |
72 |
|
|
& *_hFacW(I,J,K,bi,bj) |
73 |
|
|
aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj) + |
74 |
|
|
& faceArea*recip_dxC(I,J,bi,bj) |
75 |
|
|
faceArea = _dxG(I,J,bi,bj)*drF(K) |
76 |
|
|
& *_hFacS(I,J,K,bi,bj) |
77 |
|
|
aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj) + |
78 |
|
|
& faceArea*recip_dyC(I,J,bi,bj) |
79 |
|
|
ENDDO |
80 |
|
|
ENDDO |
81 |
|
|
ENDDO |
82 |
|
|
#ifdef ALLOW_OBCS |
83 |
|
|
IF (useOBCS) THEN |
84 |
|
|
DO I=1,sNx |
85 |
|
|
IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj),bi,bj)=0. |
86 |
|
|
IF (OB_Jn(I,bi,bj).NE.0) aS2d(I,OB_Jn(I,bi,bj)+1,bi,bj)=0. |
87 |
|
|
IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj)+1,bi,bj)=0. |
88 |
|
|
IF (OB_Js(I,bi,bj).NE.0) aS2d(I,OB_Js(I,bi,bj),bi,bj)=0. |
89 |
|
|
ENDDO |
90 |
|
|
DO J=1,sNy |
91 |
|
|
IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj),J,bi,bj)=0. |
92 |
|
|
IF (OB_Ie(J,bi,bj).NE.0) aW2d(OB_Ie(J,bi,bj)+1,J,bi,bj)=0. |
93 |
|
|
IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj)+1,J,bi,bj)=0. |
94 |
|
|
IF (OB_Iw(J,bi,bj).NE.0) aW2d(OB_Iw(J,bi,bj),J,bi,bj)=0. |
95 |
|
|
ENDDO |
96 |
|
|
ENDIF |
97 |
|
|
#endif |
98 |
jmc |
1.3 |
DO J=1,sNy+1 |
99 |
|
|
DO I=1,sNx+1 |
100 |
jmc |
1.1 |
aW2d(I,J,bi,bj) = aW2d(I,J,bi,bj)*cg2dNorm |
101 |
|
|
& *implicSurfPress*implicDiv2DFlow |
102 |
|
|
aS2d(I,J,bi,bj) = aS2d(I,J,bi,bj)*cg2dNorm |
103 |
|
|
& *implicSurfPress*implicDiv2DFlow |
104 |
|
|
ENDDO |
105 |
|
|
ENDDO |
106 |
jmc |
1.3 |
C-- Start to compute preconditioner matrix (use cg2d_q as temporary array) |
107 |
|
|
DO J=1,sNy |
108 |
|
|
DO I=1,sNx |
109 |
|
|
cg2d_q(I,J,bi,bj) = -( |
110 |
|
|
& aW2d(I,J,bi,bj) + aW2d(I+1,J ,bi,bj) |
111 |
|
|
& +aS2d(I,J,bi,bj) + aS2d(I ,J+1,bi,bj) |
112 |
|
|
& +freeSurfFac*cg2dNorm*recip_Bo(I,J,bi,bj)* |
113 |
|
|
& rA(I,J,bi,bj)/deltaTMom/deltaTMom |
114 |
|
|
& ) |
115 |
|
|
ENDDO |
116 |
|
|
ENDDO |
117 |
jmc |
1.1 |
ENDDO |
118 |
|
|
ENDDO |
119 |
|
|
|
120 |
|
|
C-- Update overlap regions |
121 |
jmc |
1.3 |
_EXCH_XY_R8(cg2d_q, myThid) |
122 |
jmc |
1.1 |
|
123 |
|
|
C-- Initialise preconditioner |
124 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
125 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
126 |
jmc |
1.3 |
DO J=1,sNy+1 |
127 |
|
|
DO I=1,sNx+1 |
128 |
|
|
IF ( cg2d_q(I,J,bi,bj) .EQ. 0. ) THEN |
129 |
jmc |
1.1 |
pC(I,J,bi,bj) = 1. _d 0 |
130 |
|
|
ELSE |
131 |
jmc |
1.3 |
pC(I,J,bi,bj) = 1. _d 0 / cg2d_q(I,J,bi,bj) |
132 |
jmc |
1.1 |
ENDIF |
133 |
jmc |
1.3 |
pW_tmp = cg2d_q(I,J,bi,bj)+cg2d_q(I-1,J,bi,bj) |
134 |
|
|
IF ( pW_tmp .EQ. 0. ) THEN |
135 |
jmc |
1.1 |
pW(I,J,bi,bj) = 0. |
136 |
|
|
ELSE |
137 |
|
|
pW(I,J,bi,bj) = |
138 |
jmc |
1.3 |
& -aW2d(I,J,bi,bj)/((cg2dpcOffDFac *pW_tmp)**2 ) |
139 |
jmc |
1.1 |
ENDIF |
140 |
jmc |
1.3 |
pS_tmp = cg2d_q(I,J,bi,bj)+cg2d_q(I,J-1,bi,bj) |
141 |
|
|
IF ( pS_tmp .EQ. 0. ) THEN |
142 |
jmc |
1.1 |
pS(I,J,bi,bj) = 0. |
143 |
|
|
ELSE |
144 |
|
|
pS(I,J,bi,bj) = |
145 |
jmc |
1.3 |
& -aS2d(I,J,bi,bj)/((cg2dpcOffDFac *pS_tmp)**2 ) |
146 |
jmc |
1.1 |
ENDIF |
147 |
|
|
ENDDO |
148 |
|
|
ENDDO |
149 |
|
|
ENDDO |
150 |
|
|
ENDDO |
151 |
|
|
|
152 |
|
|
#endif /* NONLIN_FRSURF */ |
153 |
|
|
|
154 |
|
|
RETURN |
155 |
|
|
END |