/[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.6 - (show annotations) (download)
Sat May 30 02:10:16 1998 UTC (25 years, 11 months ago) by cnh
Branch: MAIN
CVS Tags: checkpoint4
Changes since 1.5: +3 -3 lines
Further memory saving macros for particular grids

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/timestep.F,v 1.5 1998/05/25 20:05:55 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,
10 I K, pSurfX, pSurfY,
11 I myThid )
12 implicit none
13 ! Common
14 #include "SIZE.h"
15 #include "DYNVARS.h"
16 #include "PARAMS.h"
17 #include "GRID.h"
18 #include "EEPARAMS.h"
19 #include "CG2D.h"
20 C == Routine Arguments ==
21 INTEGER bi,bj,iMin,iMax,jMin,jMax
22 INTEGER K
23 _RL pSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
24 _RL pSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
25 INTEGER myThid
26 C == Local variables ==
27 INTEGER i,j
28 REAL ab15,ab05
29 _RL hxFac, hyFac
30
31 C Adams-Bashforth timestepping weights
32 ab15=1.5+abeps
33 ab05=-0.5-abeps
34
35 C On/off scaling paramters
36 hxFac = pfFacMom
37 hyFac = pfFacMom
38
39 C Step forward zonal velocity
40 DO j=jMin,jMax
41 DO i=iMin,iMax
42 uVel(i,j,k,bi,bj)=uVel(i,j,k,bi,bj)
43 & +deltaTmom*(ab15*gU(i,j,k,bi,bj)+ab05*gUNm1(i,j,k,bi,bj)
44 & -hxFac*pSurfX(i,j)/rhonil
45 & )*_maskW(i,j,k,bi,bj)
46 gUNm1(i,j,k,bi,bj)=gU(i,j,k,bi,bj)
47 ENDDO
48 ENDDO
49 C Step forward meridional velocity
50 DO j=jMin,jMax
51 DO i=iMin,iMax
52 vVel(i,j,k,bi,bj)=vVel(i,j,k,bi,bj)
53 & +deltaTmom*(ab15*gV(i,j,k,bi,bj)+ab05*gVNm1(i,j,k,bi,bj)
54 & -hyFac*pSurfY(i,j)/rhonil
55 & )*_maskS(i,j,k,bi,bj)
56 gVNm1(i,j,k,bi,bj)=gV(i,j,k,bi,bj)
57 ENDDO
58 ENDDO
59 C Step forward temperature
60 DO j=jMin,jMax
61 DO i=iMin,iMax
62 theta(i,j,k,bi,bj)=theta(i,j,k,bi,bj)
63 & +deltaTtracer*(ab15*gT(i,j,k,bi,bj)+ab05*gTNm1(i,j,k,bi,bj))
64 gTNm1(i,j,k,bi,bj)=gT(i,j,k,bi,bj)
65 ENDDO
66 ENDDO
67
68 _BARRIER
69 C CALL PLOT_FIELD_XYZR8( uVel, 'TIEMSTEP.1 uVel',Nz,1,myThid)
70 C CALL PLOT_FIELD_XYZR8( vVel, 'TIEMSTEP.1 vVel',Nz,1,myThid)
71 C CALL PLOT_FIELD_XYZR8( theta, 'TIEMSTEP.1 theta',Nz,1,myThid)
72
73
74 RETURN
75 END

  ViewVC Help
Powered by ViewVC 1.1.22