/[MITgcm]/MITgcm/model/src/timestep.F
ViewVC logotype

Diff of /MITgcm/model/src/timestep.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph | View Patch Patch

revision 1.24 by adcroft, Tue Jun 5 15:29:31 2001 UTC revision 1.27 by cnh, Wed Sep 26 18:09:16 2001 UTC
# Line 3  C $Name$ Line 3  C $Name$
3    
4  #include "CPP_OPTIONS.h"  #include "CPP_OPTIONS.h"
5    
6  C     /==========================================================\  CBOP
7  C     | S/R TIMESTEP                                             |  C     !ROUTINE: TIMESTEP
8  C     | o Step model fields forward in time                      |  C     !INTERFACE:
 C     \==========================================================/  
9        SUBROUTINE TIMESTEP( bi, bj, iMin, iMax, jMin, jMax, K,        SUBROUTINE TIMESTEP( bi, bj, iMin, iMax, jMin, jMax, K,
10       I                     phiHyd, phiSurfX, phiSurfY,       I                     phiHyd, phiSurfX, phiSurfY,
11       I                     myIter, myThid )       I                     myIter, myThid )
12        IMPLICIT NONE  C     !DESCRIPTION: \bv
13    C     *==========================================================*
14    C     | S/R TIMESTEP                                              
15    C     | o Step model fields forward in time                      
16    C     *==========================================================*
17    C     \ev
18    
19    C     !USES:
20          IMPLICIT NONE
21  C     == Global variables ==  C     == Global variables ==
22  #include "SIZE.h"  #include "SIZE.h"
23  #include "DYNVARS.h"  #include "DYNVARS.h"
24  #include "EEPARAMS.h"  #include "EEPARAMS.h"
25  #include "PARAMS.h"  #include "PARAMS.h"
26  #include "GRID.h"  #include "GRID.h"
27    #include "SURFACE.h"
28    
29    C     !INPUT/OUTPUT PARAMETERS:
30  C     == Routine Arguments ==  C     == Routine Arguments ==
31  C     phiHyd    - Hydrostatic Potential (ocean: pressure/rho)  C     phiHyd    - Hydrostatic Potential (ocean: pressure/rho)
32  C                                       (atmos: geopotentiel)  C                                       (atmos: geopotentiel)
# Line 31  C     phiSurfY             or geopotenti Line 39  C     phiSurfY             or geopotenti
39        _RL     phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL     phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
40        INTEGER myIter, myThid        INTEGER myIter, myThid
41    
42    C     !LOCAL VARIABLES:
43  C     == Local variables ==  C     == Local variables ==
44        INTEGER i,j        INTEGER i,j
45        _RL ab15,ab05        _RL ab15,ab05
46        _RL phxFac,phyFac, psFac        _RL phxFac,phyFac, psFac
47          _RL     gUtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
48          _RL     gVtmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
49    CEOP
50    
51  C     Adams-Bashforth timestepping weights  C     Adams-Bashforth timestepping weights
52  Caja  IF (myIter .EQ. 0) THEN        IF (myIter .EQ. 0) THEN
53  Caja   ab15=1.0         ab15=1.0
54  Caja   ab05=0.0         ab05=0.0
55  Caja  ELSE        ELSE
56         ab15=1.5+abeps         ab15=1.5+abeps
57         ab05=-0.5-abeps         ab05=-0.5-abeps
58  Caja  ENDIF        ENDIF
59    
60    C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
61    C-    Compute effective gU term (including Adams-Bashforth weights) :
62          DO j=jMin,jMax
63           DO i=iMin,iMax
64            gUtmp(i,j) = ab15*gU(i,j,k,bi,bj)
65         &             + ab05*gUNm1(i,j,k,bi,bj)  
66    #ifdef INCLUDE_CD_CODE
67         &             + guCD(i,j,k,bi,bj)
68    #endif
69           ENDDO
70          ENDDO
71          
72    #ifdef NONLIN_FRSURF
73          IF (.NOT. vectorInvariantMomentum
74         &    .AND. nonlinFreeSurf.GT.1) THEN
75           DO j=jMin,jMax
76            DO i=iMin,iMax
77             IF ( k.EQ.ksurfW(i,j,bi,bj) ) THEN
78               gUtmp(i,j) = gUtmp(i,j)
79         &         *hFacW(i,j,k,bi,bj)/hFac_surfW(i,j,bi,bj)
80             ENDIF
81            ENDDO
82           ENDDO
83          ENDIF
84    #endif
85    
86  C     Step forward zonal velocity (store in Gu)  C     Step forward zonal velocity (store in Gu)
87        psFac = pfFacMom*(1. _d 0 - implicSurfPress)        psFac = pfFacMom*(1. _d 0 - implicSurfPress)
88        DO j=jMin,jMax        DO j=jMin,jMax
89          DO i=iMin,iMax          DO i=iMin,iMax
90            gUNm1(i,j,k,bi,bj)=uVel(i,j,k,bi,bj) + deltaTmom * (            gUNm1(i,j,k,bi,bj) = uVel(i,j,k,bi,bj)
91       &       ab15*gU(i,j,k,bi,bj)+ab05*gUNm1(i,j,k,bi,bj)       &     +deltaTmom*(
92       &     - psFac*phiSurfX(i,j)       &         gUtmp(i,j)
93  #ifdef INCLUDE_CD_CODE       &       - psFac*phiSurfX(i,j)
      &     + guCD(i,j,k,bi,bj)  
 #endif  
94       &        )*_maskW(i,j,k,bi,bj)       &        )*_maskW(i,j,k,bi,bj)
95          ENDDO          ENDDO
96        ENDDO        ENDDO
# Line 72  C--   -grad Phi_Hyd has not been incorpo Line 108  C--   -grad Phi_Hyd has not been incorpo
108          ENDDO          ENDDO
109        ENDIF        ENDIF
110    
111    C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
112    C-    Compute effective gV term (including Adams-Bashforth weights) :
113          DO j=jMin,jMax
114           DO i=iMin,iMax
115            gVtmp(i,j) = ab15*gV(i,j,k,bi,bj)
116         &             + ab05*gVNm1(i,j,k,bi,bj)  
117    #ifdef INCLUDE_CD_CODE
118         &             + gvCD(i,j,k,bi,bj)
119    #endif
120           ENDDO
121          ENDDO
122          
123    #ifdef NONLIN_FRSURF
124          IF (.NOT. vectorInvariantMomentum
125         &    .AND. nonlinFreeSurf.GT.1) THEN
126           DO j=jMin,jMax
127            DO i=iMin,iMax
128             IF ( k.EQ.ksurfS(i,j,bi,bj) ) THEN
129               gVtmp(i,j) = gVtmp(i,j)
130         &         *hFacS(i,j,k,bi,bj)/hFac_surfS(i,j,bi,bj)
131             ENDIF
132            ENDDO
133           ENDDO
134          ENDIF
135    #endif
136    
137  C     Step forward meridional velocity (store in Gv)  C     Step forward meridional velocity (store in Gv)
138        psFac = pfFacMom*(1. _d 0 - implicSurfPress)        psFac = pfFacMom*(1. _d 0 - implicSurfPress)
139        DO j=jMin,jMax        DO j=jMin,jMax
140          DO i=iMin,iMax          DO i=iMin,iMax
141            gVNm1(i,j,k,bi,bj)=vVel(i,j,k,bi,bj) + deltaTmom * (            gVNm1(i,j,k,bi,bj) = vVel(i,j,k,bi,bj)
142       &       ab15*gV(i,j,k,bi,bj)+ab05*gVNm1(i,j,k,bi,bj)       &     +deltaTmom*(
143       &     - psFac*phiSurfY(i,j)       &         gVtmp(i,j)
144  #ifdef INCLUDE_CD_CODE       &       - psFac*phiSurfY(i,j)
      &     + gvCD(i,j,k,bi,bj)  
 #endif  
145       &        )*_maskS(i,j,k,bi,bj)       &        )*_maskS(i,j,k,bi,bj)
146          ENDDO          ENDDO
147        ENDDO        ENDDO

Legend:
Removed from v.1.24  
changed lines
  Added in v.1.27

  ViewVC Help
Powered by ViewVC 1.1.22