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

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

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


Revision 1.32 - (show annotations) (download)
Thu Apr 17 13:41:34 2003 UTC (21 years ago) by jmc
Branch: MAIN
CVS Tags: checkpoint50c_post, checkpoint50c_pre, checkpoint51, checkpoint50d_post, checkpoint51f_post, checkpoint51d_post, checkpoint51j_post, checkpoint51b_pre, checkpoint51h_pre, checkpoint50f_post, checkpoint50f_pre, branchpoint-genmake2, checkpoint51b_post, checkpoint51c_post, checkpoint50g_post, checkpoint50h_post, checkpoint50e_pre, checkpoint50i_post, checkpoint51i_pre, checkpoint50e_post, checkpoint50d_pre, checkpoint51e_post, checkpoint51f_pre, checkpoint51g_post, checkpoint50b_post, checkpoint51a_post
Branch point for: branch-genmake2
Changes since 1.31: +108 -43 lines
o store u*,v* in gU,V instead of in gu,vNm1
o to allow to put the momForcing out of the Adams-Bashforth:
  move forcing & CD-scheme calls from mom_fluxform & mom_vecinv
  to timestep.F
o new flag "useCDscheme" (default=F); replace guCD,gvCD by local arrays

1 C $Header: /u/gcmpack/MITgcm/model/src/timestep.F,v 1.31 2003/02/18 15:27:25 jmc Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: TIMESTEP
8 C !INTERFACE:
9 SUBROUTINE TIMESTEP( bi, bj, iMin, iMax, jMin, jMax, k,
10 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
11 I myTime, myIter, myThid )
12 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 IMPLICIT NONE
21 C == Global variables ==
22 #include "SIZE.h"
23 #include "DYNVARS.h"
24 #include "EEPARAMS.h"
25 #include "PARAMS.h"
26 #include "GRID.h"
27 #include "SURFACE.h"
28
29 C !INPUT/OUTPUT PARAMETERS:
30 C == Routine Arguments ==
31 C dPhiHydX,Y :: Gradient (X & Y directions) of Hydrostatic Potential
32 C phiSurfX :: gradient of Surface potential (Pressure/rho, ocean)
33 C phiSurfY :: or geopotential (atmos) in X and Y direction
34 INTEGER bi,bj,iMin,iMax,jMin,jMax
35 INTEGER k
36 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
37 _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
38 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
39 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
40 _RL myTime
41 INTEGER myIter, myThid
42
43 C !LOCAL VARIABLES:
44 C == Local variables ==
45 LOGICAL momForcing_In_AB
46 INTEGER i,j
47 _RL ab15,ab05
48 _RL phxFac,phyFac, psFac
49 _RL gUtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
50 _RL gVtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
51 #ifdef INCLUDE_CD_CODE
52 _RL guCor(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
53 _RL gvCor(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
54 #endif
55 CEOP
56
57 C Adams-Bashforth timestepping weights
58 IF (myIter .EQ. 0) THEN
59 ab15=1.0
60 ab05=0.0
61 ELSE
62 ab15=1.5+abeps
63 ab05=-0.5-abeps
64 ENDIF
65
66 C-- stagger time step: grad Phi_Hyp is not in gU,gV => add it in this S/R
67 IF (staggerTimeStep) THEN
68 phxFac = pfFacMom
69 phyFac = pfFacMom
70 ELSE
71 phxFac = 0.
72 phyFac = 0.
73 ENDIF
74
75 C-- explicit part of the surface potential gradient is added in this S/R
76 psFac = pfFacMom*(1. _d 0 - implicSurfPress)
77
78 C-- including or excluding momentum forcing from Adams-Bashforth:
79 momForcing_In_AB = forcing_In_AB
80 momForcing_In_AB = .TRUE.
81
82 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
83
84 C- Initialize local arrays (not really necessary but safer)
85 DO j=1-Oly,sNy+Oly
86 DO i=1-Olx,sNx+Olx
87 gUtmp(i,j) = 0. _d 0
88 gVtmp(i,j) = 0. _d 0
89 #ifdef INCLUDE_CD_CODE
90 guCor(i,j) = 0. _d 0
91 gvCor(i,j) = 0. _d 0
92 #endif
93 ENDDO
94 ENDDO
95
96 C-- Forcing term inside the Adams-Bashforth:
97 IF (momForcing .AND. momForcing_In_AB) THEN
98 CALL EXTERNAL_FORCING_U(
99 I iMin,iMax,jMin,jMax,bi,bj,k,
100 I myTime,myThid)
101 CALL EXTERNAL_FORCING_V(
102 I iMin,iMax,jMin,jMax,bi,bj,k,
103 I myTime,myThid)
104 ENDIF
105
106 C- Compute effective gU,gV_[n+1/2] terms (including Adams-Bashforth weights)
107 C and save gU,gV_[n] into guNm1,gvNm1 for the next time step.
108 DO j=jMin,jMax
109 DO i=iMin,iMax
110 gUtmp(i,j) = ab15*gU(i,j,k,bi,bj)
111 & + ab05*guNm1(i,j,k,bi,bj)
112 gVtmp(i,j) = ab15*gV(i,j,k,bi,bj)
113 & + ab05*gvNm1(i,j,k,bi,bj)
114 guNm1(i,j,k,bi,bj)= gU(i,j,k,bi,bj)
115 gvNm1(i,j,k,bi,bj)= gV(i,j,k,bi,bj)
116 gU(i,j,k,bi,bj) = gUtmp(i,j)
117 gV(i,j,k,bi,bj) = gVtmp(i,j)
118 ENDDO
119 ENDDO
120
121 C-- Forcing term outside the Adams-Bashforth:
122 C (not recommanded with CD-scheme ON)
123 IF (momForcing .AND. .NOT.momForcing_In_AB) THEN
124 CALL EXTERNAL_FORCING_U(
125 I iMin,iMax,jMin,jMax,bi,bj,k,
126 I myTime,myThid)
127 CALL EXTERNAL_FORCING_V(
128 I iMin,iMax,jMin,jMax,bi,bj,k,
129 I myTime,myThid)
130 IF (useCDscheme) THEN
131 C- for CD-scheme, compute gU,Vtmp = gU,V^n + forcing
132 DO j=jMin,jMax
133 DO i=iMin,iMax
134 gUtmp(i,j) = gU(i,j,k,bi,bj)-gUtmp(i,j)
135 & + guNm1(i,j,k,bi,bj)
136 gVtmp(i,j) = gV(i,j,k,bi,bj)-gVtmp(i,j)
137 & + gvNm1(i,j,k,bi,bj)
138 ENDDO
139 ENDDO
140 ELSE
141 DO j=jMin,jMax
142 DO i=iMin,iMax
143 gUtmp(i,j) = gU(i,j,k,bi,bj)
144 gVtmp(i,j) = gV(i,j,k,bi,bj)
145 ENDDO
146 ENDDO
147 ENDIF
148 ELSEIF ( useCDscheme) THEN
149 DO j=jMin,jMax
150 DO i=iMin,iMax
151 gUtmp(i,j) = guNm1(i,j,k,bi,bj)
152 gVtmp(i,j) = gvNm1(i,j,k,bi,bj)
153 ENDDO
154 ENDDO
155 ENDIF
156
157 #ifdef INCLUDE_CD_CODE
158 IF (useCDscheme) THEN
159 C- Step forward D-grid velocity using C-grid gU,Vtmp = gU,V^n + forcing
160 C and return coriolis terms on C-grid (guCor,gvCor)
161 CALL MOM_CDSCHEME(
162 I bi,bj,k, dPhiHydX,dPhiHydY, gUtmp,gVtmp,
163 O guCor,gvCor,
164 I myTime, myIter, myThid)
165 DO j=jMin,jMax
166 DO i=iMin,iMax
167 gUtmp(i,j) = gU(i,j,k,bi,bj)
168 & + guCor(i,j)
169 gVtmp(i,j) = gV(i,j,k,bi,bj)
170 & + gvCor(i,j)
171 ENDDO
172 ENDDO
173 ENDIF
174 #endif /* INCLUDE_CD_CODE */
175
176 #ifdef NONLIN_FRSURF
177 IF (.NOT. vectorInvariantMomentum
178 & .AND. nonlinFreeSurf.GT.1) THEN
179 IF (select_rStar.GT.0) THEN
180 DO j=jMin,jMax
181 DO i=iMin,iMax
182 gUtmp(i,j) = gUtmp(i,j)/rStarExpW(i,j,bi,bj)
183 gVtmp(i,j) = gVtmp(i,j)/rStarExpS(i,j,bi,bj)
184 ENDDO
185 ENDDO
186 ELSE
187 DO j=jMin,jMax
188 DO i=iMin,iMax
189 IF ( k.EQ.ksurfW(i,j,bi,bj) ) THEN
190 gUtmp(i,j) = gUtmp(i,j)
191 & *hFacW(i,j,k,bi,bj)/hFac_surfW(i,j,bi,bj)
192 ENDIF
193 IF ( k.EQ.ksurfS(i,j,bi,bj) ) THEN
194 gVtmp(i,j) = gVtmp(i,j)
195 & *hFacS(i,j,k,bi,bj)/hFac_surfS(i,j,bi,bj)
196 ENDIF
197 ENDDO
198 ENDDO
199 ENDIF
200 ENDIF
201 #endif /* NONLIN_FRSURF */
202
203 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
204
205 C Step forward zonal velocity (store in Gu)
206 DO j=jMin,jMax
207 DO i=iMin,iMax
208 gU(i,j,k,bi,bj) = uVel(i,j,k,bi,bj)
209 & +deltaTmom*(
210 & gUtmp(i,j)
211 & - psFac*phiSurfX(i,j)
212 & - phxFac*dPhiHydX(i,j)
213 & )*_maskW(i,j,k,bi,bj)
214 ENDDO
215 ENDDO
216
217 C Step forward meridional velocity (store in Gv)
218 DO j=jMin,jMax
219 DO i=iMin,iMax
220 gV(i,j,k,bi,bj) = vVel(i,j,k,bi,bj)
221 & +deltaTmom*(
222 & gVtmp(i,j)
223 & - psFac*phiSurfY(i,j)
224 & - phyFac*dPhiHydY(i,j)
225 & )*_maskS(i,j,k,bi,bj)
226 ENDDO
227 ENDDO
228
229 RETURN
230 END

  ViewVC Help
Powered by ViewVC 1.1.22