/[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.1 - (hide annotations) (download)
Sat Oct 22 19:56:33 2005 UTC (18 years, 7 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58e_post, checkpoint57y_post, checkpoint57y_pre, checkpoint58, checkpoint58f_post, checkpoint57x_post, checkpoint58d_post, checkpoint58c_post, checkpoint57w_post, checkpoint58a_post, checkpoint58g_post, checkpoint57z_post, checkpoint58b_post
add Lax-Wendroff (=DST2) & 1rst order upwind advection scheme.

1 jmc 1.1 C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_fluxlimit_adv_x.F,v 1.6 2005/08/19 22:19:35 heimbach 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, uVel,
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 uVel :: 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 uVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
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 uFld :: velocity [m/s], zonal component
53     C uCFL :: Courant-Friedrich-Levy number
54     INTEGER i,j
55     _RL uFld, uCFL, xLimit
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     c uFld = uVel(i,j,k,bi,bj)
66     uFld = uTrans(i,j)*recip_dyG(i,j,bi,bj)
67     & *recip_drF(k)*recip_hFacW(i,j,k,bi,bj)
68     uCFL = ABS(uFld*deltaTloc*recip_dxC(i,j,bi,bj))
69    
70     uT(i,j) =
71     & uTrans(i,j)*(Tracer(i-1,j)+Tracer(i,j))*0.5 _d 0
72     & + ( 1. _d 0 - xLimit*(1. _d 0 - uCFL) )*ABS(uTrans(i,j))
73     & *(tracer(i-1,j)-tracer(i,j))*0.5 _d 0
74     ENDDO
75     ENDDO
76    
77     RETURN
78     END

  ViewVC Help
Powered by ViewVC 1.1.22