/[MITgcm]/MITgcm/pkg/generic_advdiff/gad_u3_adv_r.F
ViewVC logotype

Contents of /MITgcm/pkg/generic_advdiff/gad_u3_adv_r.F

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


Revision 1.1 - (show annotations) (download)
Thu Jul 12 00:31:59 2001 UTC (22 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint40pre3, checkpoint40pre7, checkpoint40pre6, checkpoint40pre2, checkpoint40pre4, checkpoint40pre5
add a 3rd order advection scheme option

1 C $Header: $
2 C $Name: $
3
4 #include "GAD_OPTIONS.h"
5
6 SUBROUTINE GAD_U3_ADV_R(
7 I bi,bj,k,
8 I rTrans,
9 I tracer,
10 O wT,
11 I myThid )
12 C /==========================================================\
13 C | SUBROUTINE GAD_U3_ADV_R |
14 C | o Compute vertical advective Flux of Tracer using |
15 C | 3rd Order Upwind Scheme |
16 C |==========================================================|
17 IMPLICIT NONE
18
19 C == GLobal variables ==
20 #include "SIZE.h"
21 #include "GRID.h"
22 #include "GAD.h"
23
24 C == Routine arguments ==
25 INTEGER bi,bj,k
26 _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
27 _RL tracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
28 _RL wT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
29 INTEGER myThid
30
31 C == Local variables ==
32 INTEGER i,j,kp1,km1,km2
33 _RL Cr,Rjm,Rj,Rjp,Rjjm,Rjjp
34 #include "GAD_FLUX_LIMITER.h"
35
36 km2=MAX(1,k-2)
37 km1=MAX(1,k-1)
38 kp1=MIN(Nr,k+1)
39
40 IF ( k.EQ.1 .OR. k.GT.Nr ) THEN
41 DO j=1-Oly,sNy+Oly
42 DO i=1-Olx,sNx+Olx
43 wT(i,j) = 0.
44 ENDDO
45 ENDDO
46 ELSE
47 DO j=1-Oly,sNy+Oly
48 DO i=1-Olx,sNx+Olx
49 Rjp=(tracer(i,j,kp1,bi,bj)-tracer(i,j,k,bi,bj))
50 & *maskC(i,j,kp1,bi,bj)
51 Rj=(tracer(i,j,k,bi,bj)-tracer(i,j,km1,bi,bj))
52 Rjm=(tracer(i,j,km1,bi,bj)-tracer(i,j,km2,bi,bj))
53 & *maskC(i,j,km2,bi,bj)
54 Rjjp=Rjp-Rj
55 Rjjm=Rj-Rjm
56 wT(i,j) = maskC(i,j,km1,bi,bj)*(
57 & rTrans(i,j)*(
58 & (Tracer(i,j,k,bi,bj)+Tracer(i,j,km1,bi,bj))*0.5 _d 0
59 & -oneSixth*(Rjjm+Rjjp)*0.5 _d 0 )
60 & +ABS(rTrans(i,j))*
61 & oneSixth*(Rjjm-Rjjp)*0.5 _d 0
62 & )
63 ENDDO
64 ENDDO
65 ENDIF
66
67 RETURN
68 END

  ViewVC Help
Powered by ViewVC 1.1.22