[MITgcm-devel] stevens boundary conditions

Jean-Michel Campin jmc at ocean.mit.edu
Mon Dec 12 13:55:45 EST 2011


Hi Martin,

Could you leave this masking commented out (keeping #include "PACKAGES_CONFIG.h")
just in case we bump into an unknow problem one day ?
And adding a comment (related to obcs-stevens) ? 

Cheers,
Jean-Michel

On Mon, Dec 12, 2011 at 06:52:31PM +0100, Martin Losch wrote:
> Hi Jean-Michel,
> 
> I thinks I know what's going on: in timestep_tracer the gtracer is masked by maskInC. That of course spoils my idea. Do you think we can get rid of this mask? Here's your check-in message: "leave tracer unchanged outside OB interior region: This has no effect
>  on the solution but just to prevent unrealistic tracer value ouside OB."
> 
> Martin
> 
> On Dec 12, 2011, at 5:52 PM, Martin Losch wrote:
> 
> > Hi Jean-Michel,
> > 
> > before I check in my modified obcs_calc_stevens, I need to understand, why I do not see any effect  of those changes  in exp4. 
> > I expect for the western boundary of exp4 (i=1) the advective fluxes "af" (even in the diagnostic DFxE_TH) are zero, but this does not seem to be so, why? Instead, they even vary with time, what am I missing? I tried tempAdvScheme 2, 4, 33
> > 
> > Martin
> > 
> > On Dec 12, 2011, at 3:43 PM, Jean-Michel Campin wrote:
> > 
> >> Martin,
> >> 
> >> It's not a bad idea to put those masks in gmredi routines, I agree. 
> >>> so I am ready to check this in, if you do not object.
> >> Yes, go for it.
> >> 
> >> Jean-Michel
> >> 
> >> On Mon, Dec 12, 2011 at 03:34:38PM +0100, Martin Losch wrote:
> >>> Hi Jean-Michel,
> >>> 
> >>> rather than introducing another df-field, I decided to apply the maskInC directly in the gmredi routines (gmredi_rtransport and gmredi_calc_diff). All the changes that I did to implement option 2 did not change any results, so I am ready to check this in, if you do not object.
> >>> 
> >>> Martin
> >>> 
> >>> Index: generic_advdiff/gad_advection.F
> >>> ===================================================================
> >>> RCS file: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_advection.F,v
> >>> retrieving revision 1.69
> >>> diff -r1.69 gad_advection.F
> >>> 810c810
> >>> <      &        )*rkSign
> >>> ---
> >>>>    &        )*rkSign*maskInC(i,j,bi,bj)
> >>> Index: generic_advdiff/gad_calc_rhs.F
> >>> ===================================================================
> >>> RCS file: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_calc_rhs.F,v
> >>> retrieving revision 1.58
> >>> diff -r1.58 gad_calc_rhs.F
> >>> 522c522
> >>> <           fVerT(i,j,kUp) = fVerT(i,j,kUp) + af(i,j)
> >>> ---
> >>>>         fVerT(i,j,kUp) = fVerT(i,j,kUp) + af(i,j)*maskInC(i,j,bi,bj)
> >>> 718a719
> >>>> C     for Stevens OBC: keep only vertical diffusive contribution on boundaries
> >>> 724,725c725,726
> >>> <      &   *( (fZon(i+1,j)-fZon(i,j))
> >>> <      &     +(fMer(i,j+1)-fMer(i,j))
> >>> ---
> >>>>    &   *( (fZon(i+1,j)-fZon(i,j))*maskInC(i,j,bi,bj)
> >>>>    &     +(fMer(i,j+1)-fMer(i,j))*maskInC(i,j,bi,bj)
> >>> 730c731
> >>> <      &                  )
> >>> ---
> >>>>    &                  )*maskInC(i,j,bi,bj)
> >>> Index: gmredi/gmredi_calc_diff.F
> >>> ===================================================================
> >>> RCS file: /u/gcmpack/MITgcm/pkg/gmredi/gmredi_calc_diff.F,v
> >>> retrieving revision 1.10
> >>> diff -r1.10 gmredi_calc_diff.F
> >>> 58a59
> >>>>    &           *maskInC(i,j,bi,bj)
> >>> 64a66
> >>>>    &           *maskInC(i,j,bi,bj)
> >>> 75a78
> >>>>    &           *maskInC(i,j,bi,bj)
> >>> 81a85
> >>>>    &           *maskInC(i,j,bi,bj)
> >>> Index: gmredi/gmredi_rtransport.F
> >>> ===================================================================
> >>> RCS file: /u/gcmpack/MITgcm/pkg/gmredi/gmredi_rtransport.F,v
> >>> retrieving revision 1.16
> >>> diff -r1.16 gmredi_rtransport.F
> >>> 133c133
> >>> <      &      - _rA(i,j,bi,bj)
> >>> ---
> >>>>    &      - _rA(i,j,bi,bj)*maskInC(i,j,bi,bj)
> >>> 137c137
> >>> <      &      - _rA(i,j,bi,bj)
> >>> ---
> >>>>    &      - _rA(i,j,bi,bj)*maskInC(i,j,bi,bj)
> >>> 167a168
> >>>>    &        *maskInC(i,j,bi,bj)
> >>> 184a186
> >>>> c    &    *maskInC(i,j,bi,bj)
> >>> 
> >>> 
> >>> On Dec 12, 2011, at 3:03 PM, Jean-Michel Campin wrote:
> >>> 
> >>>> Hi Martin,
> >>>> 
> >>>> I did some tests to check that seaice obcs were clean
> >>>> (but don't remember how I did them !),
> >>>> and did not pass unless I put useKPP=F (would need to
> >>>> investigate within KPP code how to make diffusivity & viscosity
> >>>> independant of value from the other side). I think GGL90
> >>>> might also have similar issues but should be easier to fix.
> >>>> 
> >>>> Cheers,
> >>>> Jean-Michel
> >>>> 
> >>>> On Mon, Dec 12, 2011 at 01:34:39PM +0100, Martin Losch wrote:
> >>>>> On Dec 12, 2011, at 12:41 PM, Martin Losch wrote:
> >>>>> 
> >>>>>>> (in fact, in seaice_obcs test.exp. 
> >>>>>>> which is using KPP, some of the western OB values already 
> >>>>>>> influence the Eastern side).
> >>>>> 
> >>>>> How do you see that? Do I need to modify the exp to see that? Or is the flag OBCS_uvApplyFac enough?
> >>>>> 
> >>>>> M.
> >>>>> 
> >>>>> 
> >>>>> _______________________________________________
> >>>>> MITgcm-devel mailing list
> >>>>> MITgcm-devel at mitgcm.org
> >>>>> http://mitgcm.org/mailman/listinfo/mitgcm-devel
> >>>> 
> >>>> _______________________________________________
> >>>> MITgcm-devel mailing list
> >>>> MITgcm-devel at mitgcm.org
> >>>> http://mitgcm.org/mailman/listinfo/mitgcm-devel
> >>> 
> >>> 
> >>> _______________________________________________
> >>> MITgcm-devel mailing list
> >>> MITgcm-devel at mitgcm.org
> >>> http://mitgcm.org/mailman/listinfo/mitgcm-devel
> >> 
> >> _______________________________________________
> >> MITgcm-devel mailing list
> >> MITgcm-devel at mitgcm.org
> >> http://mitgcm.org/mailman/listinfo/mitgcm-devel
> > 
> > 
> > _______________________________________________
> > MITgcm-devel mailing list
> > MITgcm-devel at mitgcm.org
> > http://mitgcm.org/mailman/listinfo/mitgcm-devel
> 
> Martin Losch
> Martin.Losch at awi.de
> 
> 
> 
> 
> _______________________________________________
> MITgcm-devel mailing list
> MITgcm-devel at mitgcm.org
> http://mitgcm.org/mailman/listinfo/mitgcm-devel



More information about the MITgcm-devel mailing list