1 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/timestep.F,v 1.22 2001/03/06 17:10:29 jmc 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 |
|
22 |
C == Routine Arguments == |
23 |
C phiHyd - Hydrostatic Potential (ocean: pressure/rho) |
24 |
C (atmos: geopotentiel) |
25 |
C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean) |
26 |
C phiSurfY or geopotentiel (atmos) in X and Y direction |
27 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
28 |
INTEGER K |
29 |
_RL phiHyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
30 |
_RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
31 |
_RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
32 |
INTEGER myIter, myThid |
33 |
|
34 |
C == Local variables == |
35 |
INTEGER i,j |
36 |
_RL ab15,ab05 |
37 |
_RL phxFac,phyFac, psFac |
38 |
|
39 |
C Adams-Bashforth timestepping weights |
40 |
Caja IF (myIter .EQ. 0) THEN |
41 |
Caja ab15=1.0 |
42 |
Caja ab05=0.0 |
43 |
Caja ELSE |
44 |
ab15=1.5+abeps |
45 |
ab05=-0.5-abeps |
46 |
Caja ENDIF |
47 |
|
48 |
C Step forward zonal velocity (store in Gu) |
49 |
psFac = pfFacMom*(1. _d 0 - implicSurfPress) |
50 |
DO j=jMin,jMax |
51 |
DO i=iMin,iMax |
52 |
gUNm1(i,j,k,bi,bj)=uVel(i,j,k,bi,bj) + deltaTmom * ( |
53 |
& ab15*gU(i,j,k,bi,bj)+ab05*gUNm1(i,j,k,bi,bj) |
54 |
& - psFac*phiSurfX(i,j) |
55 |
#ifdef INCLUDE_CD_CODE |
56 |
& + guCD(i,j,k,bi,bj) |
57 |
#endif |
58 |
& )*_maskW(i,j,k,bi,bj) |
59 |
ENDDO |
60 |
ENDDO |
61 |
|
62 |
IF (staggerTimeStep) THEN |
63 |
C-- -grad Phi_Hyd has not been incorporated to gU and is added here: |
64 |
phxFac = pfFacMom*deltaTmom |
65 |
DO j=jMin,jMax |
66 |
DO i=iMin,iMax |
67 |
gUNm1(i,j,k,bi,bj)=gUNm1(i,j,k,bi,bj) |
68 |
& - _recip_dxC(i,j,bi,bj) |
69 |
& *(phiHyd(i,j,k)-phiHyd(i-1,j,k))*phxFac |
70 |
& *_maskW(i,j,k,bi,bj) |
71 |
ENDDO |
72 |
ENDDO |
73 |
ENDIF |
74 |
|
75 |
C Step forward meridional velocity (store in Gv) |
76 |
psFac = pfFacMom*(1. _d 0 - implicSurfPress) |
77 |
DO j=jMin,jMax |
78 |
DO i=iMin,iMax |
79 |
gVNm1(i,j,k,bi,bj)=vVel(i,j,k,bi,bj) + deltaTmom * ( |
80 |
& ab15*gV(i,j,k,bi,bj)+ab05*gVNm1(i,j,k,bi,bj) |
81 |
& - psFac*phiSurfY(i,j) |
82 |
#ifdef INCLUDE_CD_CODE |
83 |
& + gvCD(i,j,k,bi,bj) |
84 |
#endif |
85 |
& )*_maskS(i,j,k,bi,bj) |
86 |
ENDDO |
87 |
ENDDO |
88 |
|
89 |
IF (staggerTimeStep) THEN |
90 |
C-- -grad Phi_Hyd has not been incorporated to gV and is added here: |
91 |
phyFac = pfFacMom*deltaTmom |
92 |
DO j=jMin,jMax |
93 |
DO i=iMin,iMax |
94 |
gVNm1(i,j,k,bi,bj)=gVNm1(i,j,k,bi,bj) |
95 |
& - _recip_dyC(i,j,bi,bj) |
96 |
& *(phiHyd(i,j,k)-phiHyd(i,j-1,k))*phyFac |
97 |
& *_maskS(i,j,k,bi,bj) |
98 |
ENDDO |
99 |
ENDDO |
100 |
ENDIF |
101 |
|
102 |
RETURN |
103 |
END |