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

Contents 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.5 - (show annotations) (download)
Mon Jun 19 14:40:43 2006 UTC (17 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58l_post, checkpoint58r_post, checkpoint58n_post, checkpoint58q_post, checkpoint58o_post, checkpoint58k_post, checkpoint58p_post, checkpoint58m_post
Changes since 1.4: +3 -3 lines
DST advection S/R : use local copy of velocity to compute CFL

1 C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_dst2u1_adv_x.F,v 1.3 2006/06/12 16:06:27 jmc Exp $
2 C $Name: $
3
4 #include "GAD_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: GAD_DST2U1_ADV_X
8
9 C !INTERFACE: ==========================================================
10 SUBROUTINE GAD_DST2U1_ADV_X(
11 I bi,bj,k, advectionScheme,
12 I deltaTloc, uTrans, uFld,
13 I tracer,
14 O uT,
15 I myThid )
16
17 C !DESCRIPTION:
18 C Calculates the area integrated zonal flux due to advection
19 C of a tracer using second-order Direct Space and Time (DST-2)
20 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 C advectionScheme :: advection scheme to use: either 2nd Order DST
32 C or 1rst Order Upwind
33 C uTrans :: zonal volume transport
34 C uFld :: zonal flow
35 C tracer :: tracer field
36 C myThid :: thread number
37 INTEGER bi,bj,k
38 INTEGER advectionScheme
39 _RL deltaTloc
40 _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
41 _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
42 _RL tracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
43 INTEGER myThid
44
45 C !OUTPUT PARAMETERS: ==================================================
46 C uT :: zonal advective flux
47 _RL uT (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
48
49 C !LOCAL VARIABLES: ====================================================
50 C i,j :: loop indices
51 C rLimit :: centered (vs upwind) fraction
52 C uLoc :: velocity [m/s], zonal component
53 C uCFL :: Courant-Friedrich-Levy number
54 INTEGER i,j
55 _RL uLoc, uCFL, xLimit, uAbs
56 CEOP
57
58 xLimit = 0. _d 0
59 IF ( advectionScheme.EQ.ENUM_DST2 ) xLimit = 1. _d 0
60
61 DO j=1-Oly,sNy+Oly
62 uT(1-Olx,j)=0.
63 DO i=1-Olx+1,sNx+Olx
64
65 uLoc = uFld(i,j)
66 c uLoc = uTrans(i,j)*recip_dyG(i,j,bi,bj)
67 c & *recip_drF(k)*_recip_hFacW(i,j,k,bi,bj)
68 uCFL = ABS(uLoc*deltaTloc*recip_dxC(i,j,bi,bj))
69
70 c uT(i,j) =
71 c & uTrans(i,j)*(tracer(i-1,j)+tracer(i,j))*0.5 _d 0
72 c & + ( 1. _d 0 - xLimit*(1. _d 0 - uCFL) )*ABS(uTrans(i,j))
73 c & *(tracer(i-1,j)-tracer(i,j))*0.5 _d 0
74 C-- above formulation produces large truncation error when:
75 C 1rst.O upWind and u > 0 & |tracer(i-1,j)| << |tracer(i,j)|
76 C or u < 0 & |tracer(i-1,j)| >> |tracer(i,j)|
77 C-- change to a more robust expression:
78 uAbs = ABS(uTrans(i,j))
79 & *( 1. _d 0 - xLimit*(1. _d 0 - uCFL) )
80 uT(i,j) = ( uTrans(i,j)+uAbs )* 0.5 _d 0 * tracer(i-1,j)
81 & + ( uTrans(i,j)-uAbs )* 0.5 _d 0 * tracer(i,j)
82 ENDDO
83 ENDDO
84
85 RETURN
86 END

  ViewVC Help
Powered by ViewVC 1.1.22