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

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

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


Revision 1.3 - (show annotations) (download)
Fri May 30 02:50:16 2008 UTC (16 years ago) by gforget
Branch: MAIN
CVS Tags: checkpoint64i, checkpoint64h, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint64, checkpoint60, checkpoint61, checkpoint62, checkpoint63, checkpoint62c, checkpoint62b, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62w, checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint61f, checkpoint61g, checkpoint61d, checkpoint61e, checkpoint61b, checkpoint61c, checkpoint61a, checkpoint61n, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61q, checkpoint61z, checkpoint61x, checkpoint61y
Changes since 1.2: +4 -4 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_diagnostics_fill.F,v 1.2 2007/08/20 17:58:57 gforget Exp $
2 C $Name: $
3
4 #include "GMREDI_OPTIONS.h"
5
6 CStartOfInterface
7 SUBROUTINE GMREDI_DIAGNOSTICS_FILL(
8 I bi, bj, myThid )
9 C *==========================================================*
10 C | SUBROUTINE GMREDI_DIAGNOSTICS_FILL
11 C | o fill GM-Redi diagnostics
12 C *==========================================================*
13 C | note: formerly was part of S/R GMREDI_CALC_TENSOR
14 C | and was isolated in this S/R for TAF reasons
15 C *==========================================================*
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 == Routine arguments ==
26 C
27 INTEGER bi,bj
28 INTEGER myThid
29 #ifdef ALLOW_DIAGNOSTICS
30 LOGICAL DIAGNOSTICS_IS_ON
31 EXTERNAL DIAGNOSTICS_IS_ON
32 #endif /* ALLOW_DIAGNOSTICS */
33 CEndOfInterface
34
35 #ifdef ALLOW_GMREDI
36
37 C == Local variables ==
38 #ifdef ALLOW_EDDYPSI
39 INTEGER i,j,k
40 _RL tmpfld3dloc (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
41 #endif
42
43 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
44
45
46 #ifdef ALLOW_DIAGNOSTICS
47 IF ( useDiagnostics ) THEN
48
49 #ifdef GM_VISBECK_VARIABLE_K
50 IF ( GM_Visbeck_alpha.NE.0. ) THEN
51 CALL DIAGNOSTICS_FILL(VisbeckK,'GM_VisbK',0,1,1,bi,bj,myThid)
52 ENDIF
53 #endif
54 #ifdef GM_NON_UNITY_DIAGONAL
55 CALL DIAGNOSTICS_FILL(Kux,'GM_Kux ',0,Nr,1,bi,bj,myThid)
56 CALL DIAGNOSTICS_FILL(Kvy,'GM_Kvy ',0,Nr,1,bi,bj,myThid)
57 #endif
58 #ifdef GM_EXTRA_DIAGONAL
59 IF ( GM_ExtraDiag ) THEN
60 CALL DIAGNOSTICS_FILL(Kuz,'GM_Kuz ',0,Nr,1,bi,bj,myThid)
61 CALL DIAGNOSTICS_FILL(Kvz,'GM_Kvz ',0,Nr,1,bi,bj,myThid)
62 ENDIF
63 #endif
64 CALL DIAGNOSTICS_FILL(Kwx,'GM_Kwx ',0,Nr,1,bi,bj,myThid)
65 CALL DIAGNOSTICS_FILL(Kwy,'GM_Kwy ',0,Nr,1,bi,bj,myThid)
66 CALL DIAGNOSTICS_FILL(Kwz,'GM_Kwz ',0,Nr,1,bi,bj,myThid)
67 #ifdef GM_BOLUS_ADVEC
68 IF ( GM_AdvForm ) THEN
69 CALL DIAGNOSTICS_FILL(GM_PsiX,'GM_PsiX ',0,Nr,1,bi,bj,myThid)
70 CALL DIAGNOSTICS_FILL(GM_PsiY,'GM_PsiY ',0,Nr,1,bi,bj,myThid)
71 ENDIF
72 #endif
73
74 #ifdef ALLOW_EDDYPSI
75 IF ( DIAGNOSTICS_IS_ON('GMEdTauX',myThid) ) THEN
76 DO k = 1, Nr
77 DO j = 1, sny
78 DO i = 1, snx
79 tmpfld3dloc(i,j,k,bi,bj) =
80 & 0.5*rhoConst*fCori(i,j,bi,bj)*
81 & Kwy(i,j,k,bi,bj)
82 ENDDO
83 ENDDO
84 ENDDO
85 CALL DIAGNOSTICS_FILL(tmpfld3dloc,'GMEdTauX',
86 & 0,Nr,1,bi,bj,myThid)
87 ENDIF
88 c
89 IF ( DIAGNOSTICS_IS_ON('GMEdTauY',myThid) ) THEN
90 DO k = 1, Nr
91 DO j = 1, sny
92 DO i = 1, snx
93 tmpfld3dloc(i,j,k,bi,bj) =
94 & -0.5*rhoConst*fCori(i,j,bi,bj)*
95 & Kwx(i,j,k,bi,bj)
96 ENDDO
97 ENDDO
98 ENDDO
99 CALL DIAGNOSTICS_FILL(tmpfld3dloc,'GMEdTauY',
100 & 0,Nr,1,bi,bj,myThid)
101 ENDIF
102 #endif /* ALLOW_EDDYPSI */
103
104 ENDIF
105 #endif /* ALLOW_DIAGNOSTICS */
106
107 #endif /* ALLOW_GMREDI */
108
109 RETURN
110 END

  ViewVC Help
Powered by ViewVC 1.1.22