1 |
dimitri |
1.34 |
C $Header: /u/gcmpack/MITgcm/pkg/obcs/obcs_calc.F,v 1.33 2012/09/18 20:09:17 jmc Exp $ |
2 |
adcroft |
1.3 |
C $Name: $ |
3 |
adcroft |
1.2 |
|
4 |
|
|
#include "OBCS_OPTIONS.h" |
5 |
|
|
|
6 |
jmc |
1.30 |
CBOP |
7 |
|
|
C !ROUTINE: OBCS_CALC |
8 |
|
|
|
9 |
|
|
C !INTERFACE: |
10 |
jahn |
1.24 |
SUBROUTINE OBCS_CALC( futureTime, futureIter, |
11 |
jmc |
1.21 |
& uVel, vVel, wVel, theta, salt, |
12 |
adcroft |
1.2 |
& myThid ) |
13 |
jmc |
1.30 |
|
14 |
|
|
C !DESCRIPTION: |
15 |
jmc |
1.26 |
C *==========================================================* |
16 |
jmc |
1.30 |
C | SUBROUTINE OBCS_CALC |
17 |
|
|
C | o Calculate future boundary data at open boundaries |
18 |
|
|
C | at time = futureTime |
19 |
jmc |
1.26 |
C *==========================================================* |
20 |
jmc |
1.30 |
|
21 |
|
|
C !USES: |
22 |
adcroft |
1.2 |
IMPLICIT NONE |
23 |
|
|
|
24 |
|
|
C === Global variables === |
25 |
|
|
#include "SIZE.h" |
26 |
|
|
#include "EEPARAMS.h" |
27 |
|
|
#include "PARAMS.h" |
28 |
jmc |
1.21 |
#include "GRID.h" |
29 |
jmc |
1.32 |
#include "OBCS_PARAMS.h" |
30 |
|
|
#include "OBCS_GRID.h" |
31 |
|
|
#include "OBCS_FIELDS.h" |
32 |
jmc |
1.16 |
#ifdef ALLOW_PTRACERS |
33 |
jmc |
1.29 |
# include "PTRACERS_SIZE.h" |
34 |
|
|
# include "PTRACERS_PARAMS.h" |
35 |
|
|
# include "PTRACERS_FIELDS.h" |
36 |
|
|
# include "OBCS_PTRACERS.h" |
37 |
mlosch |
1.12 |
#endif /* ALLOW_PTRACERS */ |
38 |
jmc |
1.29 |
#ifdef ALLOW_NEST_CHILD |
39 |
|
|
# include "NEST_CHILD.h" |
40 |
|
|
#endif /* ALLOW_NEST_CHILD */ |
41 |
adcroft |
1.2 |
|
42 |
jmc |
1.30 |
C !INPUT/OUTPUT PARAMETERS: |
43 |
adcroft |
1.2 |
C == Routine arguments == |
44 |
jmc |
1.6 |
INTEGER futureIter |
45 |
adcroft |
1.2 |
_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 |
jmc |
1.30 |
C !LOCAL VARIABLES: |
56 |
jmc |
1.25 |
C bi, bj :: tile indices |
57 |
jmc |
1.33 |
C i,j,k :: loop indices |
58 |
jmc |
1.25 |
C I_obc, J_obc :: local index of open boundary |
59 |
|
|
C msgBuf :: Informational/error message buffer |
60 |
jahn |
1.24 |
INTEGER bi, bj |
61 |
jmc |
1.33 |
INTEGER i, j, k |
62 |
mlosch |
1.12 |
#ifdef ALLOW_PTRACERS |
63 |
jmc |
1.33 |
INTEGER I_obc, J_obc |
64 |
jmc |
1.25 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
65 |
mlosch |
1.12 |
INTEGER iTracer |
66 |
|
|
#endif /* ALLOW_PTRACERS */ |
67 |
|
|
|
68 |
jmc |
1.29 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
69 |
adcroft |
1.2 |
|
70 |
adcroft |
1.11 |
#ifdef ALLOW_DEBUG |
71 |
jmc |
1.28 |
IF (debugMode) CALL DEBUG_ENTER('OBCS_CALC',myThid) |
72 |
adcroft |
1.11 |
#endif |
73 |
|
|
|
74 |
jahn |
1.24 |
DO bj=myByLo(myThid),myByHi(myThid) |
75 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
76 |
|
|
|
77 |
jmc |
1.29 |
#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 |
heimbach |
1.7 |
#ifdef ALLOW_OBCS_EAST |
88 |
adcroft |
1.2 |
C Eastern OB |
89 |
adcroft |
1.11 |
#ifdef ALLOW_DEBUG |
90 |
jmc |
1.28 |
IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: East',myThid) |
91 |
adcroft |
1.11 |
#endif |
92 |
adcroft |
1.2 |
IF (useOrlanskiEast) THEN |
93 |
heimbach |
1.7 |
#ifdef ALLOW_ORLANSKI |
94 |
adcroft |
1.2 |
CALL ORLANSKI_EAST( |
95 |
jmc |
1.21 |
& bi, bj, futureTime, |
96 |
|
|
& uVel, vVel, wVel, theta, salt, |
97 |
adcroft |
1.2 |
& myThid ) |
98 |
heimbach |
1.7 |
#endif |
99 |
jmc |
1.29 |
#ifdef ALLOW_NEST_CHILD |
100 |
|
|
ELSEIF ( useNEST_CHILD ) THEN |
101 |
|
|
DO k=1,Nr |
102 |
jmc |
1.33 |
DO j=1-OLy,sNy+OLy |
103 |
|
|
IF ( OB_Ie(j,bi,bj).NE.OB_indexNone ) THEN |
104 |
jmc |
1.29 |
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 |
adcroft |
1.2 |
ELSE |
116 |
jmc |
1.33 |
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 |
adcroft |
1.2 |
#ifdef ALLOW_NONHYDROSTATIC |
124 |
jmc |
1.33 |
OBEw(j,k,bi,bj)=0. |
125 |
adcroft |
1.2 |
#endif |
126 |
jmc |
1.5 |
#ifdef NONLIN_FRSURF |
127 |
jmc |
1.33 |
OBEeta(j,bi,bj)=0. |
128 |
jmc |
1.5 |
#endif |
129 |
adcroft |
1.3 |
ENDIF |
130 |
adcroft |
1.2 |
ENDDO |
131 |
|
|
ENDDO |
132 |
|
|
ENDIF |
133 |
heimbach |
1.8 |
#endif /* ALLOW_OBCS_EAST */ |
134 |
heimbach |
1.7 |
|
135 |
jmc |
1.29 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
136 |
heimbach |
1.7 |
|
137 |
|
|
#ifdef ALLOW_OBCS_WEST |
138 |
adcroft |
1.2 |
C Western OB |
139 |
adcroft |
1.11 |
#ifdef ALLOW_DEBUG |
140 |
jmc |
1.28 |
IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: West',myThid) |
141 |
adcroft |
1.11 |
#endif |
142 |
adcroft |
1.2 |
IF (useOrlanskiWest) THEN |
143 |
heimbach |
1.7 |
#ifdef ALLOW_ORLANSKI |
144 |
adcroft |
1.2 |
CALL ORLANSKI_WEST( |
145 |
jmc |
1.21 |
& bi, bj, futureTime, |
146 |
|
|
& uVel, vVel, wVel, theta, salt, |
147 |
adcroft |
1.2 |
& myThid ) |
148 |
heimbach |
1.7 |
#endif |
149 |
jmc |
1.29 |
#ifdef ALLOW_NEST_CHILD |
150 |
|
|
ELSEIF ( useNEST_CHILD ) THEN |
151 |
|
|
DO k=1,Nr |
152 |
jmc |
1.33 |
DO j=1-OLy,sNy+OLy |
153 |
|
|
IF ( OB_Iw(j,bi,bj).NE.OB_indexNone ) THEN |
154 |
jmc |
1.29 |
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 |
adcroft |
1.2 |
ELSE |
166 |
jmc |
1.33 |
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 |
adcroft |
1.2 |
#ifdef ALLOW_NONHYDROSTATIC |
174 |
jmc |
1.33 |
OBWw(j,k,bi,bj)=0. |
175 |
jmc |
1.21 |
#endif |
176 |
jmc |
1.5 |
#ifdef NONLIN_FRSURF |
177 |
jmc |
1.33 |
OBWeta(j,bi,bj)=0. |
178 |
jmc |
1.5 |
#endif |
179 |
heimbach |
1.7 |
ENDIF |
180 |
|
|
ENDDO |
181 |
|
|
ENDDO |
182 |
|
|
ENDIF |
183 |
|
|
#endif /* ALLOW_OBCS_WEST */ |
184 |
|
|
|
185 |
jmc |
1.29 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
186 |
heimbach |
1.7 |
|
187 |
|
|
#ifdef ALLOW_OBCS_NORTH |
188 |
adcroft |
1.3 |
C Northern OB |
189 |
adcroft |
1.11 |
#ifdef ALLOW_DEBUG |
190 |
jmc |
1.28 |
IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: North',myThid) |
191 |
adcroft |
1.11 |
#endif |
192 |
adcroft |
1.2 |
IF (useOrlanskiNorth) THEN |
193 |
heimbach |
1.7 |
#ifdef ALLOW_ORLANSKI |
194 |
adcroft |
1.2 |
CALL ORLANSKI_NORTH( |
195 |
jmc |
1.21 |
& bi, bj, futureTime, |
196 |
|
|
& uVel, vVel, wVel, theta, salt, |
197 |
adcroft |
1.2 |
& myThid ) |
198 |
heimbach |
1.7 |
#endif |
199 |
adcroft |
1.2 |
ELSE |
200 |
jmc |
1.33 |
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 |
adcroft |
1.2 |
#ifdef ALLOW_NONHYDROSTATIC |
208 |
jmc |
1.33 |
OBNw(i,k,bi,bj)=0. |
209 |
adcroft |
1.2 |
#endif |
210 |
jmc |
1.5 |
#ifdef NONLIN_FRSURF |
211 |
jmc |
1.33 |
OBNeta(i,bi,bj)=0. |
212 |
jmc |
1.5 |
#endif |
213 |
adcroft |
1.3 |
ENDIF |
214 |
adcroft |
1.2 |
ENDDO |
215 |
|
|
ENDDO |
216 |
|
|
ENDIF |
217 |
heimbach |
1.8 |
#endif /* ALLOW_OBCS_NORTH */ |
218 |
heimbach |
1.7 |
|
219 |
jmc |
1.29 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
220 |
heimbach |
1.7 |
|
221 |
|
|
#ifdef ALLOW_OBCS_SOUTH |
222 |
adcroft |
1.3 |
C Southern OB |
223 |
adcroft |
1.11 |
#ifdef ALLOW_DEBUG |
224 |
jmc |
1.28 |
IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: South',myThid) |
225 |
adcroft |
1.11 |
#endif |
226 |
jmc |
1.21 |
IF (useOrlanskiSouth) THEN |
227 |
heimbach |
1.7 |
#ifdef ALLOW_ORLANSKI |
228 |
adcroft |
1.2 |
CALL ORLANSKI_SOUTH( |
229 |
jmc |
1.21 |
& bi, bj, futureTime, |
230 |
|
|
& uVel, vVel, wVel, theta, salt, |
231 |
adcroft |
1.2 |
& myThid ) |
232 |
heimbach |
1.7 |
#endif |
233 |
adcroft |
1.2 |
ELSE |
234 |
jmc |
1.33 |
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 |
adcroft |
1.2 |
#ifdef ALLOW_NONHYDROSTATIC |
242 |
jmc |
1.33 |
OBSw(i,k,bi,bj)=0. |
243 |
jmc |
1.5 |
#endif |
244 |
|
|
#ifdef NONLIN_FRSURF |
245 |
jmc |
1.33 |
OBSeta(i,bi,bj)=0. |
246 |
adcroft |
1.2 |
#endif |
247 |
adcroft |
1.3 |
ENDIF |
248 |
adcroft |
1.2 |
ENDDO |
249 |
|
|
ENDDO |
250 |
|
|
ENDIF |
251 |
heimbach |
1.8 |
#endif /* ALLOW_OBCS_SOUTH */ |
252 |
|
|
|
253 |
jmc |
1.29 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
254 |
|
|
|
255 |
mlosch |
1.12 |
#ifdef ALLOW_PTRACERS |
256 |
|
|
IF ( usePTRACERS ) THEN |
257 |
|
|
C |
258 |
|
|
C Calculate some default open boundary conditions for passive tracers: |
259 |
jmc |
1.21 |
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 |
mlosch |
1.12 |
C the time step during which the interior field does change; therefore |
263 |
jmc |
1.21 |
C we have to use the values from the previous time step here. If you |
264 |
mlosch |
1.12 |
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 |
jmc |
1.28 |
IF (debugMode) |
271 |
mlosch |
1.12 |
& 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 |
jmc |
1.22 |
CALL PRINT_ERROR( msgBuf, myThid ) |
277 |
mlosch |
1.12 |
WRITE(msgBuf,'(A)') |
278 |
|
|
& 'OBCS_CALC: ERROR: pTracers not yet implemented' |
279 |
jmc |
1.22 |
CALL PRINT_ERROR( msgBuf, myThid ) |
280 |
mlosch |
1.12 |
STOP 'ABNORMAL END: S/R OBCS_CALC' |
281 |
|
|
ELSE |
282 |
|
|
DO iTracer=1,PTRACERS_numInUse |
283 |
jmc |
1.33 |
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 |
mlosch |
1.12 |
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 |
jmc |
1.28 |
IF (debugMode) |
303 |
mlosch |
1.12 |
& 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 |
jmc |
1.22 |
CALL PRINT_ERROR( msgBuf, myThid ) |
309 |
mlosch |
1.12 |
WRITE(msgBuf,'(A)') |
310 |
|
|
& 'OBCS_CALC: ERROR: pTracers not yet implemented' |
311 |
jmc |
1.22 |
CALL PRINT_ERROR( msgBuf, myThid ) |
312 |
mlosch |
1.12 |
STOP 'ABNORMAL END: S/R OBCS_CALC' |
313 |
|
|
ELSE |
314 |
|
|
DO iTracer=1,PTRACERS_numInUse |
315 |
jmc |
1.33 |
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 |
mlosch |
1.12 |
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 |
jmc |
1.28 |
IF (debugMode) |
335 |
mlosch |
1.12 |
& 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 |
jmc |
1.22 |
CALL PRINT_ERROR( msgBuf, myThid ) |
341 |
mlosch |
1.12 |
WRITE(msgBuf,'(A)') |
342 |
|
|
& 'OBCS_CALC: ERROR: pTracers not yet implemented' |
343 |
jmc |
1.22 |
CALL PRINT_ERROR( msgBuf, myThid ) |
344 |
mlosch |
1.12 |
STOP 'ABNORMAL END: S/R OBCS_CALC' |
345 |
|
|
ELSE |
346 |
|
|
DO iTracer=1,PTRACERS_numInUse |
347 |
jmc |
1.33 |
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 |
mlosch |
1.12 |
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 |
jmc |
1.28 |
IF (debugMode) |
367 |
mlosch |
1.12 |
& CALL DEBUG_MSG('OBCS_CALC: South, pTracers',myThid) |
368 |
|
|
#endif |
369 |
jmc |
1.21 |
IF (useOrlanskiSouth) THEN |
370 |
mlosch |
1.12 |
WRITE(msgBuf,'(A)') |
371 |
|
|
& 'OBCS_CALC: ERROR: useOrlanskiSouth Rad OBC with' |
372 |
jmc |
1.22 |
CALL PRINT_ERROR( msgBuf, myThid ) |
373 |
mlosch |
1.12 |
WRITE(msgBuf,'(A)') |
374 |
|
|
& 'OBCS_CALC: ERROR: pTracers not yet implemented' |
375 |
jmc |
1.22 |
CALL PRINT_ERROR( msgBuf, myThid ) |
376 |
mlosch |
1.12 |
STOP 'ABNORMAL END: S/R OBCS_CALC' |
377 |
|
|
ELSE |
378 |
|
|
DO iTracer=1,PTRACERS_numInUse |
379 |
jmc |
1.33 |
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 |
mlosch |
1.12 |
ENDIF |
387 |
|
|
ENDDO |
388 |
|
|
ENDDO |
389 |
|
|
ENDDO |
390 |
|
|
ENDIF |
391 |
|
|
# endif /* ALLOW_OBCS_SOUTH */ |
392 |
|
|
C end if (usePTracers) |
393 |
jmc |
1.21 |
ENDIF |
394 |
mlosch |
1.12 |
#endif /* ALLOW_PTRACERS */ |
395 |
heimbach |
1.8 |
|
396 |
jahn |
1.24 |
C-- end bi,bj loops. |
397 |
jmc |
1.25 |
ENDDO |
398 |
jahn |
1.24 |
ENDDO |
399 |
|
|
|
400 |
jmc |
1.29 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
401 |
heimbach |
1.8 |
|
402 |
heimbach |
1.9 |
#ifdef ALLOW_OBCS_PRESCRIBE |
403 |
|
|
IF (useOBCSprescribe) THEN |
404 |
|
|
C-- Calculate future values on open boundaries |
405 |
adcroft |
1.11 |
#ifdef ALLOW_DEBUG |
406 |
jmc |
1.28 |
IF (debugMode) CALL DEBUG_CALL('OBCS_PRESCRIBE_READ',myThid) |
407 |
adcroft |
1.11 |
#endif |
408 |
jmc |
1.22 |
CALL OBCS_PRESCRIBE_READ( futureTime, futureIter, myThid ) |
409 |
heimbach |
1.9 |
ENDIF |
410 |
heimbach |
1.8 |
#endif /* ALLOW_OBCS_PRESCRIBE */ |
411 |
|
|
|
412 |
|
|
C ------------------------------------------------------------------------------ |
413 |
jmc |
1.29 |
|
414 |
mlosch |
1.27 |
#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 |
jmc |
1.28 |
IF (debugMode) CALL DEBUG_CALL('OBCS_CALC_STEVENS',myThid) |
422 |
mlosch |
1.27 |
#endif |
423 |
|
|
CALL OBCS_CALC_STEVENS( futureTime, futureIter, myThid ) |
424 |
|
|
ENDIF |
425 |
|
|
#endif /* ALLOW_OBCS_STEVENS */ |
426 |
jmc |
1.29 |
|
427 |
heimbach |
1.7 |
#ifdef ALLOW_OBCS_BALANCE |
428 |
jmc |
1.30 |
IF ( useOBCSbalance ) THEN |
429 |
|
|
CALL OBCS_BALANCE_FLOW( futureTime, futureIter, myThid ) |
430 |
heimbach |
1.7 |
ENDIF |
431 |
|
|
#endif /* ALLOW_OBCS_BALANCE */ |
432 |
|
|
|
433 |
dimitri |
1.34 |
#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 |
adcroft |
1.11 |
#ifdef ALLOW_DEBUG |
440 |
jmc |
1.28 |
IF (debugMode) CALL DEBUG_LEAVE('OBCS_CALC',myThid) |
441 |
adcroft |
1.11 |
#endif |
442 |
jmc |
1.29 |
#endif /* ALLOW_OBCS */ |
443 |
|
|
|
444 |
adcroft |
1.2 |
RETURN |
445 |
|
|
END |