[MITgcm-support] Open boundary conditions for tracers
Martin Losch
mlosch at awi-bremerhaven.de
Mon Nov 1 03:20:54 EST 2004
Hi Erik,
the ptracers-pkg does not yet work with open boundaries, or vice versa,
but I think you can make it work relatively quickly for simple
problems. Here's what I think should be done (no guarantees, but I'll
need the obcs-package for passive tracers as well, so you may see some
the stuff below---if it works for me---in the future):
1. make a copy of obcs_apply_tloc, name it, for example,
obcs_apply_ptracer. The routine routines applies the open boundaries
for temperature (OBN/S/E/Wt) to the temperature field. If you wanted to
generate very general code you'd have to define fields OBN/S/W/Eptr for
each passive tracer, but for a start and for simple boundary
conditions, you could just compute them on the fly in this routine. So
for your box with one open side (I assume for now, that it is the
east-side), it could look like this
> #ifdef ALLOW_OBCS_EAST
> I_obc = OB_Ie(J,bi,bj)
> IF (I_obc.NE.0) THEN
> obc_mask = _maskW(I_obc,J,K,bi,bj)
> tFld(I_obc,J)=tFld(I_obc-1,J)*obc_mask
> ENDIF
> #endif
or whatever you think is approriate (this I really haven't thought
about. If you know the correct boundary condition for tracers leaving
the domain, please tell me).
2. call obcs_apply_ptracer at the right places (where obcs_apply_ts or
obcs_apply_tloc are called). These places are in
a. pkg/generic_advdiff/gad_advection.F, TWO times, just add a line such
as
> ELSEIF (tracerIdentity.GE.GAD_TR1) THEN
> CALL OBCS_APPLY_PTR( bi, bj, k, localTij, myThid )
b. model/src/thermodynamics.F (after timestep_tracer and impldiff).
It's probably best to put the call into ptracers_integrate and
ptracers_impldiff into the ptracers-loop (inside a separate K-loop as
done for T/S in thermodynamics, need to define K).
That should be enough, if (AND ONLY IF) you want to treat all passive
tracers equally. If you want different open boundary conditions for
different passive tracers, the problem is that there is only one
GAD_TR1 (identifier) and so far it's not possible for gad_advection to
know which passive tracer is passed. For that, one has to call
gad_advectioin with GAD_TR1 replaced with itracer+2, also
obcs_apply_ptr should include the argument itracer for tracer
identification. In principle, it's all possible (I think) but much more
involved (something like this is done for the non-local transport of
KPP in gad_calc_rhs). You'd also have to define global fields OBNptr
etc. and define or read them in obcs_calc, apply them in
obcs_apply_ptr, etc.
Martin
> The offline tracer packages has not yet been set up to work
> with obcs (open boundaries). I don't know if anyone is planning
> on doing so? I'm not sure how much work it would entail?
>
> steph
>
>
> On Thu, 28 Oct 2004, Ed Hill wrote:
>
> >
> > Hi folks,
> >
> > Heres a forwarded message from a new MITgcm-support list member. Can
> > someone please help?
> >
> > thanks,
> > Ed
> >
> > ===
> >
> >
> _______________________________________________________________________
> _
> > > From: Erik Kvaleberg <kvaleberg at jhu.edu>
> > > To: MITgcm-support at mitgcm.org
> > > Subject: Open boundary conditions for tracers
> > > Date: Thu, 28 Oct 2004 13:14:54 -0400
> > >
> > > Hi,
> > > I'm running the MITgcm as an offline tracer model (Chris set it up
> for
> > > me in September). Currently I'm using a closed domain with a
> constant
> > > velocity field.
> > > I need to set it up so that the tracer can leave the domain at an
> open
> > > boundary, through advection and diffusion. I suspect I have to
> change
> > > something in pkg/obcs ? Hope you can help.
> > > Erik
> > >
> > > ---------------------------------------------------
> > > Erik Kvaleberg
> > >
> > > Dept. of Earth and Planetary Sciences
> > > Johns Hopkins University
> > >
> > > kvaleberg at jhu.edu
> > > ---------------------------------------------------
> >
> >
More information about the MITgcm-support
mailing list