[MITgcm-devel] seaice

chris hill cnh at mit.edu
Tue Feb 14 10:30:39 EST 2006


Martin,

  I took a quick skim through the code. You could try setting advfac to 
0.. That would mean any convergence in the ice velocity will cause ice 
to pile up (heff gets bigger) and divergence will cause it to thin (heff 
get smaller). That seems like it would make sense to me. As you have it 
now a divergent ice flow field may end up producing more ice (since with 
advfac = 1 heff will have a source term that balances the flow 
divergence/convergence). However, I think we need to get Jinlun's expert 
thoughts, as the end result is tied up with how the solver works with 
heff, uice and vice, whether uice and vice can be non-divergent etc...

Chris
Martin Losch wrote:
> GAD:
> 
> I have tried using dst3fl (withoug any explicit diffusion). It runs  
> stably for 20 years now, but now my ice (instead of being too little)  
> is growing endlessly (37 m in the western Weddell Sea), so I have the  
> suspicion that I am make severe mistakes: Could someone who  understand 
> gad_calc_rhs have a look at the attached routine and tell  me, if I have 
> done things right? In particular in the update of HEFF,  is it correct 
> to have advFac=1?:
> 
>>         DO j=1-Oly,sNy+Oly-1
>>          DO i=1-Olx,sNx+Olx-1
>>           HEFF(i,j,1,bi,bj)=HEFF(i,j,3,bi,bj) + DELTT *
>>      &         maskC(i,j,kSurface,bi,bj)*recip_rA(i,j,bi,bj)
>>      &       *( (fZon(i+1,j)-fZon(i,j))
>>      &         +(fMer(i,j+1)-fMer(i,j))
>>      &         -localT(i,j)*( (uTrans(i+1,j)-uTrans(i,j))
>>      &                       +(vTrans(i,j+1)-vTrans(i,j))
>>      &                       )*advFac
>>      &         )
> 
> 
> Martin
> 
> 
> 
> On Feb 14, 2006, at 3:57 PM, chris hill wrote:
> 
>> Martin,
>>
>>  These all sound good.
>>  Technically it should be possible to use gad to do the explicit  part 
>> of ice advection. Jinlun may have comments on whether this  makes 
>> sense algorithmically and on interactions with the implicit  parts of 
>> the ice dynamics.
>>
>> Chris
>> Martin Losch wrote:
>>
>>> Hi Dimitris,
>>> I am probably playing too much with the seaice model right now,  but  
>>> here are a few points that I'd like to make (and maybe check  them 
>>> in,  even if they break lab_sea):
>>> 1. I am now sure that there is a bug in advect.F that does NOT  
>>> affect  lat-lon-grid simulations, but WILL affect cubed-sphere  
>>> simulations  and all other irregular grid simulations. It's  
>>> basically an idexing  error (see my previous email). I think I  will 
>>> just fix that.
>>> 2. I would like to replace all DXTICE DYTICE SINEICE CSTICE etc  
>>> with  the proper combination of variables dxF,dxG, etc. from  GRID.h. 
>>> This  will --- at least as far as I can see --- make sure  that the 
>>> grid  information is correct and the same grid parameters  that are 
>>> used for  the ocean are used for seacie.
>>> Since I want to use the seaice model on a cubed sphere grid, I do   
>>> care about this. However, this will change the lab_sea and very  
>>> like  (more dramatically) any cubed sphere set-up that you may  have 
>>> (I am  currently currently playing with global_ocean.cs32x15  + 
>>> seaice). Will  I get your OK?
>>> 3. Advection schemes: for properties such as volume and  fractional  
>>> area, the advection scheme should not produce negative  (or positve)  
>>> overshoots. A 2nd order central difference scheme  does that (eg., 
>>> can  produce negative thicknesses). The scheme in  advect.F is 2nd 
>>> order  central difference, but I don't understand  the time stepping 
>>> scheme,  so it may be OK. Nevertheless, I  naively think, a positive 
>>> scheme may  be better, but it is no  longer conservative, eg. 
>>> 2n-order with flux  limiter (e.g, Hunke's  CSIM5 uses MPDATA) or 
>>> DST3FL that I use  routinely for geochemical  tracers. The nice thing 
>>> is, that all of  these schemes are there  (in generic_advdiff), one 
>>> just needs to pick  one. I have tried  dst3fl, but again, I do not 
>>> understand the time  stepping in  advect.F (nor do I understand fully 
>>> how gad_calc_rhs  works): I  have tried dst3fl and I even got it to 
>>> work, but only  halfway. If  I am not mistaken, the DST schemes look 
>>> as if they are  explicit  in time, that is, h(n+1) = h(n) + 
>>> gh(n)*deltaT. I can  compute gh (n), but for that I need to know what 
>>> the different time  levels  are, eg.,
>>> HEFF(:,:,1,:,:) = current time level?
>>> HEFF(:,:,2,:,:) = do I need these?
>>> HEFF(:,:,3,:,:) = ?
>>> Or do I just update HEFF(:,:,1,:,:) in advect.F?
>>> Martin
>>> On Feb 14, 2006, at 2:07 PM, Dimitris Menemenlis wrote:
>>>
>>>> Martin and Jinlun, I am out of my depth when it comes to  advection  
>>>> schemes.  Is there a reason for changing the scheme  that is there  
>>>> already in pkg/seaice?
>>>>
>>>> For cubed-sphere grid right now, I assumed that grid is  
>>>> rectangular  near the Poles (CS*ICE=1, TNG*ICE=0).  This was a  
>>>> quick fix to get  going but it is not exact.  So maybe that  
>>>> explains why you get  different numerical values?
>>>>
>>>> Regarding coastal sflux from seaice.  One does expect coastal   
>>>> regions around Antarctica to be ice/salt factories, but maybe  too  
>>>> much salt is being rejected.  Carl recently send me some  slides 
>>>> and  Ph.D. thesis from Dirk Notz:
>>>> http://ecco.jpl.nasa.gov/~dimitri/Notz/talk_MPI16112005.pdf
>>>> http://ecco.jpl.nasa.gov/~dimitri/Notz/PhD_thesis_Dirk.pdf
>>>> suggesting there is considerable uncertainty regarding how much   
>>>> salt is rejected during sea ice creation.
>>>>
>>>> Dimitris
>>>> _______________________________________________
>>>> 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




More information about the MITgcm-devel mailing list