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