1 |
C $Header: /u/gcmpack/MITgcm/model/src/adams_bashforth3.F,v 1.2 2005/11/06 22:19:08 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CPP_OPTIONS.h" |
5 |
|
6 |
CBOP |
7 |
C !ROUTINE: ADAMS_BASHFORTH3 |
8 |
C !INTERFACE: |
9 |
SUBROUTINE ADAMS_BASHFORTH3( |
10 |
I bi, bj, kArg, |
11 |
U gTracer, gTrNm, |
12 |
I startAB, myIter, myThid ) |
13 |
C !DESCRIPTION: \bv |
14 |
C *==========================================================* |
15 |
C | S/R ADAMS_BASHFORTH3 |
16 |
C | o Extrapolate forward in time using third order |
17 |
C | Adams-Bashforth method. |
18 |
C *==========================================================* |
19 |
C | Either apply to tendency (kArg>0) at level k=kArg, |
20 |
C | or apply to state variable (kArg=0) for all levels |
21 |
C *==========================================================* |
22 |
C \ev |
23 |
C Extrapolate forward in time using 2 A.B. parameters (alpha,beta), |
24 |
C either tendency gX : |
25 |
C \begin{equation*} |
26 |
C gX^{n+1/2} = (1 + \alpha + \beta) gX^{n} |
27 |
C - (\alpha + 2 \beta) gX^{n-1} |
28 |
C + \beta gX^{n-2} |
29 |
C \end{equation*} |
30 |
C or state variable X : |
31 |
C \begin{equation*} |
32 |
C X^{n+1/2} = (1 + \alpha + \beta) X^{n} |
33 |
C - (\alpha + 2 \beta) X^{n-1} |
34 |
C + \beta X^{n-2} |
35 |
C \end{equation*} |
36 |
C with: |
37 |
C (alpha,beta)=(1/2,5/12) : AB-3, stable until CFL = 0.724 |
38 |
C (note: beta=0.281105 give the Max stability: 0.78616) |
39 |
C (alpha,beta)=(1/2+abEps,0) : return to previous quasi AB-2. |
40 |
C (alpha,beta)=(0,0) : 1rst.O forward time stepping |
41 |
|
42 |
C !USES: |
43 |
IMPLICIT NONE |
44 |
C == Global variables === |
45 |
#include "SIZE.h" |
46 |
#include "EEPARAMS.h" |
47 |
#include "PARAMS.h" |
48 |
|
49 |
C !INPUT/OUTPUT PARAMETERS: |
50 |
C == Routine Arguments == |
51 |
C bi,bj :: Tile indices |
52 |
C kArg :: if >0: apply AB on tendency at level k=kArg |
53 |
C :: if =0: apply AB on state variable and process all levels |
54 |
C gTracer :: (kArg>0) Tendency at current time ( units of quantity/sec ) |
55 |
C :: (kArg=0) state variable at current time |
56 |
C gTrNm :: (kArg>0) Tendencies at previous times ( units of quantity/sec ) |
57 |
C :: (kArg=0) state variable at previous times |
58 |
C startAB :: number of previous time level available to start/restart AB |
59 |
C myIter :: Current time step number |
60 |
C myThid :: my Thread number Id |
61 |
INTEGER bi,bj,kArg |
62 |
_RL gTracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
63 |
_RL gTrNm (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy,2) |
64 |
INTEGER startAB |
65 |
INTEGER myIter, myThid |
66 |
|
67 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
68 |
C !LOCAL VARIABLES: |
69 |
C == Local variables == |
70 |
C k :: level index |
71 |
C i,j :: Loop counters |
72 |
C m1,m2 :: indices for the 2 previous time-step Tendency |
73 |
C ab1,ab2,ab3 :: Adams bashforth extrapolation weights. |
74 |
INTEGER i,j, k, m1,m2 |
75 |
_RL ab0, ab1, ab2 |
76 |
_RL gTrtmp |
77 |
CEOP |
78 |
|
79 |
m1 = 1 + mod(myIter+1,2) |
80 |
m2 = 1 + mod( myIter ,2) |
81 |
|
82 |
C Adams-Bashforth timestepping weights |
83 |
IF ( myIter.EQ.nIter0 .AND. startAB.EQ.0 ) THEN |
84 |
ab0 = 1. _d 0 |
85 |
ab1 = 0. _d 0 |
86 |
ab2 = 0. _d 0 |
87 |
ELSEIF ( (myIter.EQ.nIter0 .AND. startAB.EQ.1) |
88 |
& .OR. (myIter.EQ.1+nIter0 .AND. startAB.EQ.0) ) THEN |
89 |
ab0 = 1. _d 0 + alph_AB |
90 |
ab1 = -alph_AB |
91 |
ab2 = 0. _d 0 |
92 |
ELSE |
93 |
ab0 = 1. _d 0 + alph_AB + beta_AB |
94 |
ab1 = -alph_AB - 2.*beta_AB |
95 |
ab2 = beta_AB |
96 |
ENDIF |
97 |
|
98 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
99 |
|
100 |
IF ( kArg.EQ.0 ) THEN |
101 |
C- Extrapolate forward in time the state variable, with AB weights: |
102 |
DO k=1,Nr |
103 |
DO j=1-Oly,sNy+Oly |
104 |
DO i=1-Olx,sNx+Olx |
105 |
gTrtmp = ab0*gTracer(i,j,k,bi,bj) |
106 |
& + ab1*gTrNm(i,j,k,bi,bj,m1) |
107 |
& + ab2*gTrNm(i,j,k,bi,bj,m2) |
108 |
gTrNm(i,j,k,bi,bj,m2) = gTracer(i,j,k,bi,bj) |
109 |
gTracer(i,j,k,bi,bj) = gTrtmp |
110 |
ENDDO |
111 |
ENDDO |
112 |
ENDDO |
113 |
ELSE |
114 |
C- Extrapolate forward in time the tendency, with AB weights: |
115 |
k = kArg |
116 |
DO j=1-Oly,sNy+Oly |
117 |
DO i=1-Olx,sNx+Olx |
118 |
gTrtmp = ab0*gTracer(i,j,k,bi,bj) |
119 |
& + ab1*gTrNm(i,j,k,bi,bj,m1) |
120 |
& + ab2*gTrNm(i,j,k,bi,bj,m2) |
121 |
gTrNm(i,j,k,bi,bj,m2) = gTracer(i,j,k,bi,bj) |
122 |
gTracer(i,j,k,bi,bj) = gTrtmp |
123 |
ENDDO |
124 |
ENDDO |
125 |
C--- |
126 |
ENDIF |
127 |
|
128 |
#endif /* ALLOW_ADAMSBASHFORTH_3 */ |
129 |
|
130 |
RETURN |
131 |
END |