/[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.27 - (show annotations) (download)
Fri Nov 9 22:37:05 2012 UTC (11 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint65, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65o, checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, HEAD
Changes since 1.26: +10 -7 lines
- move addMass common block from DYNVARS.h to FFIELDS.h

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

  ViewVC Help
Powered by ViewVC 1.1.22