/[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.20 - (show annotations) (download)
Tue Apr 6 00:31:54 2004 UTC (20 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint57g_pre, checkpoint57s_post, checkpoint57b_post, checkpoint57g_post, checkpoint56b_post, checkpoint57y_post, checkpoint54d_post, checkpoint54e_post, checkpoint57r_post, checkpoint57d_post, checkpoint57i_post, checkpoint58, checkpoint55, checkpoint54, checkpoint57, checkpoint56, checkpoint53, checkpoint57n_post, checkpoint58a_post, checkpoint57z_post, checkpoint54f_post, checkpoint55i_post, checkpoint57l_post, checkpoint57t_post, checkpoint55c_post, checkpoint57v_post, checkpoint57f_post, checkpoint53d_post, checkpoint57a_post, checkpoint57h_pre, checkpoint54b_post, checkpoint57h_post, checkpoint52m_post, checkpoint57y_pre, checkpoint55g_post, checkpoint57c_post, checkpoint55d_post, checkpoint54a_pre, checkpoint53c_post, checkpoint55d_pre, checkpoint57c_pre, checkpoint55j_post, checkpoint54a_post, checkpoint55h_post, checkpoint57e_post, checkpoint55b_post, checkpoint53a_post, checkpoint55f_post, checkpoint53g_post, checkpoint57p_post, checkpint57u_post, checkpoint57q_post, eckpoint57e_pre, checkpoint56a_post, checkpoint53f_post, checkpoint57h_done, checkpoint57j_post, checkpoint57f_pre, checkpoint52n_post, checkpoint53b_pre, checkpoint56c_post, checkpoint57a_pre, checkpoint55a_post, checkpoint57o_post, checkpoint57k_post, checkpoint53b_post, checkpoint57w_post, checkpoint57x_post, checkpoint53d_pre, checkpoint55e_post, checkpoint54c_post
Changes since 1.19: +1 -3 lines
include FFIELDS.h not needed ;

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_div_ghat.F,v 1.19 2003/10/09 04:19:18 edhill 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, iMin, iMax, jMin, jMax - Range of points for which calculation
38 C results will be set.
39 C k - Index of layer.
40 C xA, yA - Cell face areas
41 C cg2d_b - Conjugate Gradient 2-D solver : Right-hand side vector
42 C myThid - Instance number for this call of CALC_DIV_GHAT
43 INTEGER bi,bj,iMin,iMax,jMin,jMax
44 INTEGER K
45 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
46 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
47 _RL cg2d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
48 INTEGER myThid
49
50 C !LOCAL VARIABLES:
51 C == Local variables ==
52 C i,j :: Loop counters
53 C pf :: Intermediate array for building RHS source term.
54 INTEGER i,j
55 _RL pf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
56 CEOP
57
58 C-- Pressure equation source term
59 C Term is the vertical integral of the divergence of the
60 C time tendency terms along with a relaxation term that
61 C pulls div(U) + dh/dt back toward zero.
62
63 IF (implicDiv2Dflow.EQ.1.) THEN
64 C Fully Implicit treatment of the Barotropic Flow Divergence
65 DO j=1,sNy
66 DO i=1,sNx+1
67 pf(i,j) = xA(i,j)*gU(i,j,k,bi,bj) / deltaTmom
68 ENDDO
69 ENDDO
70 ELSEIF (exactConserv) THEN
71 c ELSEIF (nonlinFreeSurf.GT.0) THEN
72 C Implicit treatment of the Barotropic Flow Divergence
73 DO j=1,sNy
74 DO i=1,sNx+1
75 pf(i,j) = implicDiv2Dflow
76 & *xA(i,j)*gU(i,j,k,bi,bj) / deltaTmom
77 ENDDO
78 ENDDO
79 ELSE
80 C Explicit+Implicit part of the Barotropic Flow Divergence
81 C => Filtering of uVel,vVel is necessary
82 C-- Now the filter are applied in the_correction_step().
83 C We have left this code here to indicate where the filters used to be
84 C in the algorithm before JMC moved them to after the pressure solver.
85 C-
86 C#ifdef ALLOW_ZONAL_FILT
87 C IF (zonal_filt_lat.LT.90.) THEN
88 C CALL ZONAL_FILTER(
89 C & uVel, hFacW, 1-1, sNy+1, k, k, bi, bj, 1, myThid)
90 C CALL ZONAL_FILTER(
91 C & vVel, hFacS, 1-1, sNy+1, k, k, bi, bj, 2, myThid)
92 C ENDIF
93 C#endif
94 DO j=1,sNy
95 DO i=1,sNx+1
96 pf(i,j) = ( implicDiv2Dflow * gU(i,j,k,bi,bj)
97 & + (1.-implicDiv2Dflow) * uVel(i,j,k,bi,bj)
98 & ) * xA(i,j) / deltaTmom
99 ENDDO
100 ENDDO
101 ENDIF
102 DO j=1,sNy
103 DO i=1,sNx
104 cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj) +
105 & pf(i+1,j)-pf(i,j)
106 ENDDO
107 ENDDO
108
109 #ifdef ALLOW_NONHYDROSTATIC
110 IF (nonHydrostatic) THEN
111 DO j=1,sNy
112 DO i=1,sNx
113 cg3d_b(i,j,k,bi,bj) = pf(i+1,j)-pf(i,j)
114 ENDDO
115 ENDDO
116 ENDIF
117 #endif
118
119 IF (implicDiv2Dflow.EQ.1.) THEN
120 C Fully Implicit treatment of the Barotropic Flow Divergence
121 DO j=1,sNy+1
122 DO i=1,sNx
123 pf(i,j) = yA(i,j)*gV(i,j,k,bi,bj) / deltatmom
124 ENDDO
125 ENDDO
126 ELSEIF (exactConserv) THEN
127 c ELSEIF (nonlinFreeSurf.GT.0) THEN
128 C Implicit treatment of the Barotropic Flow Divergence
129 DO j=1,sNy+1
130 DO i=1,sNx
131 pf(i,j) = implicDiv2Dflow
132 & *yA(i,j)*gV(i,j,k,bi,bj) / deltatmom
133 ENDDO
134 ENDDO
135 ELSE
136 C Explicit+Implicit part of the Barotropic Flow Divergence
137 DO j=1,sNy+1
138 DO i=1,sNx
139 pf(i,j) = ( implicDiv2Dflow * gV(i,j,k,bi,bj)
140 & + (1.-implicDiv2Dflow) * vVel(i,j,k,bi,bj)
141 & ) * yA(i,j) / deltaTmom
142 ENDDO
143 ENDDO
144 ENDIF
145 DO j=1,sNy
146 DO i=1,sNx
147 cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj) +
148 & pf(i,j+1)-pf(i,j)
149 ENDDO
150 ENDDO
151
152 #ifdef ALLOW_NONHYDROSTATIC
153 IF (nonHydrostatic) THEN
154 DO j=1,sNy
155 DO i=1,sNx
156 cg3d_b(i,j,k,bi,bj) = cg3d_b(i,j,k,bi,bj) +
157 & pf(i,j+1)-pf(i,j)
158 ENDDO
159 ENDDO
160 ENDIF
161 #endif
162
163 RETURN
164 END

  ViewVC Help
Powered by ViewVC 1.1.22