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

1 jmc 1.22 C $Header: /u/gcmpack/models/MITgcmUV/model/src/timestep.F,v 1.21 2001/02/20 15:06:21 jmc Exp $
2 jmc 1.21 C $Name: $
3 cnh 1.1
4 adcroft 1.10 #include "CPP_OPTIONS.h"
5 cnh 1.1
6     C /==========================================================\
7     C | S/R TIMESTEP |
8     C | o Step model fields forward in time |
9     C \==========================================================/
10 jmc 1.21 SUBROUTINE TIMESTEP( bi, bj, iMin, iMax, jMin, jMax, K,
11     I phiHyd, phiSurfX, phiSurfY,
12 adcroft 1.17 I myIter, myThid )
13 jmc 1.21 IMPLICIT NONE
14 heimbach 1.18
15     C == Global variables ==
16 cnh 1.1 #include "SIZE.h"
17     #include "DYNVARS.h"
18 cnh 1.11 #include "EEPARAMS.h"
19 cnh 1.1 #include "PARAMS.h"
20     #include "GRID.h"
21 heimbach 1.18
22 cnh 1.1 C == Routine Arguments ==
23 jmc 1.21 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 cnh 1.1 INTEGER bi,bj,iMin,iMax,jMin,jMax
27 adcroft 1.4 INTEGER K
28 adcroft 1.19 _RL phiHyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
29 jmc 1.21 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
30     _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
31 adcroft 1.17 INTEGER myIter, myThid
32 heimbach 1.18
33 cnh 1.1 C == Local variables ==
34 adcroft 1.4 INTEGER i,j
35 adcroft 1.7 _RL ab15,ab05
36 jmc 1.21 _RL phxFac,phyFac, psFac
37 cnh 1.1
38     C Adams-Bashforth timestepping weights
39 adcroft 1.17 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 cnh 1.1
47 adcroft 1.7 C Step forward zonal velocity (store in Gu)
48 jmc 1.21 psFac = pfFacMom*(1. _d 0 - implicSurfPress)
49 adcroft 1.19 DO j=jMin,jMax
50 cnh 1.1 DO i=iMin,iMax
51 jmc 1.21 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 cnh 1.16 #ifdef INCLUDE_CD_CODE
55 jmc 1.21 & + guCD(i,j,k,bi,bj)
56 cnh 1.9 #endif
57 adcroft 1.7 & )*_maskW(i,j,k,bi,bj)
58 cnh 1.1 ENDDO
59 adcroft 1.19 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 adcroft 1.12
74 adcroft 1.7 C Step forward meridional velocity (store in Gv)
75 jmc 1.21 psFac = pfFacMom*(1. _d 0 - implicSurfPress)
76 adcroft 1.19 DO j=jMin,jMax
77 cnh 1.1 DO i=iMin,iMax
78 jmc 1.21 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 cnh 1.16 #ifdef INCLUDE_CD_CODE
82 jmc 1.21 & + gvCD(i,j,k,bi,bj)
83 cnh 1.9 #endif
84 adcroft 1.7 & )*_maskS(i,j,k,bi,bj)
85 cnh 1.1 ENDDO
86 adcroft 1.19 ENDDO
87 adcroft 1.12
88 adcroft 1.19 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 adcroft 1.12 ENDDO
99     ENDIF
100 cnh 1.1
101     RETURN
102     END

  ViewVC Help
Powered by ViewVC 1.1.22