/[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.19 - (show annotations) (download)
Fri Feb 2 21:04:48 2001 UTC (23 years, 4 months ago) by adcroft
Branch: MAIN
Changes since 1.18: +34 -25 lines
Merged changes from branch "branch-atmos-merge" into MAIN (checkpoint34)
 - substantial modifications to algorithm sequence (dynamics.F)
 - packaged OBCS, Shapiro filter, Zonal filter, Atmospheric Physics

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/timestep.F,v 1.18.2.4 2001/01/17 20:58:42 jmc Exp $
2
3 #include "CPP_OPTIONS.h"
4
5 C /==========================================================\
6 C | S/R TIMESTEP |
7 C | o Step model fields forward in time |
8 C \==========================================================/
9 SUBROUTINE TIMESTEP( bi, bj, iMin, iMax, jMin, jMax,
10 I K, phiHyd,
11 I myIter, myThid )
12 implicit none
13
14 C == Global variables ==
15 #include "SIZE.h"
16 #include "DYNVARS.h"
17 #include "EEPARAMS.h"
18 #include "PARAMS.h"
19 #include "GRID.h"
20
21 C == Routine Arguments ==
22 C phiHyd - Hydrostatic pressure (ocean) or geopotentiel (atmos)
23 INTEGER bi,bj,iMin,iMax,jMin,jMax
24 INTEGER K
25 _RL phiHyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
26 INTEGER myIter, myThid
27
28 C == Local variables ==
29 INTEGER i,j
30 _RL ab15,ab05
31 _RL phxFac,phyFac
32
33 C Adams-Bashforth timestepping weights
34 Caja IF (myIter .EQ. 0) THEN
35 Caja ab15=1.0
36 Caja ab05=0.0
37 Caja ELSE
38 ab15=1.5+abeps
39 ab05=-0.5-abeps
40 Caja ENDIF
41
42 C Step forward zonal velocity (store in Gu)
43 DO j=jMin,jMax
44 DO i=iMin,iMax
45 gUNm1(i,j,k,bi,bj)=uVel(i,j,k,bi,bj)
46 & +deltaTmom*(ab15*gU(i,j,k,bi,bj)+ab05*gUNm1(i,j,k,bi,bj)
47 #ifdef INCLUDE_CD_CODE
48 & +guCD(i,j,k,bi,bj)
49 #endif
50 & )*_maskW(i,j,k,bi,bj)
51 ENDDO
52 ENDDO
53
54 IF (staggerTimeStep) THEN
55 C-- -grad Phi_Hyd has not been incorporated to gU and is added here:
56 phxFac = pfFacMom*deltaTmom*recip_rhoConst
57 DO j=jMin,jMax
58 DO i=iMin,iMax
59 gUNm1(i,j,k,bi,bj)=gUNm1(i,j,k,bi,bj)
60 & - _recip_dxC(i,j,bi,bj)
61 & *(phiHyd(i,j,k)-phiHyd(i-1,j,k))*phxFac
62 & *_maskW(i,j,k,bi,bj)
63 ENDDO
64 ENDDO
65 ENDIF
66
67 C Step forward meridional velocity (store in Gv)
68 DO j=jMin,jMax
69 DO i=iMin,iMax
70 gVNm1(i,j,k,bi,bj)=vVel(i,j,k,bi,bj)
71 & +deltaTmom*(ab15*gV(i,j,k,bi,bj)+ab05*gVNm1(i,j,k,bi,bj)
72 #ifdef INCLUDE_CD_CODE
73 & +gvCD(i,j,k,bi,bj)
74 #endif
75 & )*_maskS(i,j,k,bi,bj)
76 ENDDO
77 ENDDO
78
79 IF (staggerTimeStep) THEN
80 C-- -grad Phi_Hyd has not been incorporated to gV and is added here:
81 phyFac = pfFacMom*deltaTmom*recip_rhoConst
82 DO j=jMin,jMax
83 DO i=iMin,iMax
84 gVNm1(i,j,k,bi,bj)=gVNm1(i,j,k,bi,bj)
85 & - _recip_dyC(i,j,bi,bj)
86 & *(phiHyd(i,j,k)-phiHyd(i,j-1,k))*phyFac
87 & *_maskS(i,j,k,bi,bj)
88 ENDDO
89 ENDDO
90 ENDIF
91
92 RETURN
93 END

  ViewVC Help
Powered by ViewVC 1.1.22