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

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

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


Revision 1.11 - (show annotations) (download)
Wed Apr 4 01:39:06 2007 UTC (17 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint59, checkpoint58y_post, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59m, checkpoint59l, checkpoint59n, checkpoint59i, checkpoint59h, checkpoint59k, checkpoint59j
Changes since 1.10: +18 -16 lines
add a logical argument "calcCFL" to DST horizontal Advection S/R
(if false, assume that uFld,vFld are already CFL number in x,y dir)

1 C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_fluxlimit_adv_x.F,v 1.10 2006/12/05 22:25:41 jmc Exp $
2 C $Name: $
3
4 #include "GAD_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: GAD_FLUXLIMIT_ADV_X
8
9 C !INTERFACE: ==========================================================
10 SUBROUTINE GAD_FLUXLIMIT_ADV_X(
11 I bi,bj,k, calcCFL, deltaTloc,
12 I uTrans, uFld,
13 I maskLocW, tracer,
14 O uT,
15 I myThid )
16
17 C !DESCRIPTION:
18 C Calculates the area integrated zonal flux due to advection of a tracer
19 C using second-order interpolation with a flux limiter:
20 C \begin{equation*}
21 C F^x_{adv} = U \overline{ \theta }^i
22 C - \frac{1}{2} \left(
23 C [ 1 - \psi(C_r) ] |U|
24 C + U \frac{u \Delta t}{\Delta x_c} \psi(C_r)
25 C \right) \delta_i \theta
26 C \end{equation*}
27 C where the $\psi(C_r)$ is the limiter function and $C_r$ is
28 C the slope ratio.
29
30 C !USES: ===============================================================
31 IMPLICIT NONE
32 #include "SIZE.h"
33 #include "GRID.h"
34
35 C !INPUT PARAMETERS: ===================================================
36 C bi,bj :: tile indices
37 C k :: vertical level
38 C calcCFL :: =T: calculate CFL number ; =F: take uFld as CFL
39 C deltaTloc :: local time-step (s)
40 C uTrans :: zonal volume transport
41 C uFld :: zonal flow / CFL number
42 C tracer :: tracer field
43 C myThid :: thread number
44 INTEGER bi,bj,k
45 LOGICAL calcCFL
46 _RL deltaTloc
47 _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
48 _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
49 _RS maskLocW(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
50 _RL tracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
51 INTEGER myThid
52
53 C !OUTPUT PARAMETERS: ==================================================
54 C uT :: zonal advective flux
55 _RL uT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
56
57 C !LOCAL VARIABLES: ====================================================
58 C i,j :: loop indices
59 C Cr :: slope ratio
60 C Rjm,Rj,Rjp :: differences at i-1,i,i+1
61 INTEGER i,j
62 _RL Cr,Rjm,Rj,Rjp
63 _RL uCFL
64 C Statement function provides Limiter(Cr)
65 #include "GAD_FLUX_LIMITER.h"
66 CEOP
67
68 DO j=1-Oly,sNy+Oly
69 uT(1-Olx,j)=0.
70 uT(2-Olx,j)=0.
71 uT(sNx+Olx,j)=0.
72 DO i=1-Olx+2,sNx+Olx-1
73
74 uCFL = uFld(i,j)
75 IF ( calcCFL ) uCFL = ABS( uFld(i,j)*deltaTloc
76 & *recip_dxC(i,j,bi,bj)*recip_deepFacC(k) )
77 Rjp=(tracer(i+1,j)-tracer( i ,j))*maskLocW(i+1,j)
78 Rj =(tracer( i ,j)-tracer(i-1,j))*maskLocW( i ,j)
79 Rjm=(tracer(i-1,j)-tracer(i-2,j))*maskLocW(i-1,j)
80
81 IF (Rj.NE.0.) THEN
82 IF (uTrans(i,j).GT.0) THEN
83 Cr=Rjm/Rj
84 ELSE
85 Cr=Rjp/Rj
86 ENDIF
87 ELSE
88 IF (uTrans(i,j).GT.0) THEN
89 Cr=Rjm*1.E20
90 ELSE
91 Cr=Rjp*1.E20
92 ENDIF
93 ENDIF
94 Cr=Limiter(Cr)
95 uT(i,j) =
96 & uTrans(i,j)*(Tracer(i,j)+Tracer(i-1,j))*0.5 _d 0
97 & -ABS(uTrans(i,j))*((1.-Cr)+uCFL*Cr)
98 & *Rj*0.5 _d 0
99 ENDDO
100 ENDDO
101
102 RETURN
103 END

  ViewVC Help
Powered by ViewVC 1.1.22