[MITgcm-support] tracer's time step

Martin Losch mlosch at awi-bremerhaven.de
Thu Nov 24 02:37:07 EST 2005


Hi Eli,

from your STDOUT-files (and probably from the rest of the model output  
as well) you can see, that in both cases the model runs for 160 time  
steps. So it is not surprising, that both runs take the same time. The  
reason why increasing deltaTtracer *alone* does not help is, that  
deltaTclock sets the number of time steps =  
(endTime-startTime)/deltaTclock.
Here is a code snippet from model/src/ini_parms.F, where all deltaT's  
are set according to your namelist specifications:
>       IF ( deltaT       .EQ. 0. ) deltaT       = deltaTmom
>       IF ( deltaT       .EQ. 0. ) deltaT       = deltaTtracer
>       IF ( deltaTmom    .EQ. 0. ) deltaTmom    = deltaT
>       IF ( deltaTtracer .EQ. 0. ) deltaTtracer = deltaT
>       IF ( deltaTClock  .EQ. 0. ) deltaTClock  = deltaT
In you data files you set deltaTmom and deltaTtracer, so that  
deltaTclock=deltaT=deltaTmom.
How do you do it right? That depends on what you want to do:
1. synchronous time stepping (deltaTtracer = deltaTmom = deltaT), which  
is the only physically sensible choice: remove all  
deltaTmom/deltaTclock/deltaTtracer from your data file and use only  
deltaT.
2. tracer acceleration or asynchronos time stepping (Bryan 1984, JPO)  
to avoid CFL-probems, which you should only do for long spin-ups,  
because it changes the propagation speeds of your fast waves ("fast"  
being determined by the ratio of deltaTtracer/deltaTmom)---there's a  
paper by Danabasoglu 2004, Ocean Modelling, where he compares spinups  
with and without tracer acceleration, quite interesting. In that case  
you should do something like this:
  deltaTmom=100.0,
  deltaTtracer=800.0,
  deltaTClock =800.0,
  deltaTfreesurf=800.0,
The last term is a bit unfortunate, because it should default to  
deltaTtracer, but it doesn't for historical reasons? (Maybe it's time  
to change this default?)

Both should have the effect that you expect.

Martin


On Nov 24, 2005, at 6:49 AM, Eli Biton wrote:

> Hi Martin (and others)
> At the attached files you can find the data and STDOUT files for both
> runs.
> Thanks, Eli
>
>>>> mlosch at awi-bremerhaven.de 11/23/05 7:00 PM >>>
> Eli,
> could you show us your data file, or at least the PARM03 section in it?
> Martin
>
> On Nov 23, 2005, at 9:47 AM, eliyahu biton wrote:
>
>> Dear all
>> I have tried to save computation time by increasing my tracer's time
>> step from 100 sec to 800 sec (changing the parameter 'deltaTtracer' at
>
>> the data file under my 'build' directory). The results, however, show
>> no benefit at all in terms of comutation time.
>> Does anyone know if this is normal?
>> thanks alot, Eli
>>
>> _______________________________________________
>> 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
>
> <STDOUT_tarce800.tar.gz><data_trace800><STDOUT_tarce200.tar.gz><data_tr 
> ace200>_______________________________________________
> MITgcm-support mailing list
> MITgcm-support at mitgcm.org
> http://mitgcm.org/mailman/listinfo/mitgcm-support




More information about the MITgcm-support mailing list