/[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.22 - (show annotations) (download)
Tue Mar 6 17:10:29 2001 UTC (23 years, 2 months ago) by jmc
Branch: MAIN
Changes since 1.21: +1 -2 lines
remove "include CG2D.h"

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

  ViewVC Help
Powered by ViewVC 1.1.22