/[MITgcm]/MITgcm/pkg/obcs/obcs_calc.F
ViewVC logotype

Contents of /MITgcm/pkg/obcs/obcs_calc.F

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


Revision 1.34 - (show annotations) (download)
Thu Nov 15 15:55:42 2012 UTC (11 years, 7 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint65, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65o, 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, HEAD
Changes since 1.33: +7 -1 lines
adding tidal velocity forcing capability to obcs
 Modified Files:
  model/src/dynamics.F forward_step.F the_main_loop.F
  pkg/obcs/OBCS_FIELDS.h OBCS_OPTIONS.h OBCS_PARAMS.h
  OBCS_SEAICE.h obcs_apply_eta.F obcs_apply_r_star.F
  obcs_apply_surf_dr.F obcs_apply_ts.F obcs_apply_w.F
  obcs_calc.F obcs_check.F obcs_init_variables.F obcs_readparms.F
  verification/seaice_obcs/code/OBCS_OPTIONS.h
 Added Files:
  pkg/obcs/obcs_add_tides.F
  verification/seaice_obcs/input.tides/*
  verification/seaice_obcs/results/output.tides.txt

1 C $Header: /u/gcmpack/MITgcm/pkg/obcs/obcs_calc.F,v 1.33 2012/09/18 20:09:17 jmc Exp $
2 C $Name: $
3
4 #include "OBCS_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: OBCS_CALC
8
9 C !INTERFACE:
10 SUBROUTINE OBCS_CALC( futureTime, futureIter,
11 & uVel, vVel, wVel, theta, salt,
12 & myThid )
13
14 C !DESCRIPTION:
15 C *==========================================================*
16 C | SUBROUTINE OBCS_CALC
17 C | o Calculate future boundary data at open boundaries
18 C | at time = futureTime
19 C *==========================================================*
20
21 C !USES:
22 IMPLICIT NONE
23
24 C === Global variables ===
25 #include "SIZE.h"
26 #include "EEPARAMS.h"
27 #include "PARAMS.h"
28 #include "GRID.h"
29 #include "OBCS_PARAMS.h"
30 #include "OBCS_GRID.h"
31 #include "OBCS_FIELDS.h"
32 #ifdef ALLOW_PTRACERS
33 # include "PTRACERS_SIZE.h"
34 # include "PTRACERS_PARAMS.h"
35 # include "PTRACERS_FIELDS.h"
36 # include "OBCS_PTRACERS.h"
37 #endif /* ALLOW_PTRACERS */
38 #ifdef ALLOW_NEST_CHILD
39 # include "NEST_CHILD.h"
40 #endif /* ALLOW_NEST_CHILD */
41
42 C !INPUT/OUTPUT PARAMETERS:
43 C == Routine arguments ==
44 INTEGER futureIter
45 _RL futureTime
46 _RL uVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
47 _RL vVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
48 _RL wVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
49 _RL theta(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
50 _RL salt (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
51 INTEGER myThid
52
53 #ifdef ALLOW_OBCS
54
55 C !LOCAL VARIABLES:
56 C bi, bj :: tile indices
57 C i,j,k :: loop indices
58 C I_obc, J_obc :: local index of open boundary
59 C msgBuf :: Informational/error message buffer
60 INTEGER bi, bj
61 INTEGER i, j, k
62 #ifdef ALLOW_PTRACERS
63 INTEGER I_obc, J_obc
64 CHARACTER*(MAX_LEN_MBUF) msgBuf
65 INTEGER iTracer
66 #endif /* ALLOW_PTRACERS */
67
68 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
69
70 #ifdef ALLOW_DEBUG
71 IF (debugMode) CALL DEBUG_ENTER('OBCS_CALC',myThid)
72 #endif
73
74 DO bj=myByLo(myThid),myByHi(myThid)
75 DO bi=myBxLo(myThid),myBxHi(myThid)
76
77 #ifdef ALLOW_NEST_CHILD
78 IF ( useNEST_CHILD ) THEN
79 IF ( PASSI.LT.2 ) THEN
80 CALL NEST_CHILD_RECV ( myThid )
81 ENDIF
82 ENDIF
83 #endif /* ALLOW_NEST_CHILD */
84
85 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
86
87 #ifdef ALLOW_OBCS_EAST
88 C Eastern OB
89 #ifdef ALLOW_DEBUG
90 IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: East',myThid)
91 #endif
92 IF (useOrlanskiEast) THEN
93 #ifdef ALLOW_ORLANSKI
94 CALL ORLANSKI_EAST(
95 & bi, bj, futureTime,
96 & uVel, vVel, wVel, theta, salt,
97 & myThid )
98 #endif
99 #ifdef ALLOW_NEST_CHILD
100 ELSEIF ( useNEST_CHILD ) THEN
101 DO k=1,Nr
102 DO j=1-OLy,sNy+OLy
103 IF ( OB_Ie(j,bi,bj).NE.OB_indexNone ) THEN
104 OBEu(j,k,bi,bj)= U_F1(j,k,2)
105 OBEv(j,k,bi,bj)= V_F1(j,k,2)
106 OBEt(j,k,bi,bj)= T_F1(j,k,2)
107 OBEs(j,k,bi,bj)= S_F1(j,k,2)
108 #ifdef NONLIN_FRSURF
109 OBEeta(j,bi,bj)= ETA_F1(j,1,2)
110 #endif
111 ENDIF
112 ENDDO
113 ENDDO
114 #endif /* ALLOW_NEST_CHILD */
115 ELSE
116 DO k=1,Nr
117 DO j=1-OLy,sNy+OLy
118 IF ( OB_Ie(j,bi,bj).NE.OB_indexNone ) THEN
119 OBEu(j,k,bi,bj)=0.
120 OBEv(j,k,bi,bj)=0.
121 OBEt(j,k,bi,bj)=tRef(k)
122 OBEs(j,k,bi,bj)=sRef(k)
123 #ifdef ALLOW_NONHYDROSTATIC
124 OBEw(j,k,bi,bj)=0.
125 #endif
126 #ifdef NONLIN_FRSURF
127 OBEeta(j,bi,bj)=0.
128 #endif
129 ENDIF
130 ENDDO
131 ENDDO
132 ENDIF
133 #endif /* ALLOW_OBCS_EAST */
134
135 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
136
137 #ifdef ALLOW_OBCS_WEST
138 C Western OB
139 #ifdef ALLOW_DEBUG
140 IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: West',myThid)
141 #endif
142 IF (useOrlanskiWest) THEN
143 #ifdef ALLOW_ORLANSKI
144 CALL ORLANSKI_WEST(
145 & bi, bj, futureTime,
146 & uVel, vVel, wVel, theta, salt,
147 & myThid )
148 #endif
149 #ifdef ALLOW_NEST_CHILD
150 ELSEIF ( useNEST_CHILD ) THEN
151 DO k=1,Nr
152 DO j=1-OLy,sNy+OLy
153 IF ( OB_Iw(j,bi,bj).NE.OB_indexNone ) THEN
154 OBWu(j,k,bi,bj)= U_F1(j,k,1)
155 OBWv(j,k,bi,bj)= V_F1(j,k,1)
156 OBWt(j,k,bi,bj)= T_F1(j,k,1)
157 OBWs(j,k,bi,bj)= S_F1(j,k,1)
158 #ifdef NONLIN_FRSURF
159 OBWeta(j,bi,bj)= ETA_F1(j,1,1)
160 #endif
161 ENDIF
162 ENDDO
163 ENDDO
164 #endif /* ALLOW_NEST_CHILD */
165 ELSE
166 DO k=1,Nr
167 DO j=1-OLy,sNy+OLy
168 IF ( OB_Iw(j,bi,bj).NE.OB_indexNone ) THEN
169 OBWu(j,k,bi,bj)=0.
170 OBWv(j,k,bi,bj)=0.
171 OBWt(j,k,bi,bj)=tRef(k)
172 OBWs(j,k,bi,bj)=sRef(k)
173 #ifdef ALLOW_NONHYDROSTATIC
174 OBWw(j,k,bi,bj)=0.
175 #endif
176 #ifdef NONLIN_FRSURF
177 OBWeta(j,bi,bj)=0.
178 #endif
179 ENDIF
180 ENDDO
181 ENDDO
182 ENDIF
183 #endif /* ALLOW_OBCS_WEST */
184
185 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
186
187 #ifdef ALLOW_OBCS_NORTH
188 C Northern OB
189 #ifdef ALLOW_DEBUG
190 IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: North',myThid)
191 #endif
192 IF (useOrlanskiNorth) THEN
193 #ifdef ALLOW_ORLANSKI
194 CALL ORLANSKI_NORTH(
195 & bi, bj, futureTime,
196 & uVel, vVel, wVel, theta, salt,
197 & myThid )
198 #endif
199 ELSE
200 DO k=1,Nr
201 DO i=1-OLx,sNx+OLx
202 IF ( OB_Jn(i,bi,bj).NE.OB_indexNone ) THEN
203 OBNv(i,k,bi,bj)=0.
204 OBNu(i,k,bi,bj)=0.
205 OBNt(i,k,bi,bj)=tRef(k)
206 OBNs(i,k,bi,bj)=sRef(k)
207 #ifdef ALLOW_NONHYDROSTATIC
208 OBNw(i,k,bi,bj)=0.
209 #endif
210 #ifdef NONLIN_FRSURF
211 OBNeta(i,bi,bj)=0.
212 #endif
213 ENDIF
214 ENDDO
215 ENDDO
216 ENDIF
217 #endif /* ALLOW_OBCS_NORTH */
218
219 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
220
221 #ifdef ALLOW_OBCS_SOUTH
222 C Southern OB
223 #ifdef ALLOW_DEBUG
224 IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: South',myThid)
225 #endif
226 IF (useOrlanskiSouth) THEN
227 #ifdef ALLOW_ORLANSKI
228 CALL ORLANSKI_SOUTH(
229 & bi, bj, futureTime,
230 & uVel, vVel, wVel, theta, salt,
231 & myThid )
232 #endif
233 ELSE
234 DO k=1,Nr
235 DO i=1-OLx,sNx+OLx
236 IF ( OB_Js(i,bi,bj).NE.OB_indexNone ) THEN
237 OBSu(i,k,bi,bj)=0.
238 OBSv(i,k,bi,bj)=0.
239 OBSt(i,k,bi,bj)=tRef(k)
240 OBSs(i,k,bi,bj)=sRef(k)
241 #ifdef ALLOW_NONHYDROSTATIC
242 OBSw(i,k,bi,bj)=0.
243 #endif
244 #ifdef NONLIN_FRSURF
245 OBSeta(i,bi,bj)=0.
246 #endif
247 ENDIF
248 ENDDO
249 ENDDO
250 ENDIF
251 #endif /* ALLOW_OBCS_SOUTH */
252
253 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
254
255 #ifdef ALLOW_PTRACERS
256 IF ( usePTRACERS ) THEN
257 C
258 C Calculate some default open boundary conditions for passive tracers:
259 C The default is a homogeneous v.Neumann conditions, that is, the
260 C tracer gradient across the open boundary is nearly zero;
261 C only nearly, because the boundary conditions are applied throughout
262 C the time step during which the interior field does change; therefore
263 C we have to use the values from the previous time step here. If you
264 C really want exact v.Neumann conditions, you have to modify
265 C obcs_apply_ptracer directly.
266 C
267 # ifdef ALLOW_OBCS_EAST
268 C Eastern OB
269 # ifdef ALLOW_DEBUG
270 IF (debugMode)
271 & CALL DEBUG_MSG('OBCS_CALC: East, pTracers',myThid)
272 # endif
273 IF (useOrlanskiEast) THEN
274 WRITE(msgBuf,'(A)')
275 & 'OBCS_CALC: ERROR: useOrlanskiEast Rad OBC with'
276 CALL PRINT_ERROR( msgBuf, myThid )
277 WRITE(msgBuf,'(A)')
278 & 'OBCS_CALC: ERROR: pTracers not yet implemented'
279 CALL PRINT_ERROR( msgBuf, myThid )
280 STOP 'ABNORMAL END: S/R OBCS_CALC'
281 ELSE
282 DO iTracer=1,PTRACERS_numInUse
283 DO k=1,Nr
284 DO j=1-OLy,sNy+OLy
285 IF ( OB_Ie(j,bi,bj).NE.OB_indexNone ) THEN
286 I_obc = OB_Ie(j,bi,bj)
287 OBEptr(j,k,bi,bj,iTracer) =
288 & pTracer(I_obc-1,j,k,bi,bj,iTracer)
289 & *_maskW(I_obc,j,k,bi,bj)
290 ENDIF
291 ENDDO
292 ENDDO
293 ENDDO
294 ENDIF
295 # endif /* ALLOW_OBCS_EAST */
296
297 C ------------------------------------------------------------------------------
298
299 # ifdef ALLOW_OBCS_WEST
300 C Western OB
301 # ifdef ALLOW_DEBUG
302 IF (debugMode)
303 & CALL DEBUG_MSG('OBCS_CALC: West, pTracers',myThid)
304 # endif
305 IF (useOrlanskiWest) THEN
306 WRITE(msgBuf,'(A)')
307 & 'OBCS_CALC: ERROR: useOrlanskiWest Rad OBC with'
308 CALL PRINT_ERROR( msgBuf, myThid )
309 WRITE(msgBuf,'(A)')
310 & 'OBCS_CALC: ERROR: pTracers not yet implemented'
311 CALL PRINT_ERROR( msgBuf, myThid )
312 STOP 'ABNORMAL END: S/R OBCS_CALC'
313 ELSE
314 DO iTracer=1,PTRACERS_numInUse
315 DO k=1,Nr
316 DO j=1-OLy,sNy+OLy
317 IF ( OB_Iw(j,bi,bj).NE.OB_indexNone ) THEN
318 I_obc = OB_Iw(j,bi,bj)
319 OBWptr(j,k,bi,bj,iTracer) =
320 & pTracer(I_obc+1,j,k,bi,bj,iTracer)
321 & *_maskW(I_obc+1,j,k,bi,bj)
322 ENDIF
323 ENDDO
324 ENDDO
325 ENDDO
326 ENDIF
327 # endif /* ALLOW_OBCS_WEST */
328
329 C ------------------------------------------------------------------------------
330
331 # ifdef ALLOW_OBCS_NORTH
332 C Northern OB
333 # ifdef ALLOW_DEBUG
334 IF (debugMode)
335 & CALL DEBUG_MSG('OBCS_CALC: North, pTracers',myThid)
336 # endif
337 IF (useOrlanskiNorth) THEN
338 WRITE(msgBuf,'(A)')
339 & 'OBCS_CALC: ERROR: useOrlanskiNorth Rad OBC with'
340 CALL PRINT_ERROR( msgBuf, myThid )
341 WRITE(msgBuf,'(A)')
342 & 'OBCS_CALC: ERROR: pTracers not yet implemented'
343 CALL PRINT_ERROR( msgBuf, myThid )
344 STOP 'ABNORMAL END: S/R OBCS_CALC'
345 ELSE
346 DO iTracer=1,PTRACERS_numInUse
347 DO k=1,Nr
348 DO i=1-OLx,sNx+OLx
349 IF ( OB_Jn(i,bi,bj).NE.OB_indexNone ) THEN
350 J_obc = OB_Jn(i,bi,bj)
351 OBNptr(i,k,bi,bj,iTracer) =
352 & pTracer(i,J_obc-1,k,bi,bj,iTracer)
353 & *_maskS(i,J_obc,k,bi,bj)
354 ENDIF
355 ENDDO
356 ENDDO
357 ENDDO
358 ENDIF
359 # endif /* ALLOW_OBCS_NORTH */
360
361 C ------------------------------------------------------------------------------
362
363 # ifdef ALLOW_OBCS_SOUTH
364 C Southern OB
365 # ifdef ALLOW_DEBUG
366 IF (debugMode)
367 & CALL DEBUG_MSG('OBCS_CALC: South, pTracers',myThid)
368 #endif
369 IF (useOrlanskiSouth) THEN
370 WRITE(msgBuf,'(A)')
371 & 'OBCS_CALC: ERROR: useOrlanskiSouth Rad OBC with'
372 CALL PRINT_ERROR( msgBuf, myThid )
373 WRITE(msgBuf,'(A)')
374 & 'OBCS_CALC: ERROR: pTracers not yet implemented'
375 CALL PRINT_ERROR( msgBuf, myThid )
376 STOP 'ABNORMAL END: S/R OBCS_CALC'
377 ELSE
378 DO iTracer=1,PTRACERS_numInUse
379 DO k=1,Nr
380 DO i=1-OLx,sNx+OLx
381 IF ( OB_Js(i,bi,bj).NE.OB_indexNone ) THEN
382 J_obc = OB_Js(i,bi,bj)
383 OBSptr(i,k,bi,bj,iTracer) =
384 & pTracer(i,J_obc+1,k,bi,bj,iTracer)
385 & *_maskS(i,J_obc+1,k,bi,bj)
386 ENDIF
387 ENDDO
388 ENDDO
389 ENDDO
390 ENDIF
391 # endif /* ALLOW_OBCS_SOUTH */
392 C end if (usePTracers)
393 ENDIF
394 #endif /* ALLOW_PTRACERS */
395
396 C-- end bi,bj loops.
397 ENDDO
398 ENDDO
399
400 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
401
402 #ifdef ALLOW_OBCS_PRESCRIBE
403 IF (useOBCSprescribe) THEN
404 C-- Calculate future values on open boundaries
405 #ifdef ALLOW_DEBUG
406 IF (debugMode) CALL DEBUG_CALL('OBCS_PRESCRIBE_READ',myThid)
407 #endif
408 CALL OBCS_PRESCRIBE_READ( futureTime, futureIter, myThid )
409 ENDIF
410 #endif /* ALLOW_OBCS_PRESCRIBE */
411
412 C ------------------------------------------------------------------------------
413
414 #ifdef ALLOW_OBCS_STEVENS
415 C The Stevens (1990) boundary conditions come after reading data from
416 C files, because they use the data to compute a mix of simplified
417 C Orlanski and prescribed boundary conditions
418 IF (useStevensNorth.OR.useStevensSouth.OR.
419 & useStevensEast.OR.useStevensWest) THEN
420 #ifdef ALLOW_DEBUG
421 IF (debugMode) CALL DEBUG_CALL('OBCS_CALC_STEVENS',myThid)
422 #endif
423 CALL OBCS_CALC_STEVENS( futureTime, futureIter, myThid )
424 ENDIF
425 #endif /* ALLOW_OBCS_STEVENS */
426
427 #ifdef ALLOW_OBCS_BALANCE
428 IF ( useOBCSbalance ) THEN
429 CALL OBCS_BALANCE_FLOW( futureTime, futureIter, myThid )
430 ENDIF
431 #endif /* ALLOW_OBCS_BALANCE */
432
433 #ifdef ALLOW_OBCS_TIDES
434 IF ( useOBCStides ) THEN
435 CALL OBCS_ADD_TIDES( futureTime, futureIter, myThid )
436 ENDIF
437 #endif /* ALLOW_OBCS_TIDES */
438
439 #ifdef ALLOW_DEBUG
440 IF (debugMode) CALL DEBUG_LEAVE('OBCS_CALC',myThid)
441 #endif
442 #endif /* ALLOW_OBCS */
443
444 RETURN
445 END

  ViewVC Help
Powered by ViewVC 1.1.22