Back to home page

MITgcm

 
 

    


File indexing completed on 2019-08-24 05:11:32 UTC

view on githubraw file Latest commit abfe198b on 2019-08-23 19:59:52 UTC
00816bc2b8 Alis*0001 #include "OBCS_OPTIONS.h"
                0002 
abfe198bce Mart*0003 CBOP
                0004 C     !ROUTINE: OBCS_CALC
                0005 
                0006 C     !INTERFACE:
3a86c9b47d Oliv*0007       SUBROUTINE OBCS_CALC( futureTime, futureIter,
562221c9a5 Jean*0008      &                      uVel, vVel, wVel, theta, salt,
00816bc2b8 Alis*0009      &                      myThid )
abfe198bce Mart*0010 
                0011 C     !DESCRIPTION:
2e53480479 Jean*0012 C     *==========================================================*
2b169ecc4b Jean*0013 C     | SUBROUTINE OBCS_CALC
                0014 C     | o Calculate future boundary data at open boundaries
                0015 C     |   at time = futureTime
2e53480479 Jean*0016 C     *==========================================================*
abfe198bce Mart*0017 
                0018 C     !USES:
00816bc2b8 Alis*0019       IMPLICIT NONE
                0020 
                0021 C     === Global variables ===
                0022 #include "SIZE.h"
                0023 #include "EEPARAMS.h"
                0024 #include "PARAMS.h"
2e53480479 Jean*0025 #include "GRID.h"
65fe50e2dd Jean*0026 #include "OBCS_PARAMS.h"
                0027 #include "OBCS_GRID.h"
                0028 #include "OBCS_FIELDS.h"
00816bc2b8 Alis*0029 #include "EOS.h"
abfe198bce Mart*0030 #ifdef ALLOW_PTRACERS
                0031 # include "PTRACERS_SIZE.h"
                0032 # include "PTRACERS_PARAMS.h"
                0033 # include "PTRACERS_FIELDS.h"
                0034 # include "OBCS_PTRACERS.h"
                0035 #endif /* ALLOW_PTRACERS */
                0036 #ifdef ALLOW_NEST_CHILD
                0037 # include "NEST_CHILD.h"
                0038 #endif /* ALLOW_NEST_CHILD */
00816bc2b8 Alis*0039 
abfe198bce Mart*0040 C     !INPUT/OUTPUT PARAMETERS:
00816bc2b8 Alis*0041 C     == Routine arguments ==
                0042       INTEGER futureIter
                0043       _RL futureTime
                0044       _RL uVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
                0045       _RL vVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
                0046       _RL wVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
                0047       _RL theta(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
                0048       _RL salt (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
                0049       INTEGER myThid
                0050 
                0051 #ifdef ALLOW_OBCS
                0052 
abfe198bce Mart*0053 C     !LOCAL VARIABLES:
                0054 C     bi, bj       :: tile indices
                0055 C     i,j,k        :: loop indices
                0056 C     I_obc, J_obc :: local index of open boundary
                0057 C     msgBuf       :: Informational/error message buffer
3a86c9b47d Oliv*0058       INTEGER bi, bj
abfe198bce Mart*0059       INTEGER i, j, k
2e53480479 Jean*0060       _RL rampTime2
00816bc2b8 Alis*0061       _RL rexpt
                0062       _RL hinit, delh
                0063       _RL z(nr)
                0064       _RL Dmax,Dinf,dTemp,gp_inflow,Lrho
2e53480479 Jean*0065       _RL Width,x,Xcenter,Fz,zt,Rit
abfe198bce Mart*0066 #ifdef ALLOW_PTRACERS
                0067       INTEGER I_obc, J_obc
                0068       CHARACTER*(MAX_LEN_MBUF) msgBuf
                0069       INTEGER iTracer
                0070 #endif /* ALLOW_PTRACERS */
00816bc2b8 Alis*0071 
2b169ecc4b Jean*0072 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0073 
                0074 #ifdef ALLOW_DEBUG
                0075       IF (debugMode) CALL DEBUG_ENTER('OBCS_CALC',myThid)
                0076 #endif
                0077 
00816bc2b8 Alis*0078 c     Total depth
2e22b790fc Jean*0079       Dmax = 600. _d 0
00816bc2b8 Alis*0080 c     width of flow transition
2e22b790fc Jean*0081       delh=5. _d 0
00816bc2b8 Alis*0082 c     total depth of active layer
2e22b790fc Jean*0083       Dinf=300. _d 0
88f72205aa Jean*0084 c     g_prime of inflow
2e22b790fc Jean*0085       dTemp=(2. _d 0/rhonil)/(2. _d -4)
2b169ecc4b Jean*0086       gp_inflow=tAlpha*gravity*dTemp
00816bc2b8 Alis*0087 c     Deformation radius
                0088       Lrho=sqrt(gp_inflow*Dinf)/f0
562221c9a5 Jean*0089 c     dimensional width
2e22b790fc Jean*0090       Width=100. _d 3
00816bc2b8 Alis*0091 c     nondimensional width
                0092       Width=Width/Lrho
                0093 c      print *,'Width=',Width
                0094 c     coordinate of center of embayment
2e22b790fc Jean*0095       Xcenter = 1700. _d 3
562221c9a5 Jean*0096 c     Critical Richardson number
2e22b790fc Jean*0097       Rit = 1. _d 0/3. _d 0
00816bc2b8 Alis*0098 c      Ttarget=-5.0e6
                0099 c      areaCh=Width*Lrho*Dmax
                0100 c      relaxtime=180000.0
                0101 c      relaxConst=1.0/(areaCh*relaxtime)
                0102 
                0103 
2e22b790fc Jean*0104       z(1) = -delR(1)/2. _d 0
00816bc2b8 Alis*0105       do j=2,nr
                0106 Caja NOTE: This gives the depths of W points, not U,V,T,S,p points
                0107 C          USe instead drC but requires different indexing. See documentation
                0108 C          in AJAs head
                0109       z(j) = z(j-1) - delR(j)
                0110       enddo
                0111 
                0112 c     ramptime for velocity
2e22b790fc Jean*0113       rampTime2 = 4. _d 0*44567. _d 0
                0114       rexpt=1. _d 0/exp(futureTime/rampTime2)
00816bc2b8 Alis*0115 
3a86c9b47d Oliv*0116       DO bj=myByLo(myThid),myByHi(myThid)
                0117       DO bi=myBxLo(myThid),myBxHi(myThid)
                0118 
abfe198bce Mart*0119 #ifdef ALLOW_NEST_CHILD
                0120       IF ( useNEST_CHILD ) THEN
                0121         IF ( PASSI.LT.2 ) THEN
                0122          CALL NEST_CHILD_RECV ( myThid )
                0123         ENDIF
                0124       ENDIF
                0125 #endif /* ALLOW_NEST_CHILD */
                0126 
                0127 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0128 
                0129 #ifdef ALLOW_OBCS_EAST
00816bc2b8 Alis*0130 C     Eastern OB
abfe198bce Mart*0131 #ifdef ALLOW_DEBUG
                0132       IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: East',myThid)
                0133 #endif
00816bc2b8 Alis*0134       IF (useOrlanskiEast) THEN
abfe198bce Mart*0135 #ifdef ALLOW_ORLANSKI
00816bc2b8 Alis*0136         CALL ORLANSKI_EAST(
562221c9a5 Jean*0137      &          bi, bj, futureTime,
                0138      &          uVel, vVel, wVel, theta, salt,
00816bc2b8 Alis*0139      &          myThid )
abfe198bce Mart*0140 #endif
                0141 #ifdef ALLOW_NEST_CHILD
                0142       ELSEIF ( useNEST_CHILD ) THEN
                0143         DO k=1,Nr
                0144           DO j=1-OLy,sNy+OLy
                0145             IF ( OB_Ie(j,bi,bj).NE.OB_indexNone ) THEN
                0146               OBEu(j,k,bi,bj)= U_F1(j,k,2)
                0147               OBEv(j,k,bi,bj)= V_F1(j,k,2)
                0148               OBEt(j,k,bi,bj)= T_F1(j,k,2)
                0149               OBEs(j,k,bi,bj)= S_F1(j,k,2)
                0150 #ifdef NONLIN_FRSURF
                0151               OBEeta(j,bi,bj)= ETA_F1(j,1,2)
                0152 #endif
                0153             ENDIF
                0154           ENDDO
                0155         ENDDO
                0156 #endif /* ALLOW_NEST_CHILD */
00816bc2b8 Alis*0157       ELSE
abfe198bce Mart*0158         DO k=1,Nr
                0159           DO j=1-OLy,sNy+OLy
                0160             IF ( OB_Ie(j,bi,bj).NE.OB_indexNone ) THEN
                0161              OBEu(j,k,bi,bj)=0.
00816bc2b8 Alis*0162 c     &       +Uinflow
                0163 c     & *0.5*(1.0 - (exp((z(k)-hinit)/delh)
                0164 c     &   - exp(-(z(k)-hinit)/delh))
                0165 c     &   /(exp((z(k)-hinit)/delh)
                0166 c     &  + exp(-(z(k)-hinit)/delh)))
                0167 c     &*(1.0 - rexpt*rexpt)/(1.0 + rexpt*rexpt)
abfe198bce Mart*0168              OBEv(j,k,bi,bj)=0.
                0169              OBEt(j,k,bi,bj)=tRef(k)
00816bc2b8 Alis*0170 c     &  +Tinflow
                0171 c     & *0.5*(1.0 - (exp((z(k)-hinit)/delh)
                0172 c     &   - exp(-(z(k)-hinit)/delh))
                0173 c     &   /(exp((z(k)-hinit)/delh)
                0174 c     &  + exp(-(z(k)-hinit)/delh)))
562221c9a5 Jean*0175 c     &*(1.0 - rexpt*rexpt)/(1.0 + rexpt*rexpt)
abfe198bce Mart*0176               OBEs(j,k,bi,bj)=sRef(k)
00816bc2b8 Alis*0177 #ifdef ALLOW_NONHYDROSTATIC
abfe198bce Mart*0178               OBEw(j,k,bi,bj)=0.
                0179 #endif
                0180 #ifdef NONLIN_FRSURF
                0181               OBEeta(j,bi,bj)=0.
00816bc2b8 Alis*0182 #endif
abfe198bce Mart*0183             ENDIF
00816bc2b8 Alis*0184           ENDDO
                0185         ENDDO
                0186       ENDIF
abfe198bce Mart*0187 #endif /* ALLOW_OBCS_EAST */
                0188 
                0189 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
00816bc2b8 Alis*0190 
abfe198bce Mart*0191 #ifdef ALLOW_OBCS_WEST
00816bc2b8 Alis*0192 C     Western OB
abfe198bce Mart*0193 #ifdef ALLOW_DEBUG
                0194       IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: West',myThid)
                0195 #endif
00816bc2b8 Alis*0196       IF (useOrlanskiWest) THEN
abfe198bce Mart*0197 #ifdef ALLOW_ORLANSKI
00816bc2b8 Alis*0198         CALL ORLANSKI_WEST(
562221c9a5 Jean*0199      &          bi, bj, futureTime,
                0200      &          uVel, vVel, wVel, theta, salt,
00816bc2b8 Alis*0201      &          myThid )
abfe198bce Mart*0202 #endif
                0203 #ifdef ALLOW_NEST_CHILD
                0204       ELSEIF ( useNEST_CHILD ) THEN
                0205         DO k=1,Nr
                0206           DO j=1-OLy,sNy+OLy
                0207             IF ( OB_Iw(j,bi,bj).NE.OB_indexNone ) THEN
                0208               OBWu(j,k,bi,bj)= U_F1(j,k,1)
                0209               OBWv(j,k,bi,bj)= V_F1(j,k,1)
                0210               OBWt(j,k,bi,bj)= T_F1(j,k,1)
                0211               OBWs(j,k,bi,bj)= S_F1(j,k,1)
                0212 #ifdef NONLIN_FRSURF
                0213               OBWeta(j,bi,bj)= ETA_F1(j,1,1)
                0214 #endif
                0215            ENDIF
                0216           ENDDO
                0217         ENDDO
                0218 #endif /* ALLOW_NEST_CHILD */
00816bc2b8 Alis*0219       ELSE
abfe198bce Mart*0220         DO k=1,Nr
                0221           DO j=1-OLy,sNy+OLy
                0222             IF ( OB_Iw(j,bi,bj).NE.OB_indexNone ) THEN
                0223               OBWu(j,k,bi,bj)=0.
                0224               OBWv(j,k,bi,bj)=0.
                0225               OBWt(j,k,bi,bj)=tRef(k)
                0226               OBWs(j,k,bi,bj)=sRef(k)
00816bc2b8 Alis*0227 #ifdef ALLOW_NONHYDROSTATIC
abfe198bce Mart*0228               OBWw(j,k,bi,bj)=0.
00816bc2b8 Alis*0229 #endif
abfe198bce Mart*0230 #ifdef NONLIN_FRSURF
                0231               OBWeta(j,bi,bj)=0.
                0232 #endif
                0233            ENDIF
00816bc2b8 Alis*0234           ENDDO
                0235         ENDDO
                0236       ENDIF
abfe198bce Mart*0237 #endif /* ALLOW_OBCS_WEST */
00816bc2b8 Alis*0238 
abfe198bce Mart*0239 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0240 
                0241 #ifdef ALLOW_OBCS_NORTH
00816bc2b8 Alis*0242 C         Northern OB, template for forcing
abfe198bce Mart*0243 #ifdef ALLOW_DEBUG
                0244       IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: North',myThid)
                0245 #endif
00816bc2b8 Alis*0246       IF (useOrlanskiNorth) THEN
abfe198bce Mart*0247 #ifdef ALLOW_ORLANSKI
00816bc2b8 Alis*0248         CALL ORLANSKI_NORTH(
562221c9a5 Jean*0249      &          bi, bj, futureTime,
                0250      &          uVel, vVel, wVel, theta, salt,
00816bc2b8 Alis*0251      &          myThid )
abfe198bce Mart*0252 #endif
00816bc2b8 Alis*0253       ELSE
abfe198bce Mart*0254         DO k=1,Nr
                0255           DO i=1-OLx,sNx+OLx
                0256             IF ( OB_Jn(i,bi,bj).NE.OB_indexNone ) THEN
00816bc2b8 Alis*0257 c Make center of embayment x=0
abfe198bce Mart*0258               x=(xC(I,1,bi,bj)-Xcenter)/Lrho +Width/2. _d 0
                0259               IF ((x.GE.0.).AND.(x.LE.Width)) THEN
                0260                hinit=Dinf*(exp(-x)) - Dmax
00816bc2b8 Alis*0261 c      if (k.eq.5) then
                0262 c      print *,'x=',x,'  xC',xC(I,1,bi,bj),'  I=',I
                0263 c      print *,'hinit=',hinit
                0264 c      endif
abfe198bce Mart*0265                zt=(z(k) + Dmax - (hinit + Dmax))/(hinit+Dmax)
                0266                IF (zt.GE.(Rit/(2. _d 0 -Rit))) THEN
                0267                 Fz=1.
                0268                ELSE
                0269                 IF (zt.GE.(-Rit/(2.+Rit))) THEN
                0270                  Fz=1. _d 0/Rit*zt/(zt+1. _d 0) + 0.5 _d 0
                0271                 ELSE
                0272                  Fz=0.
                0273                 ENDIF
                0274                ENDIF
                0275               ELSE
                0276                Fz=1.
                0277               ENDIF
00816bc2b8 Alis*0278 c       if ((x.ge.0).and.(x.le.Width)) then
                0279 c        print *,'z(k)',z(k),'zt',zt,'Fz',Fz
a0318ce493 Jean*0280 c       endif
abfe198bce Mart*0281               OBNv(i,k,bi,bj)=0.
00816bc2b8 Alis*0282      &     -  sqrt(gp_inflow*Dinf)*exp(-x)
562221c9a5 Jean*0283      &                            *(1. _d 0 - Fz)
00816bc2b8 Alis*0284 c           if ((x.ge.0).and.(x.le.Width)) then
                0285 c           if (k.eq.5) then
                0286 c       print *,'V=',OBNv(I,K,bi,bj)
                0287 c           endif
                0288 c           endif
abfe198bce Mart*0289               OBNu(i,k,bi,bj)=0.
                0290               IF (tRef(K).LE. (-dTemp*(1. _d 0 - Fz))) THEN
                0291                OBNt(I,K,bi,bj) = tRef(K)
                0292               ELSE
                0293                OBNt(i,k,bi,bj) = -dTemp*(1. _d 0 - Fz)
                0294               ENDIF
00816bc2b8 Alis*0295 c           if ((x.ge.0).and.(x.le.Width)) then
abfe198bce Mart*0296 c       print *,'T=',OBNt(i,k,bi,bj)
a0318ce493 Jean*0297 c           endif
abfe198bce Mart*0298 c            OBNt(i,k,bi,bj)=tRef(k)
00816bc2b8 Alis*0299 c     & - dTemp
                0300 c     & *(1.0 - Fz)
abfe198bce Mart*0301               OBNs(i,k,bi,bj)=sRef(k) + 1. _d 0*(1. _d 0 - Fz)
00816bc2b8 Alis*0302 #ifdef ALLOW_NONHYDROSTATIC
abfe198bce Mart*0303               OBNw(i,k,bi,bj)=0.
00816bc2b8 Alis*0304 #endif
abfe198bce Mart*0305 #ifdef NONLIN_FRSURF
                0306               OBNeta(i,bi,bj)=0.
                0307 #endif
                0308             ENDIF
00816bc2b8 Alis*0309           ENDDO
                0310         ENDDO
                0311       ENDIF
abfe198bce Mart*0312 #endif /* ALLOW_OBCS_NORTH */
                0313 
                0314 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
00816bc2b8 Alis*0315 
abfe198bce Mart*0316 #ifdef ALLOW_OBCS_SOUTH
                0317 C         Southern OBC, template for forcing
                0318 #ifdef ALLOW_DEBUG
                0319       IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: South',myThid)
                0320 #endif
562221c9a5 Jean*0321       IF (useOrlanskiSouth) THEN
abfe198bce Mart*0322 #ifdef ALLOW_ORLANSKI
00816bc2b8 Alis*0323         CALL ORLANSKI_SOUTH(
562221c9a5 Jean*0324      &          bi, bj, futureTime,
                0325      &          uVel, vVel, wVel, theta, salt,
00816bc2b8 Alis*0326      &          myThid )
abfe198bce Mart*0327 #endif
00816bc2b8 Alis*0328       ELSE
abfe198bce Mart*0329         DO k=1,Nr
                0330           DO i=1-OLx,sNx+OLx
                0331             IF ( OB_Js(i,bi,bj).NE.OB_indexNone ) THEN
                0332               OBSu(i,k,bi,bj)=0.
                0333               OBSv(i,k,bi,bj)=0.
                0334               OBSt(i,k,bi,bj)=tRef(k)
                0335               OBSs(i,k,bi,bj)=sRef(k)
00816bc2b8 Alis*0336 #ifdef ALLOW_NONHYDROSTATIC
abfe198bce Mart*0337               OBSw(i,k,bi,bj)=0.
00816bc2b8 Alis*0338 #endif
abfe198bce Mart*0339 #ifdef NONLIN_FRSURF
                0340               OBSeta(i,bi,bj)=0.
                0341 #endif
                0342             ENDIF
00816bc2b8 Alis*0343           ENDDO
                0344         ENDDO
                0345       ENDIF
abfe198bce Mart*0346 #endif /* ALLOW_OBCS_SOUTH */
                0347 
                0348 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0349 
                0350 #ifdef ALLOW_PTRACERS
                0351       IF ( usePTRACERS ) THEN
                0352 C
                0353 C     Calculate some default open boundary conditions for passive tracers:
                0354 C     The default is a homogeneous v.Neumann conditions, that is, the
                0355 C     tracer gradient across the open boundary is nearly zero;
                0356 C     only nearly, because the boundary conditions are applied throughout
                0357 C     the time step during which the interior field does change; therefore
                0358 C     we have to use the values from the previous time step here. If you
                0359 C     really want exact v.Neumann conditions, you have to modify
                0360 C     obcs_apply_ptracer directly.
                0361 C
                0362 # ifdef ALLOW_OBCS_EAST
                0363 C     Eastern OB
                0364 #  ifdef ALLOW_DEBUG
                0365        IF (debugMode)
                0366      &      CALL DEBUG_MSG('OBCS_CALC: East, pTracers',myThid)
                0367 #  endif
                0368        IF (useOrlanskiEast) THEN
                0369         WRITE(msgBuf,'(A)')
                0370      &       'OBCS_CALC: ERROR: useOrlanskiEast Rad OBC with'
                0371         CALL PRINT_ERROR( msgBuf, myThid )
                0372         WRITE(msgBuf,'(A)')
                0373      &       'OBCS_CALC: ERROR: pTracers not yet implemented'
                0374         CALL PRINT_ERROR( msgBuf, myThid )
                0375         STOP 'ABNORMAL END: S/R OBCS_CALC'
                0376        ELSE
                0377         DO iTracer=1,PTRACERS_numInUse
                0378          DO k=1,Nr
                0379           DO j=1-OLy,sNy+OLy
                0380            IF ( OB_Ie(j,bi,bj).NE.OB_indexNone ) THEN
                0381             I_obc = OB_Ie(j,bi,bj)
                0382             OBEptr(j,k,bi,bj,iTracer) =
                0383      &           pTracer(I_obc-1,j,k,bi,bj,iTracer)
                0384      &           *_maskW(I_obc,j,k,bi,bj)
                0385            ENDIF
                0386           ENDDO
                0387          ENDDO
                0388         ENDDO
                0389        ENDIF
                0390 # endif /* ALLOW_OBCS_EAST */
                0391 
                0392 C ------------------------------------------------------------------------------
                0393 
                0394 # ifdef ALLOW_OBCS_WEST
                0395 C     Western OB
                0396 #  ifdef ALLOW_DEBUG
                0397        IF (debugMode)
                0398      &      CALL DEBUG_MSG('OBCS_CALC: West, pTracers',myThid)
                0399 #  endif
                0400        IF (useOrlanskiWest) THEN
                0401         WRITE(msgBuf,'(A)')
                0402      &       'OBCS_CALC: ERROR: useOrlanskiWest Rad OBC with'
                0403         CALL PRINT_ERROR( msgBuf, myThid )
                0404         WRITE(msgBuf,'(A)')
                0405      &       'OBCS_CALC: ERROR: pTracers not yet implemented'
                0406         CALL PRINT_ERROR( msgBuf, myThid )
                0407         STOP 'ABNORMAL END: S/R OBCS_CALC'
                0408        ELSE
                0409         DO iTracer=1,PTRACERS_numInUse
                0410          DO k=1,Nr
                0411           DO j=1-OLy,sNy+OLy
                0412            IF ( OB_Iw(j,bi,bj).NE.OB_indexNone ) THEN
                0413             I_obc = OB_Iw(j,bi,bj)
                0414             OBWptr(j,k,bi,bj,iTracer) =
                0415      &           pTracer(I_obc+1,j,k,bi,bj,iTracer)
                0416      &           *_maskW(I_obc+1,j,k,bi,bj)
                0417            ENDIF
                0418           ENDDO
                0419          ENDDO
                0420         ENDDO
                0421        ENDIF
                0422 # endif /* ALLOW_OBCS_WEST */
                0423 
                0424 C ------------------------------------------------------------------------------
                0425 
                0426 # ifdef ALLOW_OBCS_NORTH
                0427 C         Northern OB
                0428 #  ifdef ALLOW_DEBUG
                0429        IF (debugMode)
                0430      &     CALL DEBUG_MSG('OBCS_CALC: North, pTracers',myThid)
                0431 #  endif
                0432        IF (useOrlanskiNorth) THEN
                0433         WRITE(msgBuf,'(A)')
                0434      &       'OBCS_CALC: ERROR: useOrlanskiNorth Rad OBC with'
                0435         CALL PRINT_ERROR( msgBuf, myThid )
                0436         WRITE(msgBuf,'(A)')
                0437      &       'OBCS_CALC: ERROR: pTracers not yet implemented'
                0438         CALL PRINT_ERROR( msgBuf, myThid )
                0439         STOP 'ABNORMAL END: S/R OBCS_CALC'
                0440        ELSE
                0441         DO iTracer=1,PTRACERS_numInUse
                0442          DO k=1,Nr
                0443           DO i=1-OLx,sNx+OLx
                0444            IF ( OB_Jn(i,bi,bj).NE.OB_indexNone ) THEN
                0445             J_obc = OB_Jn(i,bi,bj)
                0446             OBNptr(i,k,bi,bj,iTracer) =
                0447      &           pTracer(i,J_obc-1,k,bi,bj,iTracer)
                0448      &           *_maskS(i,J_obc,k,bi,bj)
                0449            ENDIF
                0450           ENDDO
                0451          ENDDO
                0452         ENDDO
                0453        ENDIF
                0454 # endif /* ALLOW_OBCS_NORTH */
                0455 
                0456 C ------------------------------------------------------------------------------
                0457 
                0458 # ifdef ALLOW_OBCS_SOUTH
                0459 C         Southern OB
                0460 # ifdef ALLOW_DEBUG
                0461        IF (debugMode)
                0462      &      CALL DEBUG_MSG('OBCS_CALC: South, pTracers',myThid)
                0463 #endif
                0464        IF (useOrlanskiSouth) THEN
                0465         WRITE(msgBuf,'(A)')
                0466      &       'OBCS_CALC: ERROR: useOrlanskiSouth Rad OBC with'
                0467         CALL PRINT_ERROR( msgBuf, myThid )
                0468         WRITE(msgBuf,'(A)')
                0469      &       'OBCS_CALC: ERROR: pTracers not yet implemented'
                0470         CALL PRINT_ERROR( msgBuf, myThid )
                0471         STOP 'ABNORMAL END: S/R OBCS_CALC'
                0472        ELSE
                0473         DO iTracer=1,PTRACERS_numInUse
                0474          DO k=1,Nr
                0475           DO i=1-OLx,sNx+OLx
                0476            IF ( OB_Js(i,bi,bj).NE.OB_indexNone ) THEN
                0477             J_obc = OB_Js(i,bi,bj)
                0478             OBSptr(i,k,bi,bj,iTracer) =
                0479      &           pTracer(i,J_obc+1,k,bi,bj,iTracer)
                0480      &           *_maskS(i,J_obc+1,k,bi,bj)
                0481            ENDIF
                0482           ENDDO
                0483          ENDDO
                0484         ENDDO
                0485        ENDIF
                0486 # endif /* ALLOW_OBCS_SOUTH */
                0487 C     end if (usePTracers)
                0488       ENDIF
                0489 #endif /* ALLOW_PTRACERS */
00816bc2b8 Alis*0490 
3a86c9b47d Oliv*0491 C--   end bi,bj loops.
                0492       ENDDO
                0493       ENDDO
                0494 
abfe198bce Mart*0495 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0496 
                0497 #ifdef ALLOW_OBCS_PRESCRIBE
                0498       IF (useOBCSprescribe) THEN
                0499 C--     Calculate future values on open boundaries
                0500 #ifdef ALLOW_DEBUG
                0501         IF (debugMode) CALL DEBUG_CALL('OBCS_PRESCRIBE_READ',myThid)
                0502 #endif
                0503         CALL OBCS_PRESCRIBE_READ( futureTime, futureIter, myThid )
                0504       ENDIF
                0505 #endif /* ALLOW_OBCS_PRESCRIBE */
                0506 
                0507 C ------------------------------------------------------------------------------
                0508 
                0509 #ifdef ALLOW_OBCS_STEVENS
                0510 C     The Stevens (1990) boundary conditions come after reading data from
                0511 C     files, because they use the data to compute a mix of simplified
                0512 C     Orlanski and prescribed boundary conditions
                0513       IF (useStevensNorth.OR.useStevensSouth.OR.
                0514      &     useStevensEast.OR.useStevensWest) THEN
                0515 #ifdef ALLOW_DEBUG
                0516        IF (debugMode) CALL DEBUG_CALL('OBCS_CALC_STEVENS',myThid)
                0517 #endif
                0518        CALL OBCS_CALC_STEVENS( futureTime, futureIter, myThid )
2b169ecc4b Jean*0519       ENDIF
abfe198bce Mart*0520 #endif /* ALLOW_OBCS_STEVENS */
2b169ecc4b Jean*0521 
                0522 #ifdef ALLOW_DEBUG
                0523       IF (debugMode) CALL DEBUG_LEAVE('OBCS_CALC',myThid)
                0524 #endif
00816bc2b8 Alis*0525 #endif /* ALLOW_OBCS */
2b169ecc4b Jean*0526 
00816bc2b8 Alis*0527       RETURN
                0528       END