[MITgcm-support] Stratification creating current without forcing (Martin Losch)

Holland, Paul R. pahol at bas.ac.uk
Wed Nov 27 04:36:19 EST 2019


Hi Nadine,

Do your salinities go lower than any of those that you have set in your initial conditions?  That would help track down the problem from a dynamic issue like the one Martin suggests, to a problem with your intial scalars.

Does your hydrogSaltFile somehow contain zeros for the partial cells, or any of the cells nearby?  I would test imposing the same stratification through Sref instead.

Cheers

Paul

-----Original Message-----
From: MITgcm-support <mitgcm-support-bounces at mitgcm.org> On Behalf Of mitgcm-support-request at mitgcm.org
Sent: 27 November 2019 09:25
To: mitgcm-support at mitgcm.org
Subject: MITgcm-support Digest, Vol 197, Issue 16

Message: 3
Date: Wed, 27 Nov 2019 10:05:58 +0100
From: Martin Losch <Martin.Losch at awi.de>
To: MITgcm Support <mitgcm-support at mitgcm.org>
Subject: Re: [MITgcm-support] Stratification creating current without
forcing
Message-ID: <0C12136B-19E2-4AA5-8759-DC18C97C3E0F at awi.de>
Content-Type: text/plain; charset="utf-8"

Hi Nadine,

when you add stratification and still want to have an ocean at rest without any forcing (this includes no meltwater forcing, i.e. SHELFICEheatTransCoeff=0., SHELFICEsaltTransCoeff = 0. or similarly shiCdrag=0), you need to make sure that your iceshelf is in exact hydrostatic balance with the ocean (by specifying an appropriate SHELFICEloadAnomalyFile). This is a little tricky with a non-linear equation of state, but you can find an example in verification/isomip/input/gendata.m how to create data for the  SHELFICEloadAnomalyFile = 'phi0surf.exp1.jmd95z?, in data.shelfice. This file needs to be regenerated each time you change your stratification or your equation of state.

Martin

> On 27. Nov 2019, at 09:45, Nadine Steiger <Nadine.Steiger at uib.no> wrote:
>
> Dear all,
>
> I am using an idealized setup of a v-shaped channel in MITgcm with dimensions of 250 km x 140 km x 0.8 km.
> So far, I had no stratification in the water and added now a linear stratification in salinity (still constant temperature).
>
> To test the setup, I started without any forcing and heat/momentum source, to test if the water stays in rest. Surprisingly, it does not and I get a circulation that increases in strength with time. In the model runs without stratification everything worked fine and the ocean stayed in rest.
>
> I checked the first time steps and it looks like the salinity reduces at the bottom along the topography, where I also get non-zero velocities. It seems like the reduction in salinity happens exactly in the cells where HFacC is between 0 and 1.
>
> Does anyone know what is happening there?
>
> Thanks a lot,
> Nadine
>
>
> This is my data-file:
>
> # Continuous equation parameters
> &PARM01
> Tref = 160*-1.,
> Sref = 160*34.5,
> viscAr=1.E-4,
> viscAh=25.,
> viscA4=0.E-10,
> no_slip_sides=.FALSE.,
> no_slip_bottom=.FALSE.,
> diffKrT=1.E-5,
> diffKhT=5.,
> diffKrS=1.E-5,
> diffKhS=5.,
> bottomDragQuadratic=0.,
> eosType='JMD95Z',
> saltAdvScheme=33,
> tempAdvScheme=33,
> HeatCapacity_cp = 3974.0,
> rhoConst=1027.6,
> rhoNil=1027.6,
> gravity=9.81,
> f0=-1.398E-04,
> beta=0.E-11,
> tAlpha=3.67E-5,
> sBeta=7.89E-4,
> # tAlpha=2.E-4,
> # sBeta =0.E-4,
> convertFW2Salt = 33.4,
> rigidLid=.FALSE.,
> implicitFreeSurface=.TRUE.,
> nonlinFreeSurf = 0.,
> exactConserv=.TRUE.,
> staggerTimeStep=.TRUE.,
> useRealFreshWaterFlux = .TRUE.,
> select_rStar = 0.,
> useSingleCpuIO=.TRUE.,
> hFacMin = 0.5,
> hFacMindr = 0.0,
> hFacInf = 0.2,
> hFacSup = 2.0,
> nonHydrostatic=.FALSE.,
> readBinaryPrec=64,
> useCDScheme = .FALSE.,
> &
>
> # Elliptic solver parameters
> &PARM02
> cg2dMaxIters=1000,
> cg2dTargetResidual=1.E-13,
> &
>
> # Time stepping parameters
> &PARM03
> nIter0=0,
> nTimeSteps=60000.,
> deltaT=100.,
> abEps=0.1,
> cAdjFreq = -1.,
> #tauCD = 400000.,
> pChkptFreq=10000.0,
> chkptFreq=1000.0,
> # dumpFreq=50.,
> taveFreq=0.0,
> monitorFreq=1.,
> monitorSelect=2,
> &
>
> # Gridding parameters
> &PARM04
> usingCartesianGrid=.TRUE.,
> delX=500*500.,
> delY=240*500.,
> delR=160*5.,
> &
>
> # Input datasets
> &PARM05
> bathyFile='topog.channel',
> hydrogThetaFile='T0.bin',
> hydrogSaltFile='S0.bin',
> &
>
> --
> Nadine Steiger
> PhD Candidate
> University of Bergen | Geophysical Institute
> Bjerknes Centre for Climate Research
> _______________________________________________
> MITgcm-support mailing list
> MITgcm-support at mitgcm.org
> http://mailman.mitgcm.org/mailman/listinfo/mitgcm-support



This email and any attachments are intended solely for the use of the named recipients. If you are not the intended recipient you must not use, disclose, copy or distribute this email or any of its attachments and should notify the sender immediately and delete this email from your system.
UK Research and Innovation has taken every reasonable precaution to minimise risk of this email or any attachments containing viruses or malware but the recipient should carry out its own virus and malware checks before opening the attachments. UK Research and Innovation does not accept any liability for any losses or damages which the recipient may sustain due to presence of any viruses.
Opinions, conclusions or other information in this message and attachments that are not related directly to UK Research and Innovation business are solely those of the author and do not represent the views of UK Research and Innovation.



More information about the MITgcm-support mailing list