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

Contents of /MITgcm/model/src/timestep_tracer.F

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


Revision 1.1 - (show annotations) (download)
Fri Feb 2 21:35:19 2001 UTC (23 years, 4 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint38, c37_adj, checkpoint39, checkpoint37, checkpoint36, checkpoint35
Branch point for: pre38
Merged changes from branch "branch-atmos-merge" into MAIN (checkpoint34)
 - substantial modifications to algorithm sequence (dynamics.F)
 - packaged OBCS, Shapiro filter, Zonal filter, Atmospheric Physics

1 C $Header: $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 C /==========================================================\
7 C | S/R TIMESTEP_TRACER |
8 C | o Step model tracer field forward in time |
9 C \==========================================================/
10 SUBROUTINE TIMESTEP_TRACER(
11 I bi, bj, iMin, iMax, jMin, jMax, K,
12 U tracer, gTracer, gTrNm1,
13 I myIter, myThid )
14 implicit none
15 ! Common
16 #include "SIZE.h"
17 #include "EEPARAMS.h"
18 #include "PARAMS.h"
19 C == Routine Arguments ==
20 INTEGER bi,bj,iMin,iMax,jMin,jMax,K
21 _RL tracer (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
22 _RL gTracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
23 _RL gTrNm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
24 INTEGER myIter, myThid
25 C == Local variables ==
26 INTEGER i,j
27 _RL ab15,ab05
28
29 C Adams-Bashforth timestepping weights
30 Caja IF (myIter .EQ. 0) THEN
31 Caja ab15=1.0
32 Caja ab05=0.0
33 Caja ELSE
34 ab15=1.5+abeps
35 ab05=-0.5-abeps
36 Caja ENDIF
37
38 C Step forward temperature
39 DO j=jMin,jMax
40 DO i=iMin,iMax
41 gTrNm1(i,j,k,bi,bj)=tracer(i,j,k,bi,bj)
42 & +deltaTtracer*(
43 & ab15*gTracer(i,j,k,bi,bj)
44 & +ab05*gTrNm1(i,j,k,bi,bj) )
45 ENDDO
46 ENDDO
47
48 RETURN
49 END

  ViewVC Help
Powered by ViewVC 1.1.22