/[MITgcm]/MITgcm/model/src/calc_div_ghat.F
ViewVC logotype

Contents of /MITgcm/model/src/calc_div_ghat.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.13 - (show annotations) (download)
Tue Mar 6 16:51:02 2001 UTC (23 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint38, pre38tag1, c37_adj, pre38-close, checkpoint39, checkpoint37
Branch point for: pre38
Changes since 1.12: +6 -20 lines
separate the state variable "eta" from the 2D solver solution cg2d_x

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_div_ghat.F,v 1.12 2001/02/20 15:06:21 jmc Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 C /==========================================================\
7 C | S/R CALC_DIV_GHAT |
8 C | o Form the right hand-side of the surface pressure eqn. |
9 C |==========================================================|
10 C \==========================================================/
11 SUBROUTINE CALC_DIV_GHAT(
12 I bi,bj,iMin,iMax,jMin,jMax,K,
13 I xA,yA,
14 U cg2d_b,
15 I myThid)
16
17 IMPLICIT NONE
18
19 C == Global variables ==
20 #include "SIZE.h"
21 #include "DYNVARS.h"
22 #include "FFIELDS.h"
23 #include "EEPARAMS.h"
24 #include "PARAMS.h"
25 #include "GRID.h"
26 #ifdef ALLOW_NONHYDROSTATIC
27 #include "CG3D.h"
28 #endif
29
30 C == Routine arguments ==
31 C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
32 C results will be set.
33 C k - Index of layer.
34 C xA, yA - Cell face areas
35 C cg2d_b - Conjugate Gradient 2-D solver : Right-hand side vector
36 C myThid - Instance number for this call of CALC_DIV_GHAT
37 INTEGER bi,bj,iMin,iMax,jMin,jMax
38 INTEGER K
39 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
40 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
41 _RL cg2d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
42 INTEGER myThid
43
44 C == Local variables ==
45 INTEGER i,j
46 _RL pf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
47
48 C-- Pressure equation source term
49 C Term is the vertical integral of the divergence of the
50 C time tendency terms along with a relaxation term that
51 C pulls div(U) + dh/dt back toward zero.
52
53 IF (implicDiv2Dflow.EQ.1.) then
54 C Fully Implicit treatment of the Barotropic Flow Divergence
55 DO j=1,sNy
56 DO i=1,sNx+1
57 pf(i,j) = xA(i,j)*gUNm1(i,j,k,bi,bj) / deltaTmom
58 ENDDO
59 ENDDO
60 ELSE
61 C Explicit+Implicit part of the Barotropic Flow Divergence
62 C => Filtering of uVel,vVel is necessary
63 #ifdef ALLOW_ZONAL_FILT
64 IF (zonal_filt_lat.LT.90.) THEN
65 CALL ZONAL_FILTER(
66 & uVel, hFacW, 1-1, sNy+1, k, k, bi, bj, 1, myThid)
67 CALL ZONAL_FILTER(
68 & vVel, hFacS, 1-1, sNy+1, k, k, bi, bj, 2, myThid)
69 ENDIF
70 #endif
71 DO j=1,sNy
72 DO i=1,sNx+1
73 pf(i,j) = ( implicDiv2Dflow * gUNm1(i,j,k,bi,bj)
74 & + (1.-implicDiv2Dflow) * uVel(i,j,k,bi,bj)
75 & ) * xA(i,j) / deltaTmom
76 ENDDO
77 ENDDO
78 ENDIF
79 DO j=1,sNy
80 DO i=1,sNx
81 cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj) +
82 & pf(i+1,j)-pf(i,j)
83 ENDDO
84 ENDDO
85
86 #ifdef ALLOW_NONHYDROSTATIC
87 IF (nonHydrostatic) THEN
88 DO j=1,sNy
89 DO i=1,sNx
90 cg3d_b(i,j,k,bi,bj) = pf(i+1,j)-pf(i,j)
91 ENDDO
92 ENDDO
93 ENDIF
94 #endif
95
96 IF (implicDiv2Dflow.EQ.1.) then
97 C Fully Implicit treatment of the Barotropic Flow Divergence
98 DO j=1,sNy+1
99 DO i=1,sNx
100 pf(i,j) = yA(i,j)*gVNm1(i,j,k,bi,bj) / deltatmom
101 ENDDO
102 ENDDO
103 ELSE
104 C Explicit+Implicit part of the Barotropic Flow Divergence
105 DO j=1,sNy+1
106 DO i=1,sNx
107 pf(i,j) = ( implicDiv2Dflow * gVNm1(i,j,k,bi,bj)
108 & + (1.-implicDiv2Dflow) * vVel(i,j,k,bi,bj)
109 & ) * yA(i,j) / deltaTmom
110 ENDDO
111 ENDDO
112 ENDIF
113 DO j=1,sNy
114 DO i=1,sNx
115 cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj) +
116 & pf(i,j+1)-pf(i,j)
117 ENDDO
118 ENDDO
119
120 #ifdef ALLOW_NONHYDROSTATIC
121 IF (nonHydrostatic) THEN
122 DO j=1,sNy
123 DO i=1,sNx
124 cg3d_b(i,j,k,bi,bj) = cg3d_b(i,j,k,bi,bj) +
125 & pf(i,j+1)-pf(i,j)
126 ENDDO
127 ENDDO
128 ENDIF
129 #endif
130
131 RETURN
132 END

  ViewVC Help
Powered by ViewVC 1.1.22