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

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

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


Revision 1.3 - (hide annotations) (download)
Tue Jul 31 15:01:33 2001 UTC (22 years, 10 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint40pre7, checkpoint40pre6, checkpoint40pre5
Changes since 1.2: +8 -3 lines
Initialization of Adams-Bashforth time stepping for iter.EQ.0
- uses simple forward step
- changes trajectories so *ALL* output is affected

1 adcroft 1.3 C $Header: /u/gcmpack/models/MITgcmUV/model/src/timestep_tracer.F,v 1.2 2001/05/29 14:01:37 adcroft Exp $
2 adcroft 1.2 C $Name: $
3 adcroft 1.1
4     #include "CPP_OPTIONS.h"
5    
6 adcroft 1.2 SUBROUTINE TIMESTEP_TRACER(
7     I bi, bj, iMin, iMax, jMin, jMax, K, tauAB,
8     U tracer, gTracer, gTrNm1,
9     I myIter, myThid )
10 adcroft 1.1 C /==========================================================\
11     C | S/R TIMESTEP_TRACER |
12     C | o Step model tracer field forward in time |
13     C \==========================================================/
14 adcroft 1.2 IMPLICIT NONE
15    
16     C == Global variables ===
17 adcroft 1.1 #include "SIZE.h"
18     #include "EEPARAMS.h"
19     #include "PARAMS.h"
20 adcroft 1.2
21 adcroft 1.1 C == Routine Arguments ==
22 adcroft 1.2 C tauAB - Adams-Bashforth timestepping weight: 0=forward ; 1/2=Adams-Bashf.
23 adcroft 1.1 INTEGER bi,bj,iMin,iMax,jMin,jMax,K
24 adcroft 1.2 _RL tauAB
25 adcroft 1.1 _RL tracer (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
26     _RL gTracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
27     _RL gTrNm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
28     INTEGER myIter, myThid
29     C == Local variables ==
30     INTEGER i,j
31     _RL ab15,ab05
32    
33     C Adams-Bashforth timestepping weights
34 adcroft 1.3 IF (myIter .EQ. 0) THEN
35     ab15=1.0
36     ab05=0.0
37     ELSE
38     ab15=1.+tauAB
39     ab05=-tauAB
40     ENDIF
41 adcroft 1.1
42     C Step forward temperature
43     DO j=jMin,jMax
44     DO i=iMin,iMax
45     gTrNm1(i,j,k,bi,bj)=tracer(i,j,k,bi,bj)
46     & +deltaTtracer*(
47     & ab15*gTracer(i,j,k,bi,bj)
48     & +ab05*gTrNm1(i,j,k,bi,bj) )
49     ENDDO
50     ENDDO
51    
52     RETURN
53     END

  ViewVC Help
Powered by ViewVC 1.1.22