/[MITgcm]/MITgcm/pkg/gmredi/gmredi_calc_uvflow.F
ViewVC logotype

Contents of /MITgcm/pkg/gmredi/gmredi_calc_uvflow.F

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


Revision 1.5 - (show annotations) (download)
Fri May 30 02:50:16 2008 UTC (15 years, 11 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62c, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62w, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint63i, checkpoint63g, checkpoint64, checkpoint60, checkpoint61, checkpoint62, checkpoint63, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint62b, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint62d, checkpoint61f, checkpoint61n, checkpoint61q, checkpoint61z, checkpoint61e, checkpoint61g, checkpoint61d, checkpoint61b, checkpoint61c, checkpoint61a, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61x, checkpoint61y
Changes since 1.4: +3 -2 lines
o bridging the gap between eddy stress and GM.
  -> eddyTau is replaced with eddyPsi (eddyTau = f x rho0 x eddyPsi)
      along with a change in CPP option (now ALLOW_EDDYPSI).
  -> when using GM w/ GM_AdvForm:
      The total eddy streamfunction (Psi = eddyPsi + K x Slope)
      is applied either in the tracer Eq. or in momentum Eq.
      depending on data.gmredi (intro. GM_InMomAsStress).
  -> ALLOW_EDDYPSI_CONTROL for estimation purpose.
  The key modifications are in model/src/taueddy_external_forcing.F
  pkg/gmredi/gmredi_calc_*F pkg/gmredi/gmredi_*transport.F

1 C $Header: /u/gcmpack/MITgcm/pkg/gmredi/gmredi_calc_uvflow.F,v 1.4 2006/06/18 23:28:17 jmc Exp $
2 C $Name: $
3
4 #include "GMREDI_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: GMREDI_CALC_UVFLOW
8 C !INTERFACE:
9 SUBROUTINE GMREDI_CALC_UVFLOW(
10 U uFld, vFld, uTrans, vTrans,
11 I k, bi, bj, myThid )
12 C !DESCRIPTION:
13 C Add GM-bolus flow to Eulerian horizontal transport.
14
15 C !USES:
16 IMPLICIT NONE
17
18 C == GLobal variables ==
19 #include "SIZE.h"
20 #include "EEPARAMS.h"
21 #include "PARAMS.h"
22 #include "GRID.h"
23 #include "GMREDI.h"
24
25 C !INPUT/OUTPUT PARAMETERS:
26 C == Routine arguments ==
27 C uFld :: zonal velocity (updated)
28 C vFld :: meridional velocity (updated)
29 C uTrans :: zonal volume transport (updated)
30 C vTrans :: meridional volume transport (updated)
31 C k :: level index
32 C bi,bj :: tile indices
33 C myThid :: thread number
34 INTEGER k, bi, bj, myThid
35 _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
36 _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
37 _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
38 _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
39
40 #ifdef ALLOW_GMREDI
41 #ifdef GM_BOLUS_ADVEC
42
43 C !LOCAL VARIABLES:
44 C == Local variables ==
45 C i, j :: loop indices
46 INTEGER i, j
47 INTEGER kp1
48 _RL maskp1
49 _RL delPsi
50 CEOP
51
52 IF (GM_AdvForm .AND. .NOT.GM_AdvSeparate
53 & .AND. .NOT.GM_InMomAsStress) THEN
54
55 kp1 = MIN(k+1,Nr)
56 maskp1 = 1.
57 IF (k.GE.Nr) maskp1 = 0.
58 DO j=1-OLy,sNy+OLy
59 DO i=1-OLx,sNx+OLx
60 delPsi = GM_PsiX(i,j,kp1,bi,bj)*maskp1
61 & - GM_PsiX(i,j, k, bi,bj)
62 uFld(i,j) = uFld(i,j)
63 & + delPsi*recip_drF(k)*_recip_hFacW(i,j,k,bi,bj)
64 uTrans(i,j) = uTrans(i,j)
65 & + dyG(i,j,bi,bj)*delPsi*maskW(i,j,k,bi,bj)
66 ENDDO
67 ENDDO
68 DO j=1-OLy,sNy+OLy
69 DO i=1-OLx,sNx+OLx
70 delPsi = GM_PsiY(i,j,kp1,bi,bj)*maskp1
71 & - GM_PsiY(i,j, k, bi,bj)
72 vFld(i,j) = vFld(i,j)
73 & + delPsi*recip_drF(k)*_recip_hFacS(i,j,k,bi,bj)
74 vTrans(i,j) = vTrans(i,j)
75 & + dxG(i,j,bi,bj)*delPsi*maskS(i,j,k,bi,bj)
76 ENDDO
77 ENDDO
78
79 ENDIF
80 #endif /* GM_BOLUS_ADVEC */
81 #endif /* ALLOW_GMREDI */
82
83 RETURN
84 END

  ViewVC Help
Powered by ViewVC 1.1.22