/[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.17 - (show annotations) (download)
Wed Sep 26 18:09:13 2001 UTC (22 years, 7 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint46n_post, checkpoint47e_post, checkpoint44e_post, checkpoint46l_post, checkpoint46g_pre, checkpoint47c_post, release1_p13_pre, checkpoint46f_post, checkpoint48e_post, checkpoint44f_post, checkpoint46b_post, checkpoint43a-release1mods, ecco_c50_e32, ecco_c50_e33, ecco_c50_e30, ecco_c50_e31, release1_p13, checkpoint48i_post, checkpoint46l_pre, chkpt44d_post, checkpoint50, release1_p8, release1_p9, release1_p1, release1_p2, release1_p3, release1_p4, release1_p5, release1_p6, release1_p7, checkpoint50b_pre, checkpoint44e_pre, release1_b1, checkpoint48b_post, checkpoint43, checkpoint48c_pre, checkpoint47d_pre, release1_chkpt44d_post, checkpoint47a_post, checkpoint48d_pre, checkpoint47i_post, release1_p11, checkpoint47d_post, icebear5, icebear4, icebear3, icebear2, checkpoint46d_pre, checkpoint48d_post, release1-branch_tutorials, checkpoint48f_post, checkpoint45d_post, checkpoint46j_pre, chkpt44a_post, checkpoint44h_pre, checkpoint48h_post, ecco_c50_e29, checkpoint46a_post, checkpoint47g_post, checkpoint46j_post, checkpoint46k_post, ecco_c50_e28, chkpt44c_pre, checkpoint48a_post, checkpoint45a_post, checkpoint50a_post, ecco_c44_e19, ecco_c44_e18, ecco_c44_e17, ecco_c44_e16, release1_p12, release1_p10, release1_p16, release1_p17, release1_p14, release1_p15, checkpoint47j_post, ecco_c50_e33a, branch-exfmods-tag, checkpoint44g_post, checkpoint46e_pre, checkpoint48c_post, checkpoint45b_post, checkpoint46b_pre, release1-branch-end, release1_final_v1, checkpoint46c_pre, checkpoint46, checkpoint47b_post, checkpoint44b_post, checkpoint46h_pre, checkpoint46m_post, checkpoint46a_pre, checkpoint45c_post, ecco_ice2, ecco_ice1, checkpoint44h_post, checkpoint46g_post, release1_p12_pre, ecco_c44_e22, ecco_c44_e25, checkpoint47f_post, chkpt44a_pre, checkpoint46i_post, ecco_c44_e23, ecco_c44_e20, ecco_c44_e21, ecco_c44_e26, ecco_c44_e27, ecco_c44_e24, checkpoint46c_post, ecco-branch-mod1, ecco-branch-mod2, ecco-branch-mod3, ecco-branch-mod4, ecco-branch-mod5, checkpoint46e_post, release1_beta1, checkpoint44b_pre, checkpoint42, checkpoint41, checkpoint47, checkpoint44, checkpoint45, checkpoint48, checkpoint49, checkpoint46h_post, chkpt44c_post, checkpoint48g_post, checkpoint47h_post, checkpoint44f_pre, checkpoint46d_post, release1-branch_branchpoint
Branch point for: c24_e25_ice, branch-exfmods-curt, release1_final, release1-branch, release1, ecco-branch, release1_50yr, icebear, release1_coupled
Changes since 1.16: +21 -9 lines
Bringing comments up to data and formatting for document extraction.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_div_ghat.F,v 1.16 2001/09/19 13:58:08 jmc Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5 CBOP
6 C !ROUTINE: CALC_DIV_GHAT
7 C !INTERFACE:
8 SUBROUTINE CALC_DIV_GHAT(
9 I bi,bj,iMin,iMax,jMin,jMax,K,
10 I xA,yA,
11 U cg2d_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 "DYNVARS.h"
29 #include "FFIELDS.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)*gUNm1(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)*gUNm1(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 * gUNm1(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)*gVNm1(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)*gVNm1(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 * gVNm1(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