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

Contents of /MITgcm/pkg/gmredi/gmredi_calc_wflow.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, checkpoint64q, checkpoint64p, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, 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_wflow.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_WFLOW
8 C !INTERFACE:
9 SUBROUTINE GMREDI_CALC_WFLOW(
10 U wFld, rTrans,
11 I k, bi, bj, myThid )
12 C !DESCRIPTION:
13 C Add GM-bolus flow to Eulerian vertical 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 wFld :: vertical volume transport (updated)
28 C rTrans :: vertical volume transport (updated)
29 C k :: level index
30 C bi,bj :: tile indices
31 C myThid :: thread number
32 _RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
33 _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
34 INTEGER k, bi, bj, myThid
35
36 #ifdef ALLOW_GMREDI
37 #ifdef GM_BOLUS_ADVEC
38
39 C !LOCAL VARIABLES:
40 C == Local variables ==
41 C i, j :: loop indices
42 INTEGER i, j
43 _RL delPsi
44 CEOP
45
46 IF (GM_AdvForm .AND. .NOT.GM_AdvSeparate
47 & .AND. .NOT.GM_InMomAsStress) THEN
48
49 DO j=1-Oly,sNy+Oly-1
50 DO i=1-Olx,sNx+Olx-1
51 delPsi = ( dyG(i+1,j,bi,bj)*GM_PsiX(i+1,j,k,bi,bj)
52 & -dyG( i ,j,bi,bj)*GM_PsiX( i ,j,k,bi,bj)
53 & +dxG(i,j+1,bi,bj)*GM_PsiY(i,j+1,k,bi,bj)
54 & -dxG(i, j ,bi,bj)*GM_PsiY(i, j ,k,bi,bj)
55 & )*maskC(i,j,k,bi,bj)
56 wFld(i,j) = wFld(i,j) + delPsi*recip_rA(i,j,bi,bj)
57 rTrans(i,j) = rTrans(i,j) + delPsi
58 ENDDO
59 ENDDO
60
61 ENDIF
62 #endif /* GM_BOLUS_ADVEC */
63 #endif /* ALLOW_GMREDI */
64
65 RETURN
66 END

  ViewVC Help
Powered by ViewVC 1.1.22