1 |
C $Header: /u/gcmpack/MITgcm/model/src/timestep.F,v 1.34 2003/10/28 22:57:59 edhill Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "PACKAGES_CONFIG.h" |
5 |
#include "CPP_OPTIONS.h" |
6 |
|
7 |
CBOP |
8 |
C !ROUTINE: TIMESTEP |
9 |
C !INTERFACE: |
10 |
SUBROUTINE TIMESTEP( bi, bj, iMin, iMax, jMin, jMax, k, |
11 |
I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY, |
12 |
I myTime, myIter, myThid ) |
13 |
C !DESCRIPTION: \bv |
14 |
C *==========================================================* |
15 |
C | S/R TIMESTEP |
16 |
C | o Step model fields forward in time |
17 |
C *==========================================================* |
18 |
C \ev |
19 |
|
20 |
C !USES: |
21 |
IMPLICIT NONE |
22 |
C == Global variables == |
23 |
#include "SIZE.h" |
24 |
#include "DYNVARS.h" |
25 |
#include "EEPARAMS.h" |
26 |
#include "PARAMS.h" |
27 |
#include "GRID.h" |
28 |
#include "SURFACE.h" |
29 |
|
30 |
C !INPUT/OUTPUT PARAMETERS: |
31 |
C == Routine Arguments == |
32 |
C dPhiHydX,Y :: Gradient (X & Y directions) of Hydrostatic Potential |
33 |
C phiSurfX :: gradient of Surface potential (Pressure/rho, ocean) |
34 |
C phiSurfY :: or geopotential (atmos) in X and Y direction |
35 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
36 |
INTEGER k |
37 |
_RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
38 |
_RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
39 |
_RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
40 |
_RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
41 |
_RL myTime |
42 |
INTEGER myIter, myThid |
43 |
|
44 |
C !LOCAL VARIABLES: |
45 |
C == Local variables == |
46 |
LOGICAL momForcing_In_AB |
47 |
INTEGER i,j |
48 |
_RL ab15,ab05 |
49 |
_RL phxFac,phyFac, psFac |
50 |
_RL gUtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
51 |
_RL gVtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
52 |
#ifdef ALLOW_CD_CODE |
53 |
_RL guCor(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
54 |
_RL gvCor(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
55 |
#endif |
56 |
CEOP |
57 |
|
58 |
C Adams-Bashforth timestepping weights |
59 |
IF (myIter .EQ. 0) THEN |
60 |
ab15=1.0 |
61 |
ab05=0.0 |
62 |
ELSE |
63 |
ab15=1.5+abeps |
64 |
ab05=-0.5-abeps |
65 |
ENDIF |
66 |
|
67 |
C-- stagger time step: grad Phi_Hyp is not in gU,gV => add it in this S/R |
68 |
IF (staggerTimeStep) THEN |
69 |
phxFac = pfFacMom |
70 |
phyFac = pfFacMom |
71 |
ELSE |
72 |
phxFac = 0. |
73 |
phyFac = 0. |
74 |
ENDIF |
75 |
|
76 |
C-- explicit part of the surface potential gradient is added in this S/R |
77 |
psFac = pfFacMom*(1. _d 0 - implicSurfPress) |
78 |
|
79 |
C-- including or excluding momentum forcing from Adams-Bashforth: |
80 |
momForcing_In_AB = forcing_In_AB |
81 |
momForcing_In_AB = .TRUE. |
82 |
|
83 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
84 |
|
85 |
C- Initialize local arrays (not really necessary but safer) |
86 |
DO j=1-Oly,sNy+Oly |
87 |
DO i=1-Olx,sNx+Olx |
88 |
gUtmp(i,j) = 0. _d 0 |
89 |
gVtmp(i,j) = 0. _d 0 |
90 |
#ifdef ALLOW_CD_CODE |
91 |
guCor(i,j) = 0. _d 0 |
92 |
gvCor(i,j) = 0. _d 0 |
93 |
#endif |
94 |
ENDDO |
95 |
ENDDO |
96 |
|
97 |
C-- Forcing term inside the Adams-Bashforth: |
98 |
IF (momForcing .AND. momForcing_In_AB) THEN |
99 |
CALL EXTERNAL_FORCING_U( |
100 |
I iMin,iMax,jMin,jMax,bi,bj,k, |
101 |
I myTime,myThid) |
102 |
CALL EXTERNAL_FORCING_V( |
103 |
I iMin,iMax,jMin,jMax,bi,bj,k, |
104 |
I myTime,myThid) |
105 |
ENDIF |
106 |
|
107 |
C- Compute effective gU,gV_[n+1/2] terms (including Adams-Bashforth weights) |
108 |
C and save gU,gV_[n] into guNm1,gvNm1 for the next time step. |
109 |
DO j=jMin,jMax |
110 |
DO i=iMin,iMax |
111 |
gUtmp(i,j) = ab15*gU(i,j,k,bi,bj) |
112 |
& + ab05*guNm1(i,j,k,bi,bj) |
113 |
gVtmp(i,j) = ab15*gV(i,j,k,bi,bj) |
114 |
& + ab05*gvNm1(i,j,k,bi,bj) |
115 |
guNm1(i,j,k,bi,bj)= gU(i,j,k,bi,bj) |
116 |
gvNm1(i,j,k,bi,bj)= gV(i,j,k,bi,bj) |
117 |
gU(i,j,k,bi,bj) = gUtmp(i,j) |
118 |
gV(i,j,k,bi,bj) = gVtmp(i,j) |
119 |
ENDDO |
120 |
ENDDO |
121 |
|
122 |
C-- Forcing term outside the Adams-Bashforth: |
123 |
C (not recommanded with CD-scheme ON) |
124 |
IF (momForcing .AND. .NOT.momForcing_In_AB) THEN |
125 |
CALL EXTERNAL_FORCING_U( |
126 |
I iMin,iMax,jMin,jMax,bi,bj,k, |
127 |
I myTime,myThid) |
128 |
CALL EXTERNAL_FORCING_V( |
129 |
I iMin,iMax,jMin,jMax,bi,bj,k, |
130 |
I myTime,myThid) |
131 |
IF (useCDscheme) THEN |
132 |
C- for CD-scheme, compute gU,Vtmp = gU,V^n + forcing |
133 |
DO j=jMin,jMax |
134 |
DO i=iMin,iMax |
135 |
gUtmp(i,j) = gU(i,j,k,bi,bj)-gUtmp(i,j) |
136 |
& + guNm1(i,j,k,bi,bj) |
137 |
gVtmp(i,j) = gV(i,j,k,bi,bj)-gVtmp(i,j) |
138 |
& + gvNm1(i,j,k,bi,bj) |
139 |
ENDDO |
140 |
ENDDO |
141 |
ELSE |
142 |
DO j=jMin,jMax |
143 |
DO i=iMin,iMax |
144 |
gUtmp(i,j) = gU(i,j,k,bi,bj) |
145 |
gVtmp(i,j) = gV(i,j,k,bi,bj) |
146 |
ENDDO |
147 |
ENDDO |
148 |
ENDIF |
149 |
ELSEIF ( useCDscheme) THEN |
150 |
DO j=jMin,jMax |
151 |
DO i=iMin,iMax |
152 |
gUtmp(i,j) = guNm1(i,j,k,bi,bj) |
153 |
gVtmp(i,j) = gvNm1(i,j,k,bi,bj) |
154 |
ENDDO |
155 |
ENDDO |
156 |
ENDIF |
157 |
|
158 |
#ifdef ALLOW_CD_CODE |
159 |
IF (useCDscheme) THEN |
160 |
C- Step forward D-grid velocity using C-grid gU,Vtmp = gU,V^n + forcing |
161 |
C and return coriolis terms on C-grid (guCor,gvCor) |
162 |
CALL CD_CODE_SCHEME( |
163 |
I bi,bj,k, dPhiHydX,dPhiHydY, gUtmp,gVtmp, |
164 |
O guCor,gvCor, |
165 |
I myTime, myIter, myThid) |
166 |
DO j=jMin,jMax |
167 |
DO i=iMin,iMax |
168 |
gUtmp(i,j) = gU(i,j,k,bi,bj) |
169 |
& + guCor(i,j) |
170 |
gVtmp(i,j) = gV(i,j,k,bi,bj) |
171 |
& + gvCor(i,j) |
172 |
ENDDO |
173 |
ENDDO |
174 |
ENDIF |
175 |
#endif /* ALLOW_CD_CODE */ |
176 |
|
177 |
#ifdef NONLIN_FRSURF |
178 |
IF (.NOT. vectorInvariantMomentum |
179 |
& .AND. nonlinFreeSurf.GT.1) THEN |
180 |
IF (select_rStar.GT.0) THEN |
181 |
DO j=jMin,jMax |
182 |
DO i=iMin,iMax |
183 |
gUtmp(i,j) = gUtmp(i,j)/rStarExpW(i,j,bi,bj) |
184 |
gVtmp(i,j) = gVtmp(i,j)/rStarExpS(i,j,bi,bj) |
185 |
ENDDO |
186 |
ENDDO |
187 |
ELSE |
188 |
DO j=jMin,jMax |
189 |
DO i=iMin,iMax |
190 |
IF ( k.EQ.ksurfW(i,j,bi,bj) ) THEN |
191 |
gUtmp(i,j) = gUtmp(i,j) |
192 |
& *hFacW(i,j,k,bi,bj)/hFac_surfW(i,j,bi,bj) |
193 |
ENDIF |
194 |
IF ( k.EQ.ksurfS(i,j,bi,bj) ) THEN |
195 |
gVtmp(i,j) = gVtmp(i,j) |
196 |
& *hFacS(i,j,k,bi,bj)/hFac_surfS(i,j,bi,bj) |
197 |
ENDIF |
198 |
ENDDO |
199 |
ENDDO |
200 |
ENDIF |
201 |
ENDIF |
202 |
#endif /* NONLIN_FRSURF */ |
203 |
|
204 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
205 |
|
206 |
C Step forward zonal velocity (store in Gu) |
207 |
DO j=jMin,jMax |
208 |
DO i=iMin,iMax |
209 |
gU(i,j,k,bi,bj) = uVel(i,j,k,bi,bj) |
210 |
& +deltaTmom*( |
211 |
& gUtmp(i,j) |
212 |
& - psFac*phiSurfX(i,j) |
213 |
& - phxFac*dPhiHydX(i,j) |
214 |
& )*_maskW(i,j,k,bi,bj) |
215 |
ENDDO |
216 |
ENDDO |
217 |
|
218 |
C Step forward meridional velocity (store in Gv) |
219 |
DO j=jMin,jMax |
220 |
DO i=iMin,iMax |
221 |
gV(i,j,k,bi,bj) = vVel(i,j,k,bi,bj) |
222 |
& +deltaTmom*( |
223 |
& gVtmp(i,j) |
224 |
& - psFac*phiSurfY(i,j) |
225 |
& - phyFac*dPhiHydY(i,j) |
226 |
& )*_maskS(i,j,k,bi,bj) |
227 |
ENDDO |
228 |
ENDDO |
229 |
|
230 |
RETURN |
231 |
END |