1 |
jmc |
1.26 |
C $Header: /u/gcmpack/MITgcm/model/src/calc_div_ghat.F,v 1.25 2009/09/27 23:18:53 jmc Exp $ |
2 |
jmc |
1.12 |
C $Name: $ |
3 |
cnh |
1.1 |
|
4 |
cnh |
1.6 |
#include "CPP_OPTIONS.h" |
5 |
edhill |
1.19 |
|
6 |
cnh |
1.17 |
CBOP |
7 |
|
|
C !ROUTINE: CALC_DIV_GHAT |
8 |
|
|
C !INTERFACE: |
9 |
|
|
SUBROUTINE CALC_DIV_GHAT( |
10 |
jmc |
1.22 |
I bi,bj,k, |
11 |
|
|
U cg2d_b, cg3d_b, |
12 |
|
|
I myThid) |
13 |
cnh |
1.17 |
C !DESCRIPTION: \bv |
14 |
|
|
C *==========================================================* |
15 |
jmc |
1.21 |
C | S/R CALC_DIV_GHAT |
16 |
|
|
C | o Form the right hand-side of the surface pressure eqn. |
17 |
cnh |
1.17 |
C *==========================================================* |
18 |
|
|
C | Right hand side of pressure equation is divergence |
19 |
|
|
C | of veclocity tendency (GHAT) term along with a relaxation |
20 |
|
|
C | term equal to the barotropic flow field divergence. |
21 |
|
|
C *==========================================================* |
22 |
|
|
C \ev |
23 |
cnh |
1.1 |
|
24 |
cnh |
1.17 |
C !USES: |
25 |
cnh |
1.1 |
IMPLICIT NONE |
26 |
|
|
C == Global variables == |
27 |
|
|
#include "SIZE.h" |
28 |
|
|
#include "EEPARAMS.h" |
29 |
|
|
#include "PARAMS.h" |
30 |
|
|
#include "GRID.h" |
31 |
jmc |
1.22 |
#include "DYNVARS.h" |
32 |
cnh |
1.1 |
|
33 |
cnh |
1.17 |
C !INPUT/OUTPUT PARAMETERS: |
34 |
cnh |
1.1 |
C == Routine arguments == |
35 |
jmc |
1.21 |
C bi, bj :: tile indices |
36 |
|
|
C k :: Index of layer. |
37 |
|
|
C cg2d_b :: Conjugate Gradient 2-D solver : Right-hand side vector |
38 |
jmc |
1.22 |
C cg3d_b :: Conjugate Gradient 3-D solver : Right-hand side vector |
39 |
jmc |
1.21 |
C myThid :: Instance number for this call of CALC_DIV_GHAT |
40 |
jmc |
1.22 |
INTEGER bi,bj |
41 |
jmc |
1.21 |
INTEGER k |
42 |
jmc |
1.13 |
_RL cg2d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
43 |
jmc |
1.24 |
#ifdef ALLOW_NONHYDROSTATIC |
44 |
jmc |
1.22 |
_RL cg3d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
45 |
jmc |
1.24 |
#else |
46 |
|
|
_RL cg3d_b(1) |
47 |
|
|
#endif |
48 |
cnh |
1.1 |
INTEGER myThid |
49 |
|
|
|
50 |
cnh |
1.17 |
C !LOCAL VARIABLES: |
51 |
cnh |
1.1 |
C == Local variables == |
52 |
jmc |
1.22 |
C i,j :: Loop counters |
53 |
|
|
C xA, yA :: Cell vertical face areas |
54 |
|
|
C pf :: Intermediate array for building RHS source term. |
55 |
cnh |
1.1 |
INTEGER i,j |
56 |
jmc |
1.22 |
_RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
57 |
|
|
_RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
58 |
cnh |
1.1 |
_RL pf (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
59 |
cnh |
1.17 |
CEOP |
60 |
cnh |
1.1 |
|
61 |
jmc |
1.22 |
C Calculate vertical face areas |
62 |
|
|
DO j=1,sNy+1 |
63 |
|
|
DO i=1,sNx+1 |
64 |
|
|
xA(i,j) = _dyG(i,j,bi,bj)*deepFacC(k) |
65 |
|
|
& *drF(k)*_hFacW(i,j,k,bi,bj)*rhoFacC(k) |
66 |
|
|
yA(i,j) = _dxG(i,j,bi,bj)*deepFacC(k) |
67 |
|
|
& *drF(k)*_hFacS(i,j,k,bi,bj)*rhoFacC(k) |
68 |
|
|
ENDDO |
69 |
|
|
ENDDO |
70 |
|
|
|
71 |
cnh |
1.1 |
C-- Pressure equation source term |
72 |
jmc |
1.21 |
C Term is the vertical integral of the divergence of the |
73 |
cnh |
1.1 |
C time tendency terms along with a relaxation term that |
74 |
|
|
C pulls div(U) + dh/dt back toward zero. |
75 |
|
|
|
76 |
jmc |
1.16 |
IF (implicDiv2Dflow.EQ.1.) THEN |
77 |
jmc |
1.12 |
C Fully Implicit treatment of the Barotropic Flow Divergence |
78 |
|
|
DO j=1,sNy |
79 |
|
|
DO i=1,sNx+1 |
80 |
jmc |
1.18 |
pf(i,j) = xA(i,j)*gU(i,j,k,bi,bj) / deltaTmom |
81 |
jmc |
1.12 |
ENDDO |
82 |
|
|
ENDDO |
83 |
jmc |
1.16 |
ELSEIF (exactConserv) THEN |
84 |
|
|
c ELSEIF (nonlinFreeSurf.GT.0) THEN |
85 |
|
|
C Implicit treatment of the Barotropic Flow Divergence |
86 |
|
|
DO j=1,sNy |
87 |
|
|
DO i=1,sNx+1 |
88 |
jmc |
1.21 |
pf(i,j) = implicDiv2Dflow |
89 |
jmc |
1.18 |
& *xA(i,j)*gU(i,j,k,bi,bj) / deltaTmom |
90 |
jmc |
1.16 |
ENDDO |
91 |
|
|
ENDDO |
92 |
jmc |
1.12 |
ELSE |
93 |
|
|
C Explicit+Implicit part of the Barotropic Flow Divergence |
94 |
|
|
C => Filtering of uVel,vVel is necessary |
95 |
jmc |
1.16 |
C-- Now the filter are applied in the_correction_step(). |
96 |
|
|
C We have left this code here to indicate where the filters used to be |
97 |
jmc |
1.21 |
C in the algorithm before JMC moved them to after the pressure solver. |
98 |
jmc |
1.25 |
c#ifdef ALLOW_ZONAL_FILT |
99 |
|
|
c IF (zonal_filt_lat.LT.90.) THEN |
100 |
|
|
c CALL ZONAL_FILTER( |
101 |
|
|
c U uVel( 1-OLx,1-OLy,k,bi,bj), |
102 |
|
|
c I hFacW(1-OLx,1-OLy,k,bi,bj), |
103 |
|
|
c I 0, sNy+1, 1, bi, bj, 1, myThid ) |
104 |
|
|
c CALL ZONAL_FILTER( |
105 |
|
|
c U vVel( 1-OLx,1-OLy,k,bi,bj), |
106 |
|
|
c I hFacS(1-OLx,1-OLy,k,bi,bj), |
107 |
|
|
c I 0, sNy+1, 1, bi, bj, 2, myThid ) |
108 |
|
|
c ENDIF |
109 |
|
|
c#endif |
110 |
jmc |
1.12 |
DO j=1,sNy |
111 |
|
|
DO i=1,sNx+1 |
112 |
jmc |
1.18 |
pf(i,j) = ( implicDiv2Dflow * gU(i,j,k,bi,bj) |
113 |
jmc |
1.22 |
& + (1. _d 0-implicDiv2Dflow)* uVel(i,j,k,bi,bj) |
114 |
jmc |
1.12 |
& ) * xA(i,j) / deltaTmom |
115 |
|
|
ENDDO |
116 |
|
|
ENDDO |
117 |
|
|
ENDIF |
118 |
cnh |
1.1 |
DO j=1,sNy |
119 |
|
|
DO i=1,sNx |
120 |
adcroft |
1.8 |
cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj) + |
121 |
cnh |
1.1 |
& pf(i+1,j)-pf(i,j) |
122 |
|
|
ENDDO |
123 |
|
|
ENDDO |
124 |
|
|
|
125 |
adcroft |
1.8 |
#ifdef ALLOW_NONHYDROSTATIC |
126 |
jmc |
1.21 |
IF (use3Dsolver) THEN |
127 |
adcroft |
1.8 |
DO j=1,sNy |
128 |
|
|
DO i=1,sNx |
129 |
jmc |
1.26 |
cg3d_b(i,j,k,bi,bj) = ( pf(i+1,j)-pf(i,j) ) |
130 |
adcroft |
1.8 |
ENDDO |
131 |
|
|
ENDDO |
132 |
|
|
ENDIF |
133 |
|
|
#endif |
134 |
|
|
|
135 |
jmc |
1.16 |
IF (implicDiv2Dflow.EQ.1.) THEN |
136 |
jmc |
1.12 |
C Fully Implicit treatment of the Barotropic Flow Divergence |
137 |
|
|
DO j=1,sNy+1 |
138 |
|
|
DO i=1,sNx |
139 |
jmc |
1.18 |
pf(i,j) = yA(i,j)*gV(i,j,k,bi,bj) / deltatmom |
140 |
jmc |
1.16 |
ENDDO |
141 |
|
|
ENDDO |
142 |
|
|
ELSEIF (exactConserv) THEN |
143 |
|
|
c ELSEIF (nonlinFreeSurf.GT.0) THEN |
144 |
|
|
C Implicit treatment of the Barotropic Flow Divergence |
145 |
|
|
DO j=1,sNy+1 |
146 |
|
|
DO i=1,sNx |
147 |
|
|
pf(i,j) = implicDiv2Dflow |
148 |
jmc |
1.18 |
& *yA(i,j)*gV(i,j,k,bi,bj) / deltatmom |
149 |
jmc |
1.12 |
ENDDO |
150 |
|
|
ENDDO |
151 |
|
|
ELSE |
152 |
|
|
C Explicit+Implicit part of the Barotropic Flow Divergence |
153 |
|
|
DO j=1,sNy+1 |
154 |
|
|
DO i=1,sNx |
155 |
jmc |
1.18 |
pf(i,j) = ( implicDiv2Dflow * gV(i,j,k,bi,bj) |
156 |
jmc |
1.22 |
& + (1. _d 0-implicDiv2Dflow)* vVel(i,j,k,bi,bj) |
157 |
jmc |
1.12 |
& ) * yA(i,j) / deltaTmom |
158 |
|
|
ENDDO |
159 |
|
|
ENDDO |
160 |
|
|
ENDIF |
161 |
cnh |
1.1 |
DO j=1,sNy |
162 |
|
|
DO i=1,sNx |
163 |
adcroft |
1.8 |
cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj) + |
164 |
cnh |
1.1 |
& pf(i,j+1)-pf(i,j) |
165 |
|
|
ENDDO |
166 |
|
|
ENDDO |
167 |
cnh |
1.4 |
|
168 |
adcroft |
1.8 |
#ifdef ALLOW_NONHYDROSTATIC |
169 |
jmc |
1.21 |
IF (use3Dsolver) THEN |
170 |
adcroft |
1.8 |
DO j=1,sNy |
171 |
|
|
DO i=1,sNx |
172 |
jmc |
1.26 |
cg3d_b(i,j,k,bi,bj) = cg3d_b(i,j,k,bi,bj) |
173 |
|
|
& + ( pf(i,j+1)-pf(i,j) ) |
174 |
adcroft |
1.8 |
ENDDO |
175 |
|
|
ENDDO |
176 |
|
|
ENDIF |
177 |
|
|
#endif |
178 |
cnh |
1.1 |
|
179 |
jmc |
1.23 |
#ifdef ALLOW_ADDFLUID |
180 |
|
|
IF ( selectAddFluid.GE.1 ) THEN |
181 |
|
|
DO j=1,sNy |
182 |
|
|
DO i=1,sNx |
183 |
|
|
cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj) |
184 |
|
|
& - addMass(i,j,k,bi,bj)*mass2rUnit/deltaTmom |
185 |
|
|
ENDDO |
186 |
|
|
ENDDO |
187 |
|
|
#ifdef ALLOW_NONHYDROSTATIC |
188 |
|
|
IF (use3Dsolver) THEN |
189 |
|
|
DO j=1,sNy |
190 |
|
|
DO i=1,sNx |
191 |
|
|
cg3d_b(i,j,k,bi,bj) = cg3d_b(i,j,k,bi,bj) |
192 |
|
|
& - addMass(i,j,k,bi,bj)*mass2rUnit/deltaTmom |
193 |
|
|
ENDDO |
194 |
|
|
ENDDO |
195 |
|
|
ENDIF |
196 |
|
|
#endif |
197 |
|
|
ENDIF |
198 |
|
|
#endif /* ALLOW_ADDFLUID */ |
199 |
|
|
|
200 |
cnh |
1.1 |
RETURN |
201 |
|
|
END |