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

Annotation of /MITgcm/pkg/generic_advdiff/gad_dst2u1_adv_x.F

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


Revision 1.7 - (hide annotations) (download)
Wed Apr 4 01:39:06 2007 UTC (17 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59m, checkpoint59l, checkpoint59n, checkpoint59i, checkpoint59h, checkpoint59k, checkpoint59j, checkpoint59, checkpoint58y_post
Changes since 1.6: +9 -7 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 jmc 1.7 C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_dst2u1_adv_x.F,v 1.6 2006/12/05 22:21:50 jmc Exp $
2 jmc 1.1 C $Name: $
3    
4     #include "GAD_OPTIONS.h"
5    
6     CBOP
7     C !ROUTINE: GAD_DST2U1_ADV_X
8    
9     C !INTERFACE: ==========================================================
10 jmc 1.3 SUBROUTINE GAD_DST2U1_ADV_X(
11 jmc 1.7 I bi,bj,k, advectionScheme, calcCFL,
12 jmc 1.4 I deltaTloc, uTrans, uFld,
13 jmc 1.1 I tracer,
14     O uT,
15     I myThid )
16    
17     C !DESCRIPTION:
18     C Calculates the area integrated zonal flux due to advection
19 jmc 1.3 C of a tracer using second-order Direct Space and Time (DST-2)
20 jmc 1.1 C interpolation (=Lax-Wendroff) or simple 1rst order upwind scheme.
21    
22     C !USES: ===============================================================
23     IMPLICIT NONE
24     #include "SIZE.h"
25     #include "GRID.h"
26     #include "GAD.h"
27    
28     C !INPUT PARAMETERS: ===================================================
29     C bi,bj :: tile indices
30     C k :: vertical level
31 jmc 1.3 C advectionScheme :: advection scheme to use: either 2nd Order DST
32 jmc 1.1 C or 1rst Order Upwind
33 jmc 1.7 C calcCFL :: =T: calculate CFL number ; =F: take uFld as CFL
34     C deltaTloc :: local time-step (s)
35 jmc 1.1 C uTrans :: zonal volume transport
36 jmc 1.7 C uFld :: zonal flow / CFL number
37 jmc 1.1 C tracer :: tracer field
38     C myThid :: thread number
39     INTEGER bi,bj,k
40     INTEGER advectionScheme
41 jmc 1.7 LOGICAL calcCFL
42 jmc 1.1 _RL deltaTloc
43     _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
44 jmc 1.4 _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
45 jmc 1.1 _RL tracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
46     INTEGER myThid
47    
48     C !OUTPUT PARAMETERS: ==================================================
49     C uT :: zonal advective flux
50     _RL uT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
51    
52     C !LOCAL VARIABLES: ====================================================
53     C i,j :: loop indices
54     C rLimit :: centered (vs upwind) fraction
55     C uCFL :: Courant-Friedrich-Levy number
56     INTEGER i,j
57 jmc 1.7 _RL uCFL, xLimit, uAbs
58 jmc 1.1 CEOP
59    
60     xLimit = 0. _d 0
61     IF ( advectionScheme.EQ.ENUM_DST2 ) xLimit = 1. _d 0
62    
63     DO j=1-Oly,sNy+Oly
64     uT(1-Olx,j)=0.
65     DO i=1-Olx+1,sNx+Olx
66    
67 jmc 1.7 uCFL = uFld(i,j)
68     IF ( calcCFL ) uCFL = ABS( uFld(i,j)*deltaTloc
69 jmc 1.6 & *recip_dxC(i,j,bi,bj)*recip_deepFacC(k) )
70 jmc 1.1
71 jmc 1.3 c uT(i,j) =
72     c & uTrans(i,j)*(tracer(i-1,j)+tracer(i,j))*0.5 _d 0
73     c & + ( 1. _d 0 - xLimit*(1. _d 0 - uCFL) )*ABS(uTrans(i,j))
74     c & *(tracer(i-1,j)-tracer(i,j))*0.5 _d 0
75     C-- above formulation produces large truncation error when:
76     C 1rst.O upWind and u > 0 & |tracer(i-1,j)| << |tracer(i,j)|
77     C or u < 0 & |tracer(i-1,j)| >> |tracer(i,j)|
78     C-- change to a more robust expression:
79     uAbs = ABS(uTrans(i,j))
80     & *( 1. _d 0 - xLimit*(1. _d 0 - uCFL) )
81     uT(i,j) = ( uTrans(i,j)+uAbs )* 0.5 _d 0 * tracer(i-1,j)
82     & + ( uTrans(i,j)-uAbs )* 0.5 _d 0 * tracer(i,j)
83 jmc 1.1 ENDDO
84     ENDDO
85    
86     RETURN
87     END

  ViewVC Help
Powered by ViewVC 1.1.22