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

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

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


Revision 1.29 - (hide annotations) (download)
Sat Feb 8 02:09:20 2003 UTC (21 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint48d_pre
Changes since 1.28: +28 -8 lines
in preparation for r*:
 new S/R (calc_grad_phi_hyd.F) to compute Hydrostatic potential gradient.
 pass the 2 comp. of the grad. as arguments to momentum S/R.
for the moment, only used if it does not change the results.

1 jmc 1.29 C $Header: /u/gcmpack/MITgcm/model/src/timestep.F,v 1.28 2003/01/26 21:06:11 jmc Exp $
2 jmc 1.21 C $Name: $
3 cnh 1.1
4 adcroft 1.10 #include "CPP_OPTIONS.h"
5 cnh 1.1
6 cnh 1.27 CBOP
7     C !ROUTINE: TIMESTEP
8     C !INTERFACE:
9 jmc 1.21 SUBROUTINE TIMESTEP( bi, bj, iMin, iMax, jMin, jMax, K,
10 jmc 1.29 I phiHyd, dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
11 adcroft 1.17 I myIter, myThid )
12 cnh 1.27 C !DESCRIPTION: \bv
13     C *==========================================================*
14     C | S/R TIMESTEP
15     C | o Step model fields forward in time
16     C *==========================================================*
17     C \ev
18    
19     C !USES:
20 jmc 1.21 IMPLICIT NONE
21 heimbach 1.18 C == Global variables ==
22 cnh 1.1 #include "SIZE.h"
23     #include "DYNVARS.h"
24 cnh 1.11 #include "EEPARAMS.h"
25 cnh 1.1 #include "PARAMS.h"
26     #include "GRID.h"
27 jmc 1.26 #include "SURFACE.h"
28 heimbach 1.18
29 cnh 1.27 C !INPUT/OUTPUT PARAMETERS:
30 cnh 1.1 C == Routine Arguments ==
31 jmc 1.29 C phiHyd :: Hydrostatic Potential (ocean: pressure/rho)
32     C (atmos: geopotential)
33     C dPhiHydX,Y :: Gradient (X & Y directions) of Hydrostatic Potential
34     C phiSurfX :: gradient of Surface potential (Pressure/rho, ocean)
35     C phiSurfY :: or geopotential (atmos) in X and Y direction
36 cnh 1.1 INTEGER bi,bj,iMin,iMax,jMin,jMax
37 adcroft 1.4 INTEGER K
38 adcroft 1.19 _RL phiHyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
39 jmc 1.29 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
40     _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
41 jmc 1.21 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
42     _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
43 adcroft 1.17 INTEGER myIter, myThid
44 heimbach 1.18
45 cnh 1.27 C !LOCAL VARIABLES:
46 cnh 1.1 C == Local variables ==
47 adcroft 1.4 INTEGER i,j
48 adcroft 1.7 _RL ab15,ab05
49 jmc 1.21 _RL phxFac,phyFac, psFac
50 jmc 1.26 _RL gUtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
51     _RL gVtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
52 cnh 1.27 CEOP
53 cnh 1.1
54     C Adams-Bashforth timestepping weights
55 adcroft 1.25 IF (myIter .EQ. 0) THEN
56     ab15=1.0
57     ab05=0.0
58     ELSE
59 adcroft 1.17 ab15=1.5+abeps
60     ab05=-0.5-abeps
61 adcroft 1.25 ENDIF
62 cnh 1.1
63 jmc 1.29 C-- stagger time step: grad Phi_Hyp is not in gU,gV => add it in this S/R
64     c IF (staggerTimeStep) THEN
65     IF (.FALSE.) THEN
66     phxFac = pfFacMom
67     phyFac = pfFacMom
68     ELSE
69     phxFac = 0.
70     phyFac = 0.
71     ENDIF
72    
73     C-- explicit part of the surface potential gradient is added in this S/R
74     psFac = pfFacMom*(1. _d 0 - implicSurfPress)
75    
76 jmc 1.26 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
77     C- Compute effective gU term (including Adams-Bashforth weights) :
78     DO j=jMin,jMax
79     DO i=iMin,iMax
80     gUtmp(i,j) = ab15*gU(i,j,k,bi,bj)
81     & + ab05*gUNm1(i,j,k,bi,bj)
82     #ifdef INCLUDE_CD_CODE
83     & + guCD(i,j,k,bi,bj)
84     #endif
85     ENDDO
86     ENDDO
87    
88     #ifdef NONLIN_FRSURF
89     IF (.NOT. vectorInvariantMomentum
90     & .AND. nonlinFreeSurf.GT.1) THEN
91 jmc 1.28 IF (select_rStar.GT.0) THEN
92     DO j=jMin,jMax
93     DO i=iMin,iMax
94     gUtmp(i,j) = gUtmp(i,j)/rStarExpW(i,j,bi,bj)
95     ENDDO
96     ENDDO
97     ELSE
98     DO j=jMin,jMax
99     DO i=iMin,iMax
100     IF ( k.EQ.ksurfW(i,j,bi,bj) ) THEN
101 jmc 1.26 gUtmp(i,j) = gUtmp(i,j)
102     & *hFacW(i,j,k,bi,bj)/hFac_surfW(i,j,bi,bj)
103 jmc 1.28 ENDIF
104     ENDDO
105 jmc 1.26 ENDDO
106 jmc 1.28 ENDIF
107 jmc 1.26 ENDIF
108     #endif
109    
110 adcroft 1.7 C Step forward zonal velocity (store in Gu)
111 adcroft 1.19 DO j=jMin,jMax
112 cnh 1.1 DO i=iMin,iMax
113 jmc 1.26 gUNm1(i,j,k,bi,bj) = uVel(i,j,k,bi,bj)
114     & +deltaTmom*(
115     & gUtmp(i,j)
116     & - psFac*phiSurfX(i,j)
117 jmc 1.29 & - phxFac*dPhiHydX(i,j)
118 adcroft 1.7 & )*_maskW(i,j,k,bi,bj)
119 cnh 1.1 ENDDO
120 adcroft 1.19 ENDDO
121    
122 jmc 1.29 c IF (.FALSE.) THEN
123 adcroft 1.19 IF (staggerTimeStep) THEN
124     C-- -grad Phi_Hyd has not been incorporated to gU and is added here:
125 jmc 1.23 phxFac = pfFacMom*deltaTmom
126 adcroft 1.19 DO j=jMin,jMax
127     DO i=iMin,iMax
128     gUNm1(i,j,k,bi,bj)=gUNm1(i,j,k,bi,bj)
129 jmc 1.29 c & - phxFac*dPhiHydX(i,j)*_maskW(i,j,k,bi,bj)
130 adcroft 1.19 & - _recip_dxC(i,j,bi,bj)
131     & *(phiHyd(i,j,k)-phiHyd(i-1,j,k))*phxFac
132     & *_maskW(i,j,k,bi,bj)
133     ENDDO
134     ENDDO
135     ENDIF
136 adcroft 1.12
137 jmc 1.26 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
138     C- Compute effective gV term (including Adams-Bashforth weights) :
139     DO j=jMin,jMax
140     DO i=iMin,iMax
141     gVtmp(i,j) = ab15*gV(i,j,k,bi,bj)
142     & + ab05*gVNm1(i,j,k,bi,bj)
143     #ifdef INCLUDE_CD_CODE
144     & + gvCD(i,j,k,bi,bj)
145     #endif
146     ENDDO
147     ENDDO
148    
149     #ifdef NONLIN_FRSURF
150     IF (.NOT. vectorInvariantMomentum
151     & .AND. nonlinFreeSurf.GT.1) THEN
152 jmc 1.28 IF (select_rStar.GT.0) THEN
153     DO j=jMin,jMax
154     DO i=iMin,iMax
155     gVtmp(i,j) = gVtmp(i,j)/rStarExpS(i,j,bi,bj)
156     ENDDO
157     ENDDO
158     ELSE
159     DO j=jMin,jMax
160     DO i=iMin,iMax
161     IF ( k.EQ.ksurfS(i,j,bi,bj) ) THEN
162 jmc 1.26 gVtmp(i,j) = gVtmp(i,j)
163     & *hFacS(i,j,k,bi,bj)/hFac_surfS(i,j,bi,bj)
164 jmc 1.28 ENDIF
165     ENDDO
166 jmc 1.26 ENDDO
167 jmc 1.28 ENDIF
168 jmc 1.26 ENDIF
169     #endif
170    
171 adcroft 1.7 C Step forward meridional velocity (store in Gv)
172 adcroft 1.19 DO j=jMin,jMax
173 cnh 1.1 DO i=iMin,iMax
174 jmc 1.26 gVNm1(i,j,k,bi,bj) = vVel(i,j,k,bi,bj)
175     & +deltaTmom*(
176     & gVtmp(i,j)
177     & - psFac*phiSurfY(i,j)
178 jmc 1.29 & - phyFac*dPhiHydY(i,j)
179 adcroft 1.7 & )*_maskS(i,j,k,bi,bj)
180 cnh 1.1 ENDDO
181 adcroft 1.19 ENDDO
182 adcroft 1.12
183 jmc 1.29 c IF (.FALSE.) THEN
184 adcroft 1.19 IF (staggerTimeStep) THEN
185     C-- -grad Phi_Hyd has not been incorporated to gV and is added here:
186 jmc 1.23 phyFac = pfFacMom*deltaTmom
187 adcroft 1.19 DO j=jMin,jMax
188     DO i=iMin,iMax
189     gVNm1(i,j,k,bi,bj)=gVNm1(i,j,k,bi,bj)
190 jmc 1.29 c & - phyFac*dPhiHydY(i,j)*_maskS(i,j,k,bi,bj)
191 adcroft 1.19 & - _recip_dyC(i,j,bi,bj)
192 adcroft 1.24 & *(phiHyd(i,j,k)-phiHyd(i,j-1,k))*phyFac
193 adcroft 1.19 & *_maskS(i,j,k,bi,bj)
194     ENDDO
195 adcroft 1.12 ENDDO
196     ENDIF
197 cnh 1.1
198     RETURN
199     END

  ViewVC Help
Powered by ViewVC 1.1.22