/[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.3 - (show annotations) (download)
Fri Apr 24 03:07:12 1998 UTC (26 years, 1 month ago) by cnh
Branch: MAIN
CVS Tags: kloop1
Changes since 1.2: +1 -9 lines
Removed development debuggin

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/timestep.F,v 1.2 1998/04/24 02:05:42 cnh Exp $
2
3 #include "CPP_EEOPTIONS.h"
4
5 C /==========================================================\
6 C | S/R TIMESTEP |
7 C | o Step model fields forward in time |
8 C \==========================================================/
9 SUBROUTINE TIMESTEP( bi, bj, iMin, iMax, jMin, jMax, myThid )
10 implicit none
11 ! Common
12 #include "SIZE.h"
13 #include "DYNVARS.h"
14 #include "PARAMS.h"
15 #include "GRID.h"
16 #include "EEPARAMS.h"
17 #include "CG2D.h"
18 C == Routine Arguments ==
19 INTEGER bi,bj,iMin,iMax,jMin,jMax
20 INTEGER myThid
21 C == Local variables ==
22 C pg - Pressure gradient terms. Note cg2d_x
23 C holds term in units so that lateral
24 C gradient is all that is needed.
25 INTEGER i,j,k
26 REAL ab15,ab05
27 _RL pg(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
28
29 C Adams-Bashforth timestepping weights
30 ab15=1.5+abeps
31 ab05=-0.5-abeps
32
33 C Zonal pressure term
34 DO j=jMin,jMax
35 DO i=iMin,iMax
36 pg(i,j)=rDxC(i,j,bi,bj)*
37 & (cg2d_x(i,j,bi,bj)-cg2d_x(i-1,j,bi,bj))
38 & *gravity*rhonil
39 ENDDO
40 ENDDO
41 C Step forward zonal velocity
42 DO k=1,Nz
43 DO j=jMin,jMax
44 DO i=iMin,iMax
45 uVel(i,j,k,bi,bj)=uVel(i,j,k,bi,bj)
46 & +deltaTmom*(ab15*gU(i,j,k,bi,bj)+ab05*gUNm1(i,j,k,bi,bj)
47 & -pg(i,j)/rhonil
48 & )*maskW(i,j,k,bi,bj)
49 gUNm1(i,j,k,bi,bj)=gU(i,j,k,bi,bj)
50 ENDDO
51 ENDDO
52 ENDDO
53 C Meridional pressure term
54 DO j=jMin,jMax
55 DO i=iMin,iMax
56 pg(i,j)=rDyC(i,j,bi,bj)*
57 & (cg2d_x(i,j,bi,bj)-cg2d_x(i,j-1,bi,bj))
58 & *gravity*rhonil
59 ENDDO
60 ENDDO
61 C Step forward meridional velocity
62 DO k=1,Nz
63 DO j=jMin,jMax
64 DO i=iMin,iMax
65 vVel(i,j,k,bi,bj)=vVel(i,j,k,bi,bj)
66 & +deltaTmom*(ab15*gV(i,j,k,bi,bj)+ab05*gVNm1(i,j,k,bi,bj)
67 & -pg(i,j)/rhonil
68 & )*maskS(i,j,k,bi,bj)
69 gVNm1(i,j,k,bi,bj)=gV(i,j,k,bi,bj)
70 ENDDO
71 ENDDO
72 ENDDO
73 C Step forward temperature
74 DO k=1,Nz
75 DO j=jMin,jMax
76 DO i=iMin,iMax
77 theta(i,j,k,bi,bj)=theta(i,j,k,bi,bj)
78 & +deltaTtracer*(ab15*gT(i,j,k,bi,bj)+ab05*gTNm1(i,j,k,bi,bj))
79 gTNm1(i,j,k,bi,bj)=gT(i,j,k,bi,bj)
80 ENDDO
81 ENDDO
82 ENDDO
83
84 _BARRIER
85 C CALL PLOT_FIELD_XYZR8( uVel, 'TIEMSTEP.1 uVel',Nz,1,myThid)
86 C CALL PLOT_FIELD_XYZR8( vVel, 'TIEMSTEP.1 vVel',Nz,1,myThid)
87 C CALL PLOT_FIELD_XYZR8( theta, 'TIEMSTEP.1 theta',Nz,1,myThid)
88
89
90 RETURN
91 END

  ViewVC Help
Powered by ViewVC 1.1.22