/[MITgcm]/MITgcm/model/src/taueddy_external_forcing.F
ViewVC logotype

Contents of /MITgcm/model/src/taueddy_external_forcing.F

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


Revision 1.4 - (show annotations) (download)
Sat May 31 16:50:35 2008 UTC (16 years ago) by gforget
Branch: MAIN
CVS Tags: checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64n, 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.3: +6 -2 lines
o bridging the gap between eddy stress and GM.
... had missed this bit.

1 C $Header: /u/gcmpack/MITgcm/model/src/taueddy_external_forcing.F,v 1.3 2008/05/30 02:45:43 gforget Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 #ifdef ALLOW_GMREDI
8 # include "GMREDI_OPTIONS.h"
9 #endif
10
11 CBOP
12 C !ROUTINE: TAUEDDY_EXTERNAL_FORCING_U
13 C !INTERFACE:
14 SUBROUTINE TAUEDDY_EXTERNAL_FORCING_U(
15 I iMin, iMax, jMin, jMax,bi,bj,kLev,
16 I myCurrentTime,myThid)
17 C !DESCRIPTION: \bv
18 C *==========================================================*
19 C | S/R TAUEDDY_EXTERNAL_FORCING_U
20 C | o Contains problem specific forcing for zonal velocity.
21 C *==========================================================*
22 C | Adds terms to gU for forcing by external sources
23 C | e.g. wind stress, bottom friction etc..................
24 C *==========================================================*
25 C \ev
26
27 C !USES:
28 IMPLICIT NONE
29 C == Global data ==
30 #include "SIZE.h"
31 #include "EEPARAMS.h"
32 #include "PARAMS.h"
33 #include "GRID.h"
34 #include "DYNVARS.h"
35 #include "FFIELDS.h"
36 #ifdef ALLOW_GMREDI
37 #include "GMREDI.h"
38 #endif
39
40 C !INPUT/OUTPUT PARAMETERS:
41 C == Routine arguments ==
42 C iMin - Working range of tile for applying forcing.
43 C iMax
44 C jMin
45 C jMax
46 C kLev
47 INTEGER iMin, iMax, jMin, jMax, kLev, bi, bj
48 _RL myCurrentTime
49 INTEGER myThid
50 CEOP
51
52 #ifdef ALLOW_EDDYPSI
53 C !LOCAL VARIABLES:
54 C == Local variables ==
55 C Loop counters
56 INTEGER I, J
57 C number of surface interface layer
58 INTEGER kSurface, Kp1
59 _RL maskm1, maskp1
60
61 IF ( fluidIsAir ) THEN
62 kSurface = 0
63 ELSEIF ( usingPCoords ) THEN
64 kSurface = Nr
65 ELSE
66 kSurface = 1
67 ENDIF
68
69 C Add zonal eddy momentum impulse into the layer
70 #ifdef ALLOW_GMREDI
71 if ( GM_InMomAsStress ) then
72 #endif
73 Kp1=min(klev+1,Nr)
74 maskp1=1.
75 maskm1=1.
76 IF (klev.EQ.Nr) maskp1=0.
77 IF (klev.EQ.1) maskm1=0.
78 DO j=jMin,jMax
79 DO i=iMin,iMax
80 gU(i,j,kLev,bi,bj) = gU(i,j,kLev,bi,bj)
81 & +foFacMom*
82 #ifdef ALLOW_GMREDI
83 & (maskm1*_maskW(i,j,klev,bi,bj)*GM_PsiX(i,j,klev,bi,bj)
84 & -maskp1*_maskW(i,j,kp1,bi,bj)*GM_PsiX(i,j,kp1,bi,bj))
85 #else
86 & (maskm1*_maskW(i,j,klev,bi,bj)*eddyPsiX(i,j,klev,bi,bj)
87 & -maskp1*_maskW(i,j,kp1,bi,bj)*eddyPsiX(i,j,kp1,bi,bj))
88 #endif
89 & *recip_drF(klev)*_recip_hFacW(i,j,klev,bi,bj)
90 & *0.5*(_fCori(i,j,bi,bj)+_fCori(i-1,j,bi,bj))
91 ENDDO
92 ENDDO
93 #ifdef ALLOW_GMREDI
94 endif
95 #endif
96
97 #endif
98
99 RETURN
100 END
101 CBOP
102 C !ROUTINE: TAUEDDY_EXTERNAL_FORCING_V
103 C !INTERFACE:
104 SUBROUTINE TAUEDDY_EXTERNAL_FORCING_V(
105 I iMin, iMax, jMin, jMax,bi,bj,kLev,
106 I myCurrentTime,myThid)
107 C !DESCRIPTION: \bv
108 C *==========================================================*
109 C | S/R TAUEDDY_EXTERNAL_FORCING_V
110 C | o Contains problem specific forcing for merid velocity.
111 C *==========================================================*
112 C | Adds terms to gV for forcing by external sources
113 C | e.g. wind stress, bottom friction etc..................
114 C *==========================================================*
115 C \ev
116
117 C !USES:
118 IMPLICIT NONE
119 C == Global data ==
120 #include "SIZE.h"
121 #include "EEPARAMS.h"
122 #include "PARAMS.h"
123 #include "GRID.h"
124 #include "DYNVARS.h"
125 #include "FFIELDS.h"
126 #ifdef ALLOW_GMREDI
127 #include "GMREDI.h"
128 #endif
129
130 C !INPUT/OUTPUT PARAMETERS:
131 C == Routine arguments ==
132 C iMin - Working range of tile for applying forcing.
133 C iMax
134 C jMin
135 C jMax
136 C kLev
137 INTEGER iMin, iMax, jMin, jMax, kLev, bi, bj
138 _RL myCurrentTime
139 INTEGER myThid
140 CEOP
141
142 #ifdef ALLOW_EDDYPSI
143 C !LOCAL VARIABLES:
144 C == Local variables ==
145 C Loop counters
146 INTEGER I, J
147 C number of surface interface layer
148 INTEGER kSurface, Kp1
149 _RL maskm1, maskp1
150
151 IF ( fluidIsAir ) THEN
152 kSurface = 0
153 ELSEIF ( usingPCoords ) THEN
154 kSurface = Nr
155 ELSE
156 kSurface = 1
157 ENDIF
158
159 C Add meridional eddy momentum impulse into the layer
160 #ifdef ALLOW_GMREDI
161 if ( GM_InMomAsStress ) then
162 #endif
163 Kp1=min(klev+1,Nr)
164 maskp1=1.
165 maskm1=1.
166 IF (klev.EQ.Nr) maskp1=0.
167 IF (klev.EQ.1) maskm1=0.
168 DO j=jMin,jMax
169 DO i=iMin,iMax
170 gV(i,j,kLev,bi,bj) = gV(i,j,kLev,bi,bj)
171 & +foFacMom*
172 #ifdef ALLOW_GMREDI
173 & (maskm1*_maskS(i,j,klev,bi,bj)*GM_PsiY(i,j,klev,bi,bj)
174 & -maskp1*_maskS(i,j,kp1,bi,bj)*GM_PsiY(i,j,kp1,bi,bj))
175 #else
176 & (maskm1*_maskS(i,j,klev,bi,bj)*eddyPsiY(i,j,klev,bi,bj)
177 & -maskp1*_maskS(i,j,kp1,bi,bj)*eddyPsiY(i,j,kp1,bi,bj))
178 #endif
179 & *recip_drF(klev)*_recip_hFacS(i,j,klev,bi,bj)
180 & *0.5*(_fCori(i,j,bi,bj)+_fCori(i,j-1,bi,bj))
181 ENDDO
182 ENDDO
183 #ifdef ALLOW_GMREDI
184 endif
185 #endif
186
187 #endif
188
189 RETURN
190 END

  ViewVC Help
Powered by ViewVC 1.1.22