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

Annotation of /MITgcm/pkg/gmredi/gmredi_calc_psi_b.F

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


Revision 1.10 - (hide annotations) (download)
Fri May 30 02:50:16 2008 UTC (16 years ago) by gforget
Branch: MAIN
CVS Tags: checkpoint62c, checkpoint62q, checkpoint62p, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint60, checkpoint61, checkpoint62, checkpoint62b, checkpoint61f, checkpoint61n, checkpoint61q, checkpoint61e, checkpoint61g, checkpoint61d, checkpoint61b, checkpoint61c, checkpoint61a, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61z, checkpoint61x, checkpoint61y
Changes since 1.9: +8 -1 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 gforget 1.10 C $Header: /u/gcmpack/MITgcm/pkg/gmredi/gmredi_calc_psi_b.F,v 1.9 2008/03/28 18:48:05 heimbach Exp $
2 jmc 1.1 C $Name: $
3    
4     #include "GMREDI_OPTIONS.h"
5    
6     CStartOfInterface
7     SUBROUTINE GMREDI_CALC_PSI_B(
8     I bi, bj, iMin, iMax, jMin, jMax,
9     I sigmaX, sigmaY, sigmaR,
10 jmc 1.6 I ldd97_LrhoW, ldd97_LrhoS,
11 jmc 1.1 I myThid )
12     C /==========================================================\
13     C | SUBROUTINE GMREDI_CALC_PSI_B |
14     C | o Calculate stream-functions for GM bolus velocity |
15     C |==========================================================|
16     C \==========================================================/
17     IMPLICIT NONE
18    
19     C == Global variables ==
20     #include "SIZE.h"
21     #include "GRID.h"
22     #include "DYNVARS.h"
23     #include "EEPARAMS.h"
24     #include "PARAMS.h"
25     #include "GMREDI.h"
26 gforget 1.10 #include "FFIELDS.h"
27 jmc 1.1
28 heimbach 1.2 #ifdef ALLOW_AUTODIFF_TAMC
29     #include "tamc.h"
30     #include "tamc_keys.h"
31     #endif /* ALLOW_AUTODIFF_TAMC */
32    
33 jmc 1.1 C == Routine arguments ==
34     C
35     _RL sigmaX(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
36     _RL sigmaY(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
37     _RL sigmaR(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
38 jmc 1.6 _RL ldd97_LrhoW(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
39     _RL ldd97_LrhoS(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
40 jmc 1.1 INTEGER bi,bj,iMin,iMax,jMin,jMax
41     INTEGER myThid
42     CEndOfInterface
43    
44     #ifdef ALLOW_GMREDI
45     #ifdef GM_BOLUS_ADVEC
46    
47     C == Local variables ==
48     INTEGER i,j,k, km1
49     _RL SlopeX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
50     _RL SlopeY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
51     _RL dSigmaDrW(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
52     _RL dSigmaDrS(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
53     _RL taperX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
54     _RL taperY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
55    
56     C- Initialization : <= done in S/R gmredi_init
57    
58 heimbach 1.2 #ifdef ALLOW_AUTODIFF_TAMC
59     act1 = bi - myBxLo(myThid)
60     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
61     act2 = bj - myByLo(myThid)
62     max2 = myByHi(myThid) - myByLo(myThid) + 1
63     act3 = myThid - 1
64     max3 = nTx*nTy
65     act4 = ikey_dynamics - 1
66 heimbach 1.5 igmkey = (act1 + 1) + act2*max1
67     & + act3*max1*max2
68     & + act4*max1*max2*max3
69 heimbach 1.2 #endif /* ALLOW_AUTODIFF_TAMC */
70    
71 heimbach 1.3 #ifdef ALLOW_AUTODIFF_TAMC
72     # ifdef GM_VISBECK_VARIABLE_K
73 heimbach 1.5 CADJ STORE VisbeckK(:,:,bi,bj) = comlev1_bibj, key=igmkey, byte=isbyte
74 heimbach 1.3 # endif
75     #endif
76 jmc 1.1 IF (GM_AdvForm) THEN
77     DO k=2,Nr
78 heimbach 1.2 km1 = k-1
79    
80     #ifdef ALLOW_AUTODIFF_TAMC
81 heimbach 1.5 kkey = (igmkey-1)*Nr + k
82 heimbach 1.2 DO j=1-Oly,sNy+Oly
83     DO i=1-Olx,sNx+Olx
84     SlopeX(i,j) = 0. _d 0
85     SlopeY(i,j) = 0. _d 0
86     dSigmaDrW(i,j) = 0. _d 0
87     dSigmaDrS(i,j) = 0. _d 0
88     ENDDO
89     ENDDO
90     #endif
91 jmc 1.1
92 jmc 1.6 C Gradient of Sigma below U and V points
93     DO j=1-Oly,sNy+Oly
94     DO i=1-Olx+1,sNx+Olx
95     SlopeX(i,j)=op5*( sigmaX(i,j,km1)+sigmaX(i,j,k) )
96     & *maskW(i,j,k,bi,bj)
97     dSigmaDrW(i,j)=op5*( sigmaR(i-1,j,k)+sigmaR(i,j,k) )
98     & *maskW(i,j,k,bi,bj)
99     ENDDO
100     ENDDO
101     DO j=1-Oly+1,sNy+Oly
102     DO i=1-Olx,sNx+Olx
103     SlopeY(i,j)=op5*( sigmaY(i,j,km1)+sigmaY(i,j,k) )
104     & *maskS(i,j,k,bi,bj)
105     dSigmaDrS(i,j)=op5*( sigmaR(i,j-1,k)+sigmaR(i,j,k) )
106     & *maskS(i,j,k,bi,bj)
107     ENDDO
108     ENDDO
109 jmc 1.1
110 jmc 1.6 C Calculate slopes , taper and/or clip
111     CALL GMREDI_SLOPE_PSI(
112     O taperX, taperY,
113 jmc 1.1 U SlopeX, SlopeY,
114 jmc 1.6 U dSigmaDrW, dSigmaDrS,
115     I ldd97_LrhoW, ldd97_LrhoS, rF(k), k,
116 jmc 1.1 I bi, bj, myThid )
117 heimbach 1.2
118     #ifdef ALLOW_AUTODIFF_TAMC
119     CADJ STORE SlopeX(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
120     CADJ STORE SlopeY(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
121     CADJ STORE taperX(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
122     CADJ STORE taperY(:,:) = comlev1_bibj_k, key=kkey, byte=isbyte
123     #endif /* ALLOW_AUTODIFF_TAMC */
124 jmc 1.1
125     C- Compute the 2 stream-function Components ( GM bolus vel.)
126 jmc 1.6 DO j=1-Oly,sNy+Oly
127     DO i=1-Olx+1,sNx+Olx
128     GM_PsiX(i,j,k,bi,bj) = SlopeX(i,j)*taperX(i,j)
129 heimbach 1.9 #if (defined (ALLOW_AUTODIFF) && defined (ALLOW_KAPGM_CONTROL))
130 heimbach 1.8 & *( kapgm(i,j,k,bi,bj)
131     #else
132 jmc 1.6 & *( GM_background_K
133 heimbach 1.8 #endif
134 jmc 1.1 #ifdef GM_VISBECK_VARIABLE_K
135 jmc 1.6 & +op5*(VisbeckK(i-1,j,bi,bj)+VisbeckK(i,j,bi,bj))
136 jmc 1.1 #endif
137 jmc 1.6 & )*maskW(i,j,k,bi,bj)
138 gforget 1.10 #ifdef ALLOW_EDDYPSI
139     & +eddyPsiX(i,j,k,bi,bj)*maskW(i,j,k,bi,bj)
140     #endif
141 jmc 1.6 ENDDO
142     ENDDO
143     DO j=1-Oly+1,sNy+Oly
144     DO i=1-Olx,sNx+Olx
145     GM_PsiY(i,j,k,bi,bj) = SlopeY(i,j)*taperY(i,j)
146 heimbach 1.9 #if (defined (ALLOW_AUTODIFF) && defined (ALLOW_KAPGM_CONTROL))
147 heimbach 1.8 & *( kapgm(i,j,k,bi,bj)
148     #else
149 jmc 1.6 & *( GM_background_K
150 heimbach 1.8 #endif
151 jmc 1.1 #ifdef GM_VISBECK_VARIABLE_K
152 jmc 1.6 & +op5*(VisbeckK(i,j-1,bi,bj)+VisbeckK(i,j,bi,bj))
153 jmc 1.1 #endif
154 jmc 1.6 & )*maskS(i,j,k,bi,bj)
155 gforget 1.10 #ifdef ALLOW_EDDYPSI
156     & +eddyPsiY(i,j,k,bi,bj)*maskS(i,j,k,bi,bj)
157     #endif
158 jmc 1.6 ENDDO
159     ENDDO
160 jmc 1.1
161 jmc 1.6 C----- end of loop on level k
162 jmc 1.1 ENDDO
163    
164     ENDIF
165     #endif /* GM_BOLUS_ADVEC */
166     #endif /* ALLOW_GMREDI */
167    
168     RETURN
169     END

  ViewVC Help
Powered by ViewVC 1.1.22