1 |
mlosch |
1.8 |
C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_dst2u1_adv_x.F,v 1.7 2007/04/04 01:39:06 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 |
mlosch |
1.8 |
ENDDO |
66 |
|
|
DO j=1-Oly,sNy+Oly |
67 |
jmc |
1.1 |
DO i=1-Olx+1,sNx+Olx |
68 |
|
|
|
69 |
jmc |
1.7 |
uCFL = uFld(i,j) |
70 |
|
|
IF ( calcCFL ) uCFL = ABS( uFld(i,j)*deltaTloc |
71 |
jmc |
1.6 |
& *recip_dxC(i,j,bi,bj)*recip_deepFacC(k) ) |
72 |
jmc |
1.1 |
|
73 |
jmc |
1.3 |
c uT(i,j) = |
74 |
|
|
c & uTrans(i,j)*(tracer(i-1,j)+tracer(i,j))*0.5 _d 0 |
75 |
|
|
c & + ( 1. _d 0 - xLimit*(1. _d 0 - uCFL) )*ABS(uTrans(i,j)) |
76 |
|
|
c & *(tracer(i-1,j)-tracer(i,j))*0.5 _d 0 |
77 |
|
|
C-- above formulation produces large truncation error when: |
78 |
|
|
C 1rst.O upWind and u > 0 & |tracer(i-1,j)| << |tracer(i,j)| |
79 |
|
|
C or u < 0 & |tracer(i-1,j)| >> |tracer(i,j)| |
80 |
|
|
C-- change to a more robust expression: |
81 |
|
|
uAbs = ABS(uTrans(i,j)) |
82 |
|
|
& *( 1. _d 0 - xLimit*(1. _d 0 - uCFL) ) |
83 |
|
|
uT(i,j) = ( uTrans(i,j)+uAbs )* 0.5 _d 0 * tracer(i-1,j) |
84 |
|
|
& + ( uTrans(i,j)-uAbs )* 0.5 _d 0 * tracer(i,j) |
85 |
jmc |
1.1 |
ENDDO |
86 |
|
|
ENDDO |
87 |
|
|
|
88 |
|
|
RETURN |
89 |
|
|
END |