/[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.26 - (show annotations) (download)
Mon Aug 27 18:48:57 2001 UTC (22 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint40pre9, checkpoint40
Changes since 1.25: +64 -13 lines
modified for NonLin-FreeSurf : This affects the truncation error
==> all output.txt need to be updated.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/timestep.F,v 1.25 2001/07/31 15:01:33 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, 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 #include "SURFACE.h"
22
23 C == Routine Arguments ==
24 C phiHyd - Hydrostatic Potential (ocean: pressure/rho)
25 C (atmos: geopotentiel)
26 C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean)
27 C phiSurfY or geopotentiel (atmos) in X and Y direction
28 INTEGER bi,bj,iMin,iMax,jMin,jMax
29 INTEGER K
30 _RL phiHyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
31 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
32 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
33 INTEGER myIter, myThid
34
35 C == Local variables ==
36 INTEGER i,j
37 _RL ab15,ab05
38 _RL phxFac,phyFac, psFac
39 _RL gUtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
40 _RL gVtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
41
42 C Adams-Bashforth timestepping weights
43 IF (myIter .EQ. 0) THEN
44 ab15=1.0
45 ab05=0.0
46 ELSE
47 ab15=1.5+abeps
48 ab05=-0.5-abeps
49 ENDIF
50
51 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
52 C- Compute effective gU term (including Adams-Bashforth weights) :
53 DO j=jMin,jMax
54 DO i=iMin,iMax
55 gUtmp(i,j) = ab15*gU(i,j,k,bi,bj)
56 & + ab05*gUNm1(i,j,k,bi,bj)
57 #ifdef INCLUDE_CD_CODE
58 & + guCD(i,j,k,bi,bj)
59 #endif
60 ENDDO
61 ENDDO
62
63 #ifdef NONLIN_FRSURF
64 IF (.NOT. vectorInvariantMomentum
65 & .AND. nonlinFreeSurf.GT.1) THEN
66 DO j=jMin,jMax
67 DO i=iMin,iMax
68 IF ( k.EQ.ksurfW(i,j,bi,bj) ) THEN
69 gUtmp(i,j) = gUtmp(i,j)
70 & *hFacW(i,j,k,bi,bj)/hFac_surfW(i,j,bi,bj)
71 ENDIF
72 ENDDO
73 ENDDO
74 ENDIF
75 #endif
76
77 C Step forward zonal velocity (store in Gu)
78 psFac = pfFacMom*(1. _d 0 - implicSurfPress)
79 DO j=jMin,jMax
80 DO i=iMin,iMax
81 gUNm1(i,j,k,bi,bj) = uVel(i,j,k,bi,bj)
82 & +deltaTmom*(
83 & gUtmp(i,j)
84 & - psFac*phiSurfX(i,j)
85 & )*_maskW(i,j,k,bi,bj)
86 ENDDO
87 ENDDO
88
89 IF (staggerTimeStep) THEN
90 C-- -grad Phi_Hyd has not been incorporated to gU and is added here:
91 phxFac = pfFacMom*deltaTmom
92 DO j=jMin,jMax
93 DO i=iMin,iMax
94 gUNm1(i,j,k,bi,bj)=gUNm1(i,j,k,bi,bj)
95 & - _recip_dxC(i,j,bi,bj)
96 & *(phiHyd(i,j,k)-phiHyd(i-1,j,k))*phxFac
97 & *_maskW(i,j,k,bi,bj)
98 ENDDO
99 ENDDO
100 ENDIF
101
102 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
103 C- Compute effective gV term (including Adams-Bashforth weights) :
104 DO j=jMin,jMax
105 DO i=iMin,iMax
106 gVtmp(i,j) = ab15*gV(i,j,k,bi,bj)
107 & + ab05*gVNm1(i,j,k,bi,bj)
108 #ifdef INCLUDE_CD_CODE
109 & + gvCD(i,j,k,bi,bj)
110 #endif
111 ENDDO
112 ENDDO
113
114 #ifdef NONLIN_FRSURF
115 IF (.NOT. vectorInvariantMomentum
116 & .AND. nonlinFreeSurf.GT.1) THEN
117 DO j=jMin,jMax
118 DO i=iMin,iMax
119 IF ( k.EQ.ksurfS(i,j,bi,bj) ) THEN
120 gVtmp(i,j) = gVtmp(i,j)
121 & *hFacS(i,j,k,bi,bj)/hFac_surfS(i,j,bi,bj)
122 ENDIF
123 ENDDO
124 ENDDO
125 ENDIF
126 #endif
127
128 C Step forward meridional velocity (store in Gv)
129 psFac = pfFacMom*(1. _d 0 - implicSurfPress)
130 DO j=jMin,jMax
131 DO i=iMin,iMax
132 gVNm1(i,j,k,bi,bj) = vVel(i,j,k,bi,bj)
133 & +deltaTmom*(
134 & gVtmp(i,j)
135 & - psFac*phiSurfY(i,j)
136 & )*_maskS(i,j,k,bi,bj)
137 ENDDO
138 ENDDO
139
140 IF (staggerTimeStep) THEN
141 C-- -grad Phi_Hyd has not been incorporated to gV and is added here:
142 phyFac = pfFacMom*deltaTmom
143 DO j=jMin,jMax
144 DO i=iMin,iMax
145 gVNm1(i,j,k,bi,bj)=gVNm1(i,j,k,bi,bj)
146 & - _recip_dyC(i,j,bi,bj)
147 & *(phiHyd(i,j,k)-phiHyd(i,j-1,k))*phyFac
148 & *_maskS(i,j,k,bi,bj)
149 ENDDO
150 ENDDO
151 ENDIF
152
153 RETURN
154 END

  ViewVC Help
Powered by ViewVC 1.1.22