[petsc-users] Turning off TSADAPT still adapts time step

Mark Lohry mlohry at princeton.edu
Mon Apr 27 11:12:47 CDT 2015


TSAdaptChoose_None is being called, yes. Output with -ts_adapt monitor 
looks like this:

       TSAdapt 'none': step  17 accepted t=680        + 4.000e+01 
family='beuler' scheme=0:'(null)' dt=4.000e+01
       TSAdapt 'none': step  18 accepted t=720        + 4.000e+01 
family='beuler' scheme=0:'(null)' dt=4.000e+01
       TSAdapt 'none': step  19 accepted t=760        + 4.000e+01 
family='beuler' scheme=0:'(null)' dt=4.000e+01
       TSAdapt 'none': step  20 stage rejected t=800        + 4.000e+01 
retrying with dt=1.000e+01
       TSAdapt 'none': step  20 accepted t=800        + 1.000e+01 
family='beuler' scheme=0:'(null)' dt=1.000e+01
       TSAdapt 'none': step  21 stage rejected t=810        + 4.000e+01 
retrying with dt=1.000e+01
       TSAdapt 'none': step  21 accepted t=810        + 1.000e+01 
family='beuler' scheme=0:'(null)' dt=1.000e+01
       TSAdapt 'none': step  22 stage rejected t=820        + 4.000e+01 
retrying with dt=1.000e+01
       TSAdapt 'none': step  22 accepted t=820        + 1.000e+01 
family='beuler' scheme=0:'(null)' dt=1.000e+01
       TSAdapt 'none': step  23 stage rejected t=830        + 4.000e+01 
retrying with dt=1.000e+01
       TSAdapt 'none': step  23 accepted t=830        + 1.000e+01 
family='beuler' scheme=0:'(null)' dt=1.000e+01


I would have expected with "none" that no stages would be rejected, at 
the very least dt would be constant. Is there an error tolerance 
somewhere I need to set to force it to accept all stages?


On 04/27/2015 11:58 AM, Jed Brown wrote:
> Mark Lohry <mlohry at princeton.edu> writes:
>
>> The TSAdaptSetType is being called. Looking at TSStep_Theta, is it
>> possible that here:
>>
>> /* Register only the current method as a candidate because we're not
>> supporting multiple candidates yet. */
>>       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
>>       ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
>>       ierr =
>> TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);CHKERRQ(ierr);
>>       ierr =
>> TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
>>
>> TSAdaptCandidatesClear is unsetting whatever method I'm setting?
> No.  It does not affect the implementation, it just clears the candidate
> list that the implementation may choose from.
>
> PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
> {
>    PetscErrorCode ierr;
>
>    PetscFunctionBegin;
>    PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
>    ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr);
>    PetscFunctionReturn(0);
> }
>
>
> When you set -ts_adapt_type none, the function below should be called.
> Verify that in a debugger.  Also, always run with -ts_adapt_monitor when
> trying to debug something like this.
>
> static PetscErrorCode TSAdaptChoose_None(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte)
> {
>
>    PetscFunctionBegin;
>    *accept  = PETSC_TRUE;
>    *next_sc = 0;                 /* Reuse the same order scheme */
>    *next_h  = h;                 /* Reuse the old step */
>    *wlte    = -1;                /* Weighted local truncation error was not evaluated */
>    PetscFunctionReturn(0);
> }



More information about the petsc-users mailing list