[MITgcm-support] Spurious internal buoyancy minima in long runs
Jody Klymak
jklymak at uvic.ca
Tue May 3 11:58:44 EDT 2011
Hi Christopher.
Are you sure all your overlap specifications are correct for your
advection scheme (i.e. OLy and OLx), and that you don't have pkg/obcs
turned on, and maybe sponging towards Tref? Did you try changing Tref
to something high to see if something is trying to drive your model to
Tref?
Just guessing.
Cheers, Jody
On May 2, 2011, at 16:43 PM, Christopher L. Wolfe wrote:
>
> Hi all,
>
> Here are some interesting new wrinkles to the saga of the spurious
> internal buoyancy minima:
>
> 1. The same thing occurs if I use a single component (linear)
> equation of state with fixed flux boundary conditions, although the
> effect is smaller (~15% of the domain denser than the surface rather
> than 40-50%).
>
> 2. The effect persists if I use a linear centered-difference scheme
> with horizontal diffusion, but it's smaller again than with the
> single component equation of state (~6% of the domain denser than
> the surface). However, the improvement may just be due to the
> addition of horizontal diffusion mixing away the overdense water.
>
> 3. My default domain has a zonally-reentrant channel occupying the
> southern-most 1/8th of the domain. If I block off the channel
> (effectively building a wall across my "Drake Passage"), the problem
> goes away entirely.
>
> So, the problem doesn't seem to be (directly) related to the
> advection scheme, but only appears if I have a fixed flux boundary
> condition *and* a zonally-reentrant channel. Hopefully these clues
> point to the source of the problem ...
>
> Jean-Michel: Have you had a chance to look at the model setup I sent
> you last week?
>
> Thanks,
> Christopher
>
> On Apr 26, 2011, at 4:10 PM, Jean-Michel Campin wrote:
>
>> Hi Christopher,
>>
>> Could you put the set-up in a tar file and put it some-where
>> where I can get it ? or if it's not too big, send it to me by
>> email (there is a strict limit-size on mitgcm-support).
>> When I have some time, can have a look at it.
>>
>> Thanks,
>> Jean-Michel
>>
>> On Tue, Apr 26, 2011 at 03:59:24PM -0700, Christopher L. Wolfe wrote:
>>>
>>> Hi Jean-Michel,
>>>
>>> cg2d converges to within 1e-7 after about 40 iterations to 1e-11
>>> within about 80 iterations, so I don't think it has to do with
>>> poor convergence of the surface pressure. Also, decreasing the
>>> tolerance to 1e-11 didn't seem to have any impact on the spurious
>>> internal minima.
>>>
>>> Cheers,
>>> Christopher
>>>
>>> On Apr 25, 2011, at 4:23 PM, Jean-Michel Campin wrote:
>>>
>>>> Hi Christopher,
>>>>
>>>> So we know it's not coming from GM.
>>>>
>>>> I am not sure about using only the default params for PARAM02:
>>>>>>> # Elliptic solver parameters
>>>>>>> &PARM02
>>>>>>> &
>>>>
>>>> Sometimes (but it's only scpecial cases), with simple domains
>>>> (cartesian grid with flat bottom) the 2-d solver does not
>>>> converge well.
>>>> The default we have are:
>>>> cg2dMaxIters = 150
>>>> cg2dTargetResidual = 1. _d -7
>>>> And may be you reach the MaxIters number without good convergence ?
>>>> ( using some monitor output and default value for debugLevel
>>>> should help to figure out if this is the case).
>>>>
>>>> An other thing related to the 1rst one: If the solver does not
>>>> converge well (could also be because the Target-residual is too
>>>> large), setting "exactConserv=.TRUE.," might make a difference
>>>> (whereas in general, it does not make much difference).
>>>>
>>>> Cheers,
>>>> Jean-Michel
>>>>
>>>> On Wed, Apr 20, 2011 at 08:52:07PM -0700, Christopher L.P. Wolfe
>>>> wrote:
>>>>>
>>>>> Hi Jean-Michel,
>>>>>
>>>>> Here's my data.gmredi. I've already got GM_AdvForm=.TRUE.
>>>>>
>>>>> &GM_PARM01
>>>>> GM_background_K = 400,
>>>>> GM_taper_scheme = 'dm95',
>>>>> GM_AdvForm = .TRUE.,
>>>>> GM_UseBVP = .TRUE.,
>>>>> GM_BVP_ModeNumber = 2,
>>>>> GM_BVP_cMin = .1,
>>>>> &
>>>>>
>>>>> Thanks,
>>>>> Christopher
>>>>>
>>>>>
>>>>>
>>>>> On Apr 20, 2011, at 6:34 PM, Jean-Michel Campin wrote:
>>>>>
>>>>>> Hi Christopher,
>>>>>>
>>>>>> Might be usefull to also send your data.gmredi file.
>>>>>> Sometime turning on "GM_AdvForm=.TRUE.," helps.
>>>>>>
>>>>>> Cheers,
>>>>>> Jean-Michel
>>>>>>
>>>>>> On Wed, Apr 20, 2011 at 06:04:15PM -0700, Christopher L. Wolfe
>>>>>> wrote:
>>>>>>> support] Spurious internal buoyancy minima in long runs
>>>>>>> X-SA-Exim-Version: 4.2.1 (built Sun, 08 Nov 2009 07:31:22 +0000)
>>>>>>> X-SA-Exim-Scanned: Yes (on ocean.mit.edu)
>>>>>>> Status: O
>>>>>>> Content-Length: 4534
>>>>>>> Lines: 150
>>>>>>>
>>>>>>> Hi Modelers,
>>>>>>>
>>>>>>> I'm having a bit of a strange problem the MITgcm. I'm running
>>>>>>> at coarse resolution (80 km) in a simple domain (a 2000 km by
>>>>>>> 8000 km by 2000 m box) with a linear equation of state, GM/
>>>>>>> Redi eddies, and a "mixed layer" modeled by specifying a large
>>>>>>> diffusivity near the surface which smoothly decays into the
>>>>>>> interior. I'm forcing it with smooth winds, temperature
>>>>>>> restoring (9 days) to a fixed profile, and a fixed and
>>>>>>> balanced surface salt flux.
>>>>>>>
>>>>>>> The problem is that I keep ending up with bottom water which
>>>>>>> is denser than any water found at the surface. Both the
>>>>>>> temperature and salinity fields are bounded by their surface
>>>>>>> values, which is expected since both tracers satisfy advection-
>>>>>>> diffusion equations with extremum principles which state that
>>>>>>> any extrema are found at the surface. With a linear equation
>>>>>>> of state, the buoyancy satisfies a similar advection-diffusion
>>>>>>> equation with a similar extremum principle, yet the buoyancy
>>>>>>> of the bottom water is persistently lower than any water found
>>>>>>> at the surface. This is a serious problem, since I'm trying to
>>>>>>> study deep stratification and it's hard to draw any useful
>>>>>>> conclusions with a huge volume of my domain filled with an
>>>>>>> "impossible" water mass.
>>>>>>>
>>>>>>> I don't think this is a problem with my initial conditions,
>>>>>>> since I've integrated the model for almost 17,000 years with
>>>>>>> an interior diffusivity of 8e-5. I've tried using different
>>>>>>> advection schemes (2, 7, and 81) and en/disabling the
>>>>>>> "Smolarkiewicz hack" without any change in the results. For
>>>>>>> the nonlinear advection schemes, I've also varied the
>>>>>>> horizontal diffusivity from 0 to 100 with no effect.
>>>>>>>
>>>>>>> I'm very much bothered by this problem and have run out of
>>>>>>> things I can think of trying. If anyone here has any
>>>>>>> suggestions/explanations, I'd be extremely grateful to hear
>>>>>>> them. I'm enclosing my data file below for reference.
>>>>>>>
>>>>>>> Thanks,
>>>>>>> Christopher
>>>>>>>
>>>>>>> # ====================
>>>>>>> # | Model parameters |
>>>>>>> # ====================
>>>>>>> #
>>>>>>> # Continuous equation parameters
>>>>>>> &PARM01
>>>>>>> sRef=20*35,
>>>>>>> tRef=20*0.0,
>>>>>>> viscAh=10E3,
>>>>>>> viscAz=0.25E-3,
>>>>>>> no_slip_sides=.TRUE.,
>>>>>>> no_slip_bottom=.FALSE.,
>>>>>>> diffK4S=0.0,
>>>>>>> diffKhS=0.0,
>>>>>>> diffKzS=0.0,
>>>>>>> diffK4T=0.0,
>>>>>>> diffKhT=0.0,
>>>>>>> diffKzT=0.0,
>>>>>>> f0=-1.3683e-04,
>>>>>>> beta=3.4208e-11,
>>>>>>> tAlpha=2.E-4,
>>>>>>> sBeta =7.4e-4,
>>>>>>> gravity=10.,
>>>>>>> rhoConst=1000.,
>>>>>>> rhoNil=1000.,
>>>>>>> rigidLid=.FALSE.,
>>>>>>> implicitFreeSurface=.TRUE.,
>>>>>>> saltAdvection=.TRUE.,
>>>>>>> saltForcing=.TRUE.,
>>>>>>> saltStepping=.TRUE.,
>>>>>>> tempAdvection=.TRUE.,
>>>>>>> tempForcing=.TRUE.,
>>>>>>> tempStepping=.TRUE.,
>>>>>>> eosType='LINEAR',
>>>>>>> nonHydrostatic=.FALSE.,
>>>>>>> momAdvection=.TRUE.,
>>>>>>> implicitViscosity=.TRUE.,
>>>>>>> implicitDiffusion=.TRUE.,
>>>>>>> ivdc_kappa=10.,
>>>>>>> readBinaryPrec=64,
>>>>>>> writeBinaryPrec=32,
>>>>>>> tempAdvScheme=81,
>>>>>>> saltAdvScheme=81,
>>>>>>> staggerTimeStep=.TRUE.,
>>>>>>> bottomDragLinear=1.1135E-3,
>>>>>>> debugLevel=-1,
>>>>>>> useJamartWetPoints=.TRUE.,
>>>>>>> useSingleCpuIo=.TRUE.,
>>>>>>> diffKrNrS=4.069423e-03,
>>>>>>> 3.808919e-03,
>>>>>>> 2.924517e-03,
>>>>>>> 1.608878e-03,
>>>>>>> 5.410003e-04,
>>>>>>> 1.345410e-04,
>>>>>>> 8.147859e-05,
>>>>>>> 8.000402e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> diffKrNrT=4.069423e-03,
>>>>>>> 3.808919e-03,
>>>>>>> 2.924517e-03,
>>>>>>> 1.608878e-03,
>>>>>>> 5.410003e-04,
>>>>>>> 1.345410e-04,
>>>>>>> 8.147859e-05,
>>>>>>> 8.000402e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> 8.000000e-05,
>>>>>>> &
>>>>>>>
>>>>>>> # Elliptic solver parameters
>>>>>>> &PARM02
>>>>>>> &
>>>>>>>
>>>>>>> # Time stepping parameters
>>>>>>> &PARM03
>>>>>>> nIter0=0006205000,
>>>>>>> # 5000 years
>>>>>>> nTimeSteps=1825000,
>>>>>>> deltaTMom=3600,
>>>>>>> deltaTtracer=86400,
>>>>>>> abEps=0.1,
>>>>>>> # write out every 50 years
>>>>>>> pChkptFreq=1.5768e10,
>>>>>>> chkptFreq=1.5768e9,
>>>>>>> dumpFreq=1.5768e9,
>>>>>>> monitorFreq=1.5768e7,
>>>>>>> tauThetaClimRelax=764400.,
>>>>>>> cAdjFreq=0,
>>>>>>> pickupStrictlyMatch=.FALSE.,
>>>>>>> &
>>>>>>>
>>>>>>> # Gridding parameters
>>>>>>> &PARM04
>>>>>>> usingCartesianGrid=.TRUE.,
>>>>>>> usingSphericalPolarGrid=.FALSE.,
>>>>>>> dXspacing=80.e3,
>>>>>>> dYspacing=80.e3,
>>>>>>> delZ= 12.0505, 14.9159, 18.4328, 22.7337, 27.9692,
>>>>>>> 34.3067, 41.9245, 51.0025, 61.7057, 74.1598,
>>>>>>> 88.4181, 104.4189, 121.9376, 140.5412, 159.5568,
>>>>>>> 178.0710, 194.9769, 209.0770, 219.2372, 224.5642,
>>>>>>> &
>>>>>>>
>>>>>>> # Input datasets
>>>>>>> &PARM05
>>>>>>> bathyFile='topo_two.bin',
>>>>>>> thetaClimFile='TSurf.bin',
>>>>>>> zonalWin
>>>>>>
>>>>>> _______________________________________________
>>>>>> MITgcm-support mailing list
>>>>>> MITgcm-support at mitgcm.org
>>>>>> http://mitgcm.org/mailman/listinfo/mitgcm-support
>>>>>
>>>>> -----------------------------------------------------------
>>>>> Dr. Christopher L. Wolfe 858-534-4560
>>>>> Physical Oceanography Research Division OAR 357
>>>>> Scripps Institution of Oceanography, UCSD clwolfe at ucsd.edu
>>>>> -----------------------------------------------------------
>>>>>
>>>>>
>>>>>
>>>>>
>>>>> _______________________________________________
>>>>> MITgcm-support mailing list
>>>>> MITgcm-support at mitgcm.org
>>>>> http://mitgcm.org/mailman/listinfo/mitgcm-support
>>>>
>>>> _______________________________________________
>>>> MITgcm-support mailing list
>>>> MITgcm-support at mitgcm.org
>>>> http://mitgcm.org/mailman/listinfo/mitgcm-support
>>>
>>>
>>> _______________________________________________
>>> MITgcm-support mailing list
>>> MITgcm-support at mitgcm.org
>>> http://mitgcm.org/mailman/listinfo/mitgcm-support
>>
>> _______________________________________________
>> MITgcm-support mailing list
>> MITgcm-support at mitgcm.org
>> http://mitgcm.org/mailman/listinfo/mitgcm-support
>
>
> _______________________________________________
> MITgcm-support mailing list
> MITgcm-support at mitgcm.org
> http://mitgcm.org/mailman/listinfo/mitgcm-support
--
Jody Klymak
http://web.uvic.ca/~jklymak/
More information about the MITgcm-support
mailing list