1 |
C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_dst2u1_adv_y.F,v 1.6 2006/12/05 22:21:50 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "GAD_OPTIONS.h" |
5 |
|
6 |
CBOP |
7 |
C !ROUTINE: GAD_DST2U1_ADV_Y |
8 |
|
9 |
C !INTERFACE: ========================================================== |
10 |
SUBROUTINE GAD_DST2U1_ADV_Y( |
11 |
I bi,bj,k, advectionScheme, calcCFL, |
12 |
I deltaTloc, vTrans, vFld, |
13 |
I tracer, |
14 |
O vT, |
15 |
I myThid ) |
16 |
|
17 |
C !DESCRIPTION: |
18 |
C Calculates the area integrated meridional 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 calcCFL :: =T: calculate CFL number ; =F: take vFld as CFL |
34 |
C deltaTloc :: local time-step (s) |
35 |
C vTrans :: meridional volume transport |
36 |
C vFld :: meridional flow / CFL number |
37 |
C tracer :: tracer field |
38 |
C myThid :: thread number |
39 |
INTEGER bi,bj, k, advectionScheme |
40 |
LOGICAL calcCFL |
41 |
_RL deltaTloc |
42 |
_RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
43 |
_RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
44 |
_RL tracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
45 |
INTEGER myThid |
46 |
|
47 |
C !OUTPUT PARAMETERS: ================================================== |
48 |
C vT :: meridional advective flux |
49 |
_RL vT (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
50 |
|
51 |
C !LOCAL VARIABLES: ==================================================== |
52 |
C i,j :: loop indices |
53 |
C yLimit :: centered (vs upwind) fraction |
54 |
C vCFL :: Courant-Friedrich-Levy number |
55 |
INTEGER i,j |
56 |
_RL vCFL, yLimit, vAbs |
57 |
CEOP |
58 |
|
59 |
yLimit = 0. _d 0 |
60 |
IF ( advectionScheme.EQ.ENUM_DST2 ) yLimit = 1. _d 0 |
61 |
|
62 |
DO i=1-Olx,sNx+Olx |
63 |
vT(i,1-Oly)=0. |
64 |
ENDDO |
65 |
DO j=1-Oly+1,sNy+Oly |
66 |
DO i=1-Olx,sNx+Olx |
67 |
|
68 |
vCFL = vFld(i,j) |
69 |
IF ( calcCFL ) vCFL = ABS( vFld(i,j)*deltaTloc |
70 |
& *recip_dyC(i,j,bi,bj)*recip_deepFacC(k) ) |
71 |
|
72 |
c vT(i,j) = |
73 |
c & vTrans(i,j)*(tracer(i,j-1)+tracer(i,j))*0.5 _d 0 |
74 |
c & + ( 1. _d 0 - yLimit*(1. _d 0 - vCFL) )*ABS(vTrans(i,j)) |
75 |
c & *(tracer(i,j-1)-tracer(i,j))*0.5 _d 0 |
76 |
C-- above formulation produces large truncation error when: |
77 |
C 1rst.O upWind and v > 0 & |tracer(i,j-1)| << |tracer(i,j)| |
78 |
C or v < 0 & |tracer(i,j-1)| >> |tracer(i,j)| |
79 |
C-- change to a more robust expression: |
80 |
vAbs = ABS(vTrans(i,j)) |
81 |
& *( 1. _d 0 - yLimit*(1. _d 0 - vCFL) ) |
82 |
vT(i,j) = ( vTrans(i,j)+vAbs )* 0.5 _d 0 * tracer(i,j-1) |
83 |
& + ( vTrans(i,j)-vAbs )* 0.5 _d 0 * tracer(i,j) |
84 |
|
85 |
ENDDO |
86 |
ENDDO |
87 |
|
88 |
RETURN |
89 |
END |