[MITgcm-support] Model crashing due to vertical velocities

cimatori cimatori at knmi.nl
Thu Aug 18 07:55:09 EDT 2011


 Indeed switching from a domain along y to a domain along x gives quite
a substantial performance improvement, thanks a lot for the suggestion.
For the stability issue, it seems that the use of an explicit vertical
diffusion was the problem (in combination with the vertical walls).
Removing diffkzT solves it. Anyway, I'll have to experiment with
different advection schemes as I do want an explicit vertical diffusion,
and I would prefer a box domain.

Andrea


On 08/18/2011 10:17 AM, Martin Losch wrote:
> Hi Andrea,
>
> Your vertical diffusion/viscosity appears a little large (in the ocean I would use order 1e-5 and 1e-3 for viscosity, viscAz), but I would expect with your configuration, that nothing moves at all (stable stratification, no forcing), so that's puzzling to begin with.
>
> I'd try
> viscAz = 1.e-3,
> viscAhGrid = 0.1, (that scales the viscosity with grid cell size, interesting, when your grid size varies)
>
> You could simplify your problem even more by having constant dy/dx and dz, remove the boundaries (make it periodic), see what happens, and gradually re-introduce complexity to see where the problem stems from (I am guessing no slip on boundaries).
>
> What happens, when you turn the domain around so that is has ny=1 and nx=1000? (On a vector machine this will be much faster, but even on a single cpu with some vectorization, it will be faster).
>
> Martin
>
> PS. Just some remarks:
>
> With 
> tempAdvScheme=33
> you do not need any explicit diffusion for stability, just viscosity.
>
> You do not need these flags, they have the default values:
> momViscosity=.TRUE.,
> tempAdvection=.TRUE.,
> tempDiffusion=.TRUE.,
> saltStepping=.TRUE.,
> saltAdvection=.TRUE.,
> saltDiffusion=.TRUE.,
>
>
>
>
> On Aug 18, 2011, at 10:35 AM, Andrea Cimatoribus wrote:
>
>>
>> Hi everybody,
>> I'm doing some tests with a 2D model (1 grid cell in the x direction, 1000 in the y direction, non hydrostatic). The resolution goes from approximately 700m to 7km, in the vertical it's from 10m to 40m (the low vertical resolution is to try to keep computational cost low). Just to start, I have no forcing, and an initial stable stratification in temperature (salt is constant).
>> The model is crashing due to very strong vertical velocities that develop close to the north/south boundaries. I did lots of testing with different values for viscosity and timestep (down to 1s), with no success. Anyway, if I change the bathimetry from flat bottom (3000m) to a sufficiently gentle slope, the model runs with no problem. I attach my data file. Any idea of what I'm missing? Could it be connected with the low resolution in the vertical? Before flooding you with figures from the model, I'd like to know if there's something simple that I'm missing.
>> Thanks a lot,
>> AC
>>
>>
>> # ====================
>> # | Model parameters |
>> # ====================
>> #
>> # Continuous equation parameters
>> &PARM01
>> tRef= 20.0, 19.7, 19.4, 19.1, 18.8, 18.4, 18.0, 17.7,
>> 17.3, 16.9, 16.5, 16.0, 15.6, 15.2, 14.8, 14.3,
>> 13.8, 13.4, 13.0, 12.6, 12.2, 11.7, 11.3, 11.0,
>> 10.6, 10.1, 9.7, 9.3, 8.9, 8.6, 8.2, 8.0,
>> 7.7, 7.5, 7.2, 7.0, 6.9, 6.7, 6.6, 6.4,
>> 6.3, 6.2, 6.1, 6.0, 5.9, 5.8, 5.7, 5.68,
>> 5.6, 5.56, 5.5, 5.47, 5.43, 5.39, 5.36, 5.32,
>> 5.30, 5.27, 5.25, 5.22, 5.21, 5.19, 5.17, 5.16,
>> 5.14, 5.13, 5.12, 5.11, 5.10, 5.09, 5.08, 5.07,
>> 5.06, 5.05, 5.04, 5.04, 5.04, 5.03, 5.02, 5.01,
>> sRef=80*30.,
>> viscAh=12,
>> viscAz=3.E-1,
>> momViscosity=.TRUE.,
>> no_slip_sides=.TRUE.,
>> no_slip_bottom=.TRUE.,
>> diffKhT=1.E-1,
>> diffKzT=1.E-1,
>> # diffK4T=1.E-3,
>> # diffK4S=1.E-3,
>> f0=0.E-4,
>> beta=0.E-11,
>> tAlpha=2.0E-4,
>> sBeta =7.4E-4,
>> gravity=10.,
>> rhoConst=1030.,
>> rhoNil=1030.,
>> heatCapacity_Cp=4000.,
>> #rigidLid=.TRUE.,
>> implicitFreeSurface=.TRUE.,
>> # exactConserv=.TRUE.,
>> eosType='LINEAR',
>> nonHydrostatic=.TRUE.,
>> # implicitDiffusion=.TRUE.,
>> tempAdvection=.TRUE.,
>> tempDiffusion=.TRUE.,
>> tempForcing=.FALSE.,
>> saltStepping=.TRUE.,
>> saltAdvection=.TRUE.,
>> saltDiffusion=.TRUE.,
>> saltForcing=.FALSE.,
>> tempAdvScheme=33,
>> # saltAdvScheme=33,
>> staggerTimeStep=.TRUE.,
>> &
>>
>> # Elliptic solver parameters
>> &PARM02
>> cg2dMaxIters=1000,
>> cg2dTargetResidual=1.E-9,
>> cg3dMaxIters=100,
>> cg3dTargetResidual=1.E-9,
>> &
>>
>> # Time stepping parameters
>> &PARM03
>> nIter0=0,
>> endTime=86400.,
>> #nTimeSteps=3,
>> deltaT=10.,
>> pChkptFreq=86400.,
>> chkptFreq=86400.,
>> dumpFreq=30.,
>> monitorFreq=1800.,
>> monitorSelect=1800,
>> tauSaltClimRelax = 864000.,
>> tauThetaClimRelax = 864000.,
>> &
>>
>> # Gridding parameters
>> &PARM04
>> usingCartesianGrid=.TRUE.,
>> dXspacing=1000.,
>> # dYspacing=5000.,
>> delZ=5*10.,5*15.,5*20.,5*25,5*30,55*46
>> delYfile = 'delY.bin',
>> &
>>
>> # Input datasets
>> &PARM05
>> thetaClimFile = 'SurfTemp.bin',
>> saltClimFile = 'SurfSalt.bin',
>> bathyFile = 'bathy.bin',
>> &
>>
>>
>> _______________________________________________
>> 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




More information about the MITgcm-support mailing list