[MITgcm-devel] some idea for pkg/shelfice
Martin Losch
Martin.Losch at awi.de
Tue Apr 9 05:32:49 EDT 2013
Dan, you are right about this.
Jean-Michel, great and thanks for fixing this. If that was the only problem "erroneously related to the NLFS", then that sheds (again) a poor light on me. I probably should have found that problem years ago.
Martin
On Apr 8, 2013, at 11:59 PM, Daniel Goldberg <dngoldberg at gmail.com> wrote:
> OK. In general, one cannot use an AB update on the first timestep, right? Because the tendency from the previous timestep (or two timesteps) is needed?
>
> Dan
>
>
> On Mon, Apr 8, 2013 at 5:46 PM, Jean-Michel Campin <jmc at ocean.mit.edu> wrote:
> Hi Martin,
>
> The issue with the CD-scheme is fixed, I think.
> I was doing a test with very large gradient of Eta (balanced by
> phi0surf) and the 1rst iteration Adams-Bashforth in CD scheme
> was not consistent with the rest of the momentum (no AB on 1rst iter).
>
> Cheers,
> Jean-Michel
>
> On Mon, Apr 08, 2013 at 09:49:08AM -0400, Jean-Michel Campin wrote:
> > Hi Martin,
> >
> > I have started to do (A), and have 2 issues:
> > 1) it does change results, because of machine truncation
> > (conversion shelficeLoadAnomaly --> shelficeMass --> shelficeLoadAnomaly)
> > can always put a test " IF ( SHELFICEloadAnomalyFile .EQ. ' ' ) THEN"
> > but I would prefer not to, and update the isomip results.
> > 2) I am getting many "TAF recomputation warnings" (10 more than before)
> > and not sure how to fix this.
> >
> > Regarding the other things, need to look into the details of realFreshWater
> > code, to see if it has a good chance to works when kSurf <> 1.
> > And will check also the CD-code with NonLin-FreeSurf (should not be too
> > difficult since it happens at the 1rst iteration), but will be easier once
> > I am able to specify the ice-shelf mass.
> >
> > Cheers,
> > Jean-Michel
> >
> > On Thu, Apr 04, 2013 at 04:24:33PM +0200, Martin Losch wrote:
> > > Hi Jean-Michel,
> > >
> > > some useful suggestions!
> > >
> > > I totally agree with A. Do you want to do that?
> > >
> > > B: I found this
> > > > C with "real fresh water flux" (affecting ETAN),
> > > > C there is more to modify
> > > > rFac = 1. _d 0
> > > > IF ( SHELFICEconserve .AND. useRealFreshWaterFlux ) rFac = 0. _d 0
> > > in the code, but rFac is never used. Not sure what happened here. I think we postponed the real freshwater flux business here. If you know how to do that, I am all for it.
> > >
> > > C: I use the CD scheme, because I needed to get all the smoothing I needed (have a look at section 3 of Losch (2008)). I am sure the CD-scheme can be replaced by biharmonic viscosity. I never got the shelfice package to work with the non-linear free surface. even without r*, but I cannot remember if I tried it only with CD turned on. It's probably better to try without.
> > >
> > > The noise that "called for" the CD scheme is due to the step-wise topography (and partial cells). I always planned to try r* to get a (numerically) "smooth" surface, by depressing the surface layer according the SHELFICEloadAnomaly (or better your SHELFICE_Mass). That might solve many numerical problems, e.g. noise, dynamic coupling, and introduce new ones (e.g., layers that get too thin, pressure gradient errors). I just never got around to trying that out.
> > >
> > > Martin
> > >
> > > On Apr 3, 2013, at 3:47 AM, Jean-Michel Campin <jmc at ocean.mit.edu> wrote:
> > >
> > > > Hi Martin,
> > > >
> > > > Some here have been using pkg/shelfice and, in order to make
> > > > more progress, there are few things that could be changed:
> > > >
> > > > A) add a new variable for the mass of ice (per unit area):
> > > > I propose to do this, and to allow to specify (through a file)
> > > > the ice-shelf mass.
> > > > To start with, fully backward compatible with no result changes,
> > > > will allow either to specify:
> > > > a) R_shelfice & SHELFICEloadAnomaly (same as now)
> > > > or b) R_shelfice & SHELFICE_Mass
> > > > and with the relation between the 2:
> > > > SHELFICE_Mass = SHELFICEloadAnomaly/gravity - Ro_surf*rhoConst
> > > > (since Ro_surf is negative).
> > > > At the end of shelfice_thermodynamics.F, will compute SHELFICEloadAnomaly,
> > > > and in case of (a), SHELFICEloadAnomaly from the file will be converted
> > > > (in shelfice_init_varia.F) to SHELFICE_Mass using the same relation.
> > > > Advantages:
> > > > 1) mass of ice is (much) more natural than SHELFICEloadAnomaly;
> > > > + it does not dependent on hFacMin truncation.
> > > > 2) should allow in the future to step forward SHELFICE_Mass.
> > > > 3) they are few places where SHELFICE_Mass can be used
> > > > (e.g., pressure @ the base of the ice-shelf = SHELFICE_Mass*gravity,
> > > > instead of R_shelfice as it is now,
> > > > also total thickness of the ice = SHELFICE_Mass/rho_ice ...)
> > > > but will let you decide what should be modified or not.
> > > >
> > > > B) implement real-fresh-water flux with pkg/shelfice:
> > > > I think it would not be too hard to allow the surface forcing
> > > > to be applied at the surface level (kSurf) even if kSurf <> 1
> > > > This way, if we put the SHELFICE melting into PmE, we will get
> > > > the real-fresh-water formulation with little more modif.
> > > >
> > > > C) current isomip experiment: found a problem when turning on
> > > > the non-lin free-surf (without z*), got some strange (and wrong)
> > > > currents because of the CD-scheme. Would need confirmation,
> > > > Likely due to a Pb in CD-scheme with non-lin free-surf and
> > > > large variations of hFac (+ kSurf <> 1 ?), will need to be confirmed.
> > > > Can you remember us why CD-scheme is used (since the resolution
> > > > is high enough, ~10.km, to run without) ?
> > > >
> > > > Cheers,
> > > > Jean-Michel
> > > >
> > > >
> > > > _______________________________________________
> > > > 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
>
>
>
> --
>
> Daniel Goldberg, PhD
> NSF Postdoctoral Fellow
> Earth, Atmosphere and Planetary Sciences, MIT
> Cambridge, MA 02139
>
> em: dgoldber at mit.edu
> web: http://ocean.mit.edu/~dgoldberg
> _______________________________________________
> MITgcm-devel mailing list
> MITgcm-devel at mitgcm.org
> http://mitgcm.org/mailman/listinfo/mitgcm-devel
More information about the MITgcm-devel
mailing list