[MITgcm-support] Way around ptracers & Orlanski OBCS?
Martin Losch
Martin.Losch at awi.de
Thu Oct 4 04:59:37 EDT 2012
Hi Jessica,
I think you have two options:
1. implement the orlanski bc for passive tracers (by copying what it done for temperature and salinity)
2. ignore this issue and comment out the corresponding code in obcs_calc.F that stops the model (but keep the default ptracers, e.g. for the northern boundary:
# ifdef ALLOW_OBCS_NORTH
C Northern OB
# ifdef ALLOW_DEBUG
IF (debugMode)
& CALL DEBUG_MSG('OBCS_CALC: North, pTracers',myThid)
# endif
C IF (useOrlanskiNorth) THEN
C WRITE(msgBuf,'(A)')
C & 'OBCS_CALC: ERROR: useOrlanskiNorth Rad OBC with'
C CALL PRINT_ERROR( msgBuf, myThid )
C WRITE(msgBuf,'(A)')
C & 'OBCS_CALC: ERROR: pTracers not yet implemented'
C CALL PRINT_ERROR( msgBuf, myThid )
C STOP 'ABNORMAL END: S/R OBCS_CALC'
C ELSE
DO iTracer=1,PTRACERS_numInUse
DO k=1,Nr
DO i=1-OLx,sNx+OLx
IF ( OB_Jn(i,bi,bj).NE.OB_indexNone ) THEN
J_obc = OB_Jn(i,bi,bj)
OBNptr(i,k,bi,bj,iTracer) =
& pTracer(i,J_obc-1,k,bi,bj,iTracer)
& *_maskS(i,J_obc,k,bi,bj)
ENDIF
ENDDO
ENDDO
ENDDO
C ENDIF
# endif /* ALLOW_OBCS_NORTH */
) and keep your fingers crossed. As the tracers are passive, you'll have not dynamical problems, but you will not have the same behavior of ptracers and T/S near the boundary. If you can live with this ...
Martin
On Oct 3, 2012, at 1:58 AM, Jessica Spurgin wrote:
> Hello there,
>
> I am trying to use the mitgcm to look at vertical movement of water around a submarine canyon. My domain is set such that I have closed boundaries to the east and west (onshore and offshore), and Orlanksi open boundary conditions to the north & south (in the alongshore). I am using a wind stress to drive circulation in the alongshore direction of the upper layer.
>
> I would also like to apply the ptracers package while using the Orlanski condition. I know the code is such that these two do not work together. I was just wondering if anybody has found a way around this, or have an idea of another way to track movement of certain water properties while using Orlanski boundaries?
>
> Thanks so much,
> Jessica
>
> _______________________________________________
> MITgcm-support mailing list
> MITgcm-support at mitgcm.org
> http://mitgcm.org/mailman/listinfo/mitgcm-support
More information about the MITgcm-support
mailing list