/[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.19 - (hide 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 adcroft 1.19 C $Header: /u/gcmpack/models/MITgcmUV/model/src/timestep.F,v 1.18.2.4 2001/01/17 20:58:42 jmc Exp $
2 cnh 1.1
3 adcroft 1.10 #include "CPP_OPTIONS.h"
4 cnh 1.1
5     C /==========================================================\
6     C | S/R TIMESTEP |
7     C | o Step model fields forward in time |
8     C \==========================================================/
9 adcroft 1.4 SUBROUTINE TIMESTEP( bi, bj, iMin, iMax, jMin, jMax,
10 adcroft 1.19 I K, phiHyd,
11 adcroft 1.17 I myIter, myThid )
12 cnh 1.1 implicit none
13 heimbach 1.18
14     C == Global variables ==
15 cnh 1.1 #include "SIZE.h"
16     #include "DYNVARS.h"
17 cnh 1.11 #include "EEPARAMS.h"
18 cnh 1.1 #include "PARAMS.h"
19     #include "GRID.h"
20 heimbach 1.18
21 cnh 1.1 C == Routine Arguments ==
22 adcroft 1.19 C phiHyd - Hydrostatic pressure (ocean) or geopotentiel (atmos)
23 cnh 1.1 INTEGER bi,bj,iMin,iMax,jMin,jMax
24 adcroft 1.4 INTEGER K
25 adcroft 1.19 _RL phiHyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
26 adcroft 1.17 INTEGER myIter, myThid
27 heimbach 1.18
28 cnh 1.1 C == Local variables ==
29 adcroft 1.4 INTEGER i,j
30 adcroft 1.7 _RL ab15,ab05
31 adcroft 1.19 _RL phxFac,phyFac
32 cnh 1.1
33     C Adams-Bashforth timestepping weights
34 adcroft 1.17 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 cnh 1.1
42 adcroft 1.7 C Step forward zonal velocity (store in Gu)
43 adcroft 1.19 DO j=jMin,jMax
44 cnh 1.1 DO i=iMin,iMax
45 adcroft 1.19 gUNm1(i,j,k,bi,bj)=uVel(i,j,k,bi,bj)
46 cnh 1.1 & +deltaTmom*(ab15*gU(i,j,k,bi,bj)+ab05*gUNm1(i,j,k,bi,bj)
47 cnh 1.16 #ifdef INCLUDE_CD_CODE
48 cnh 1.9 & +guCD(i,j,k,bi,bj)
49     #endif
50 adcroft 1.7 & )*_maskW(i,j,k,bi,bj)
51 cnh 1.1 ENDDO
52 adcroft 1.19 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 adcroft 1.12
67 adcroft 1.7 C Step forward meridional velocity (store in Gv)
68 adcroft 1.19 DO j=jMin,jMax
69 cnh 1.1 DO i=iMin,iMax
70 adcroft 1.19 gVNm1(i,j,k,bi,bj)=vVel(i,j,k,bi,bj)
71 cnh 1.1 & +deltaTmom*(ab15*gV(i,j,k,bi,bj)+ab05*gVNm1(i,j,k,bi,bj)
72 cnh 1.16 #ifdef INCLUDE_CD_CODE
73 cnh 1.9 & +gvCD(i,j,k,bi,bj)
74     #endif
75 adcroft 1.7 & )*_maskS(i,j,k,bi,bj)
76 cnh 1.1 ENDDO
77 adcroft 1.19 ENDDO
78 adcroft 1.12
79 adcroft 1.19 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 adcroft 1.12 ENDDO
90     ENDIF
91 cnh 1.1
92     RETURN
93     END

  ViewVC Help
Powered by ViewVC 1.1.22