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

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

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


Revision 1.16 - (hide annotations) (download)
Wed Sep 19 13:58:08 2001 UTC (22 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint40
Changes since 1.15: +33 -13 lines
"Volume exact-Conservation" modified for
non-linear free-surface + Crank-Nickelson

1 jmc 1.16 C $Header: /u/gcmpack/models/MITgcmUV/model/src/calc_div_ghat.F,v 1.15 2001/06/29 17:14:49 adcroft Exp $
2 jmc 1.12 C $Name: $
3 cnh 1.1
4 cnh 1.6 #include "CPP_OPTIONS.h"
5 cnh 1.1
6     C /==========================================================\
7     C | S/R CALC_DIV_GHAT |
8     C | o Form the right hand-side of the surface pressure eqn. |
9 adcroft 1.9 C |==========================================================|
10 cnh 1.1 C \==========================================================/
11     SUBROUTINE CALC_DIV_GHAT(
12 jmc 1.13 I bi,bj,iMin,iMax,jMin,jMax,K,
13 cnh 1.1 I xA,yA,
14 jmc 1.13 U cg2d_b,
15 cnh 1.1 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 adcroft 1.15 #include "SOLVE_FOR_PRESSURE3D.h"
27 cnh 1.1
28     C == Routine arguments ==
29     C bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
30     C results will be set.
31 adcroft 1.9 C k - Index of layer.
32     C xA, yA - Cell face areas
33 jmc 1.13 C cg2d_b - Conjugate Gradient 2-D solver : Right-hand side vector
34     C myThid - Instance number for this call of CALC_DIV_GHAT
35 cnh 1.1 INTEGER bi,bj,iMin,iMax,jMin,jMax
36     INTEGER K
37     _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
38     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
39 jmc 1.13 _RL cg2d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
40 cnh 1.1 INTEGER myThid
41    
42     C == Local variables ==
43     INTEGER i,j
44     _RL pf (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
45    
46     C-- Pressure equation source term
47     C Term is the vertical integral of the divergence of the
48     C time tendency terms along with a relaxation term that
49     C pulls div(U) + dh/dt back toward zero.
50    
51 jmc 1.16 IF (implicDiv2Dflow.EQ.1.) THEN
52 jmc 1.12 C Fully Implicit treatment of the Barotropic Flow Divergence
53     DO j=1,sNy
54     DO i=1,sNx+1
55     pf(i,j) = xA(i,j)*gUNm1(i,j,k,bi,bj) / deltaTmom
56     ENDDO
57     ENDDO
58 jmc 1.16 ELSEIF (exactConserv) THEN
59     c ELSEIF (nonlinFreeSurf.GT.0) THEN
60     C Implicit treatment of the Barotropic Flow Divergence
61     DO j=1,sNy
62     DO i=1,sNx+1
63     pf(i,j) = implicDiv2Dflow
64     & *xA(i,j)*gUNm1(i,j,k,bi,bj) / deltaTmom
65     ENDDO
66     ENDDO
67 jmc 1.12 ELSE
68     C Explicit+Implicit part of the Barotropic Flow Divergence
69     C => Filtering of uVel,vVel is necessary
70 jmc 1.16 C-- Now the filter are applied in the_correction_step().
71     C We have left this code here to indicate where the filters used to be
72     C in the algorithm before JMC moved them to after the pressure solver.
73     C-
74     C#ifdef ALLOW_ZONAL_FILT
75     C IF (zonal_filt_lat.LT.90.) THEN
76     C CALL ZONAL_FILTER(
77     C & uVel, hFacW, 1-1, sNy+1, k, k, bi, bj, 1, myThid)
78     C CALL ZONAL_FILTER(
79     C & vVel, hFacS, 1-1, sNy+1, k, k, bi, bj, 2, myThid)
80     C ENDIF
81     C#endif
82 jmc 1.12 DO j=1,sNy
83     DO i=1,sNx+1
84     pf(i,j) = ( implicDiv2Dflow * gUNm1(i,j,k,bi,bj)
85     & + (1.-implicDiv2Dflow) * uVel(i,j,k,bi,bj)
86     & ) * xA(i,j) / deltaTmom
87     ENDDO
88     ENDDO
89     ENDIF
90 cnh 1.1 DO j=1,sNy
91     DO i=1,sNx
92 adcroft 1.8 cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj) +
93 cnh 1.1 & pf(i+1,j)-pf(i,j)
94     ENDDO
95     ENDDO
96    
97 adcroft 1.8 #ifdef ALLOW_NONHYDROSTATIC
98     IF (nonHydrostatic) THEN
99     DO j=1,sNy
100     DO i=1,sNx
101     cg3d_b(i,j,k,bi,bj) = pf(i+1,j)-pf(i,j)
102     ENDDO
103     ENDDO
104     ENDIF
105     #endif
106    
107 jmc 1.16 IF (implicDiv2Dflow.EQ.1.) THEN
108 jmc 1.12 C Fully Implicit treatment of the Barotropic Flow Divergence
109     DO j=1,sNy+1
110     DO i=1,sNx
111     pf(i,j) = yA(i,j)*gVNm1(i,j,k,bi,bj) / deltatmom
112 jmc 1.16 ENDDO
113     ENDDO
114     ELSEIF (exactConserv) THEN
115     c ELSEIF (nonlinFreeSurf.GT.0) THEN
116     C Implicit treatment of the Barotropic Flow Divergence
117     DO j=1,sNy+1
118     DO i=1,sNx
119     pf(i,j) = implicDiv2Dflow
120     & *yA(i,j)*gVNm1(i,j,k,bi,bj) / deltatmom
121 jmc 1.12 ENDDO
122     ENDDO
123     ELSE
124     C Explicit+Implicit part of the Barotropic Flow Divergence
125     DO j=1,sNy+1
126     DO i=1,sNx
127     pf(i,j) = ( implicDiv2Dflow * gVNm1(i,j,k,bi,bj)
128     & + (1.-implicDiv2Dflow) * vVel(i,j,k,bi,bj)
129     & ) * yA(i,j) / deltaTmom
130     ENDDO
131     ENDDO
132     ENDIF
133 cnh 1.1 DO j=1,sNy
134     DO i=1,sNx
135 adcroft 1.8 cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj) +
136 cnh 1.1 & pf(i,j+1)-pf(i,j)
137     ENDDO
138     ENDDO
139 cnh 1.4
140 adcroft 1.8 #ifdef ALLOW_NONHYDROSTATIC
141     IF (nonHydrostatic) THEN
142     DO j=1,sNy
143     DO i=1,sNx
144     cg3d_b(i,j,k,bi,bj) = cg3d_b(i,j,k,bi,bj) +
145     & pf(i,j+1)-pf(i,j)
146     ENDDO
147     ENDDO
148     ENDIF
149     #endif
150 cnh 1.1
151     RETURN
152     END

  ViewVC Help
Powered by ViewVC 1.1.22