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