/[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.20 - (show annotations) (download)
Sun Feb 4 14:38:48 2001 UTC (23 years, 4 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint35
Changes since 1.19: +2 -1 lines
Made sure each .F and .h file had
the CVS keywords Header and Name at its start.
Most had header but very few currently have Name, so
lots of changes!

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/timestep.F,v 1.19 2001/02/02 21:04:48 adcroft 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,
11 I K, phiHyd,
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 INTEGER bi,bj,iMin,iMax,jMin,jMax
25 INTEGER K
26 _RL phiHyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
27 INTEGER myIter, myThid
28
29 C == Local variables ==
30 INTEGER i,j
31 _RL ab15,ab05
32 _RL phxFac,phyFac
33
34 C Adams-Bashforth timestepping weights
35 Caja IF (myIter .EQ. 0) THEN
36 Caja ab15=1.0
37 Caja ab05=0.0
38 Caja ELSE
39 ab15=1.5+abeps
40 ab05=-0.5-abeps
41 Caja ENDIF
42
43 C Step forward zonal velocity (store in Gu)
44 DO j=jMin,jMax
45 DO i=iMin,iMax
46 gUNm1(i,j,k,bi,bj)=uVel(i,j,k,bi,bj)
47 & +deltaTmom*(ab15*gU(i,j,k,bi,bj)+ab05*gUNm1(i,j,k,bi,bj)
48 #ifdef INCLUDE_CD_CODE
49 & +guCD(i,j,k,bi,bj)
50 #endif
51 & )*_maskW(i,j,k,bi,bj)
52 ENDDO
53 ENDDO
54
55 IF (staggerTimeStep) THEN
56 C-- -grad Phi_Hyd has not been incorporated to gU and is added here:
57 phxFac = pfFacMom*deltaTmom*recip_rhoConst
58 DO j=jMin,jMax
59 DO i=iMin,iMax
60 gUNm1(i,j,k,bi,bj)=gUNm1(i,j,k,bi,bj)
61 & - _recip_dxC(i,j,bi,bj)
62 & *(phiHyd(i,j,k)-phiHyd(i-1,j,k))*phxFac
63 & *_maskW(i,j,k,bi,bj)
64 ENDDO
65 ENDDO
66 ENDIF
67
68 C Step forward meridional velocity (store in Gv)
69 DO j=jMin,jMax
70 DO i=iMin,iMax
71 gVNm1(i,j,k,bi,bj)=vVel(i,j,k,bi,bj)
72 & +deltaTmom*(ab15*gV(i,j,k,bi,bj)+ab05*gVNm1(i,j,k,bi,bj)
73 #ifdef INCLUDE_CD_CODE
74 & +gvCD(i,j,k,bi,bj)
75 #endif
76 & )*_maskS(i,j,k,bi,bj)
77 ENDDO
78 ENDDO
79
80 IF (staggerTimeStep) THEN
81 C-- -grad Phi_Hyd has not been incorporated to gV and is added here:
82 phyFac = pfFacMom*deltaTmom*recip_rhoConst
83 DO j=jMin,jMax
84 DO i=iMin,iMax
85 gVNm1(i,j,k,bi,bj)=gVNm1(i,j,k,bi,bj)
86 & - _recip_dyC(i,j,bi,bj)
87 & *(phiHyd(i,j,k)-phiHyd(i,j-1,k))*phyFac
88 & *_maskS(i,j,k,bi,bj)
89 ENDDO
90 ENDDO
91 ENDIF
92
93 RETURN
94 END

  ViewVC Help
Powered by ViewVC 1.1.22