/[MITgcm]/MITgcm/pkg/diagnostics/diag_calc_psivel.F
ViewVC logotype

Contents of /MITgcm/pkg/diagnostics/diag_calc_psivel.F

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


Revision 1.4 - (show annotations) (download)
Wed Aug 24 02:24:12 2011 UTC (12 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g, checkpoint63b, checkpoint63c, checkpoint64, checkpoint65, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e
Changes since 1.3: +9 -1 lines
add DEBUG calls

1 C $Header: /u/gcmpack/MITgcm/pkg/diagnostics/diag_calc_psivel.F,v 1.3 2011/08/01 20:40:14 jmc Exp $
2 C $Name: $
3
4 #include "DIAG_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: DIAG_CALC_PSIVEL
8 C !INTERFACE:
9 SUBROUTINE DIAG_CALC_PSIVEL(
10 I k, iPsi0, jPsi0, uTrans, vTrans,
11 O psiVel, psiLoc,
12 I prtMsg, myTime, myIter, myThid )
13 C !DESCRIPTION: \bv
14 C *==========================================================*
15 C | SUBROUTINE DIAG_CALC_PSIVEL
16 C | o Calculate horizontal transport stream-function
17 C | from non-divergent horizontal transport.
18 C *==========================================================*
19 C *==========================================================*
20 C \ev
21
22 C !USES:
23 IMPLICIT NONE
24 C === Global data ===
25 #include "SIZE.h"
26 #include "EEPARAMS.h"
27 #include "PARAMS.h"
28 #ifdef ALLOW_OBCS
29 # include "GRID.h"
30 # include "OBCS_GRID.h"
31 #endif /* ALLOW_OBCS */
32
33 C !INPUT/OUTPUT PARAMETERS:
34 C === Routine arguments ===
35 C k :: current level
36 C i,jPsi0 :: indices of grid-point location where Psi == 0
37 C uTrans :: horizontal transport, u-component
38 C vTrans :: horizontal transport, u-component
39 C psiVel :: horizontal stream-function
40 C psiLoc :: horizontal stream-function at special location
41 C prtMsg :: do print message to standard-output
42 C myTime :: current time of simulation (s)
43 C myIter :: current iteration number
44 C myThid :: my Thread Id number
45 INTEGER k
46 INTEGER iPsi0(nSx,nSy)
47 INTEGER jPsi0(nSx,nSy)
48 _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
49 _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
50 _RL psiVel(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
51 _RL psiLoc(2)
52 LOGICAL prtMsg
53 _RL myTime
54 INTEGER myIter
55 INTEGER myThid
56
57 C !LOCAL VARIABLES:
58 C === Local variables ====
59 C bi, bj :: tile indices
60 C i, j :: loop indices
61 C dPsiX :: tile stream-function increment along X-dir
62 C dPsiY :: tile stream-function increment along Y-dir
63 C psiOri :: stream-function value at tile origin
64 INTEGER bi, bj
65 INTEGER i, j
66 _RL dPsiX (nSx,nSy)
67 _RL dPsiY (nSx,nSy)
68 _RL psiOri(nSx,nSy)
69 _RL offSet
70 LOGICAL zeroPsi
71 c CHARACTER*(MAX_LEN_MBUF) msgBuf
72 #ifdef ALLOW_OBCS
73 INTEGER is, js, ix, jx, iy, jy, ijCnt
74 INTEGER npass, nPts, prev_nPts
75 LOGICAL kPsi(1:sNx+1,1:sNy+1)
76 #endif /* ALLOW_OBCS */
77 CEOP
78
79 #ifdef ALLOW_DEBUG
80 IF (debugMode) CALL DEBUG_ENTER('DIAG_CALC_PSIVEL',myThid)
81 #endif
82
83 C-- Initialise
84 zeroPsi = iPsi0(1,1).GE.0
85 psiLoc(1) = 0.
86 psiLoc(2) = 0.
87
88 C-- step.1 : compute Psi over each tile separately
89 DO bj=myByLo(myThid),myByHi(myThid)
90 DO bi=myBxLo(myThid),myBxHi(myThid)
91 dPsiX (bi,bj) = 0.
92 dPsiY (bi,bj) = 0.
93 psiOri(bi,bj) = 0.
94 IF ( useOBCS ) THEN
95 C- Case with OBC:
96 #ifdef ALLOW_OBCS
97 C- Note: OB may introduce discontinuity in domain & tile stream-function map;
98 C within a tile: define a local "is-set" mask (=kPsi) and
99 C propagate stream-function value without assumption.
100 C between tiles: present code is not "general", likely to work with
101 C simple OB setting and/or simple tile connection (no exch2).
102 C A truly general algorithm requires to change CUMULSUM_Z_TILE (adding 1
103 C more input dPsi/tile) and to account for disabled tile-connection due
104 C to OB when setting cumsum tile-mapping (W2_SET_MAP_CUMSUM).
105 DO j=1,sNy+1
106 DO i=1,sNx+1
107 kPsi(i,j) = .FALSE.
108 psiVel(i,j,bi,bj) = 0.
109 ENDDO
110 ENDDO
111 C- select starting point
112 ijCnt = sNx+sNy+1
113 is = 0
114 js = 0
115 DO j=1,sNy
116 DO i=1,sNx
117 IF ( OBCS_insideMask(i,j,bi,bj).EQ.1.
118 & .AND. (i+j).LE.ijCnt ) THEN
119 is = i
120 js = j
121 ijCnt = i+j
122 ENDIF
123 ENDDO
124 ENDDO
125 IF ( is.EQ.0 ) THEN
126 nPts = 0
127 ELSE
128 kPsi(is,js) = .TRUE.
129 nPts = 1
130 ENDIF
131 npass = 0
132 prev_nPts = 0
133 DO WHILE ( nPts.GT.prev_nPts )
134 prev_nPts = nPts
135 npass = npass + 1
136 DO j=1,sNy+1
137 DO i=1,sNx
138 IF ( OBCS_insideMask(i,j-1,bi,bj).EQ.1. .OR.
139 & OBCS_insideMask(i, j ,bi,bj).EQ.1. ) THEN
140 IF ( kPsi(i,j) .AND. .NOT.kPsi(i+1,j) ) THEN
141 nPts = nPts + 1
142 kPsi(i+1,j) = .TRUE.
143 psiVel(i+1,j,bi,bj) = psiVel(i,j,bi,bj)
144 & + vTrans(i,j,bi,bj)
145 ENDIF
146 IF ( .NOT.kPsi(i,j) .AND. kPsi(i+1,j) ) THEN
147 nPts = nPts + 1
148 kPsi(i,j) = .TRUE.
149 psiVel(i,j,bi,bj) = psiVel(i+1,j,bi,bj)
150 & - vTrans(i,j,bi,bj)
151 ENDIF
152 ENDIF
153 ENDDO
154 ENDDO
155 DO j=1,sNy
156 DO i=1,sNx+1
157 IF ( OBCS_insideMask(i-1,j,bi,bj).EQ.1. .OR.
158 & OBCS_insideMask( i ,j,bi,bj).EQ.1. ) THEN
159 IF ( kPsi(i,j) .AND. .NOT.kPsi(i,j+1) ) THEN
160 nPts = nPts + 1
161 kPsi(i,j+1) = .TRUE.
162 psiVel(i,j+1,bi,bj) = psiVel(i,j,bi,bj)
163 & - uTrans(i,j,bi,bj)
164 ENDIF
165 IF ( .NOT.kPsi(i,j) .AND. kPsi(i,j+1) ) THEN
166 nPts = nPts + 1
167 kPsi(i,j) = .TRUE.
168 psiVel(i,j,bi,bj) = psiVel(i,j+1,bi,bj)
169 & + uTrans(i,j,bi,bj)
170 ENDIF
171 ENDIF
172 ENDDO
173 ENDDO
174 IF ( prtMsg .AND. nPts.GT.prev_nPts ) THEN
175 _BEGIN_MASTER( myThid )
176 WRITE(standardMessageUnit,'(A,2I4,A,I6,I8)')
177 & ' diag_calc_psivel: bi,bj=', bi, bj,
178 & ' : npass,nPts=', npass, nPts
179 _END_MASTER( myThid )
180 ENDIF
181 C- end do while (npass count)
182 ENDDO
183 C- set tile increments
184 ix = 0
185 jx = 0
186 DO i=sNx+1,1,-1
187 DO j=1,sNy
188 IF ( kPsi(i,j) .AND. jx.EQ.0 ) THEN
189 ix = i
190 jx = j
191 ENDIF
192 ENDDO
193 ENDDO
194 IF ( jx.NE.0 ) dPsiX (bi,bj) = psiVel(ix,jx,bi,bj)
195 iy = 0
196 jy = 0
197 DO j=sNy+1,1,-1
198 DO i=1,sNx
199 IF ( kPsi(i,j) .AND. iy.EQ.0 ) THEN
200 iy = i
201 jy = j
202 ENDIF
203 ENDDO
204 ENDDO
205 IF ( iy.NE.0 ) dPsiY (bi,bj) = psiVel(iy,jy,bi,bj)
206 IF ( prtMsg ) THEN
207 _BEGIN_MASTER( myThid )
208 WRITE(standardMessageUnit,'(3(A,2I4))')
209 & ' diag_calc_psivel: is,js=', is,js,
210 & ' ; ix,jx=', ix,jx, ' ; iy,jy=', iy,jy
211 _END_MASTER( myThid )
212 ENDIF
213 #endif /* ALLOW_OBCS */
214 ELSE
215 C- Case without OBC:
216 psiVel(1,1,bi,bj) = psiOri(bi,bj)
217 j = 1
218 DO i=1,sNx
219 psiVel(i+1,j,bi,bj) = psiVel(i,j,bi,bj)
220 & + vTrans(i,j,bi,bj)
221 ENDDO
222 C- note: can vectorise inner loop
223 DO j=1,sNy
224 DO i=1,sNx+1
225 psiVel(i,j+1,bi,bj) = psiVel(i,j,bi,bj)
226 & - uTrans(i,j,bi,bj)
227 ENDDO
228 ENDDO
229 dPsiX (bi,bj) = psiVel(1+sNx,1,bi,bj)
230 dPsiY (bi,bj) = psiVel(1,1+sNy,bi,bj)
231 C- end with/without OBC cases
232 ENDIF
233 ENDDO
234 ENDDO
235
236 CALL CUMULSUM_Z_TILE_RL(
237 O psiOri, psiLoc,
238 I dPsiX, dPsiY, myThid )
239
240 C-- step.2 : account for Psi @ tile origin
241 offSet = 0.
242 DO bj=myByLo(myThid),myByHi(myThid)
243 DO bi=myBxLo(myThid),myBxHi(myThid)
244 DO j=1,sNy+1
245 DO i=1,sNx+1
246 psiVel(i,j,bi,bj) = psiVel(i,j,bi,bj) + psiOri(bi,bj)
247 ENDDO
248 ENDDO
249 IF ( iPsi0(bi,bj)*jPsi0(bi,bj).GT.0 )
250 & offSet = -psiVel(iPsi0(bi,bj),jPsi0(bi,bj),bi,bj)
251 ENDDO
252 ENDDO
253
254 IF ( zeroPsi ) THEN
255 C-- step.3 : shift stream-function to satisfy Psi == 0 @ a particular location
256 _GLOBAL_SUM_RL( offSet, myThid )
257 psiLoc(1) = psiLoc(1) + offSet
258 psiLoc(2) = psiLoc(2) + offSet
259 DO bj=myByLo(myThid),myByHi(myThid)
260 DO bi=myBxLo(myThid),myBxHi(myThid)
261 DO j=1,sNy+1
262 DO i=1,sNx+1
263 psiVel(i,j,bi,bj) = psiVel(i,j,bi,bj) + offSet
264 ENDDO
265 ENDDO
266 ENDDO
267 ENDDO
268 ENDIF
269
270 #ifdef ALLOW_DEBUG
271 IF (debugMode) CALL DEBUG_LEAVE('DIAG_CALC_PSIVEL',myThid)
272 #endif
273
274 RETURN
275 END

  ViewVC Help
Powered by ViewVC 1.1.22