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

Mark Lohry mlohry at princeton.edu
Mon Apr 27 10:50:45 CDT 2015


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?



On 04/24/2015 06:41 PM, Barry Smith wrote:
>> On Apr 24, 2015, at 4:08 PM, Mark Lohry <mlohry at princeton.edu> wrote:
>>
>> I'm new to PETSc and not sure how to go about debugging. Running with -start_in_debugger runs in gdb with the normal output (where would i put a breakpoint for this?)
>    You could start by putting a break point in TSAdaptSetType() and see if your call is made and then if any call is made after that by other code changing it to something else. You could put a breakpoint in TSStep_Theta (which supports beuler) and "next" through that to see why it is using an adapter that actually does something.
>
>    Barry
>
>
>
>> and running with -info dumps a lot of information but only related to SNES and KSP, nothing about time stepping or adapt.
>>
>>
>>
>> On 04/24/2015 04:29 PM, Barry Smith wrote:
>>>    Run in the debugger (with one process obviously) to see why it is still calling an adaptor.
>>>
>>>    Barry
>>>
>>>> On Apr 24, 2015, at 3:16 PM, Mark Lohry <mlohry at princeton.edu> wrote:
>>>>
>>>> Trying a simple stiff test case (y'=y^2*(1-y)) using TSBEULER, to get a feel for implicit stepping with JFNK.
>>>>
>>>> I'm giving it a pretty large timestep, and the default appears to start adapting the time step around the stiff part of the solution. So I tried turning it off:
>>>>
>>>>   TSAdapt adapt;
>>>>   TSGetAdapt(its,&adapt);
>>>>   TSAdaptSetType(adapt,TSADAPTNONE);
>>>>
>>>> and even fixing the limits:
>>>>
>>>>   TSAdaptSetStepLimits(adapt,40.0,40.0);
>>>>
>>>> And it still ignores it and adapts anyway. At t=800 (where the solution gets very stiff, so one would normally see adaptive stepping there if it was desired) it's dropping dt from 40 to 10.
>>>>
>>>> Printing (t,y) here from the ts monitor:
>>>>
>>>> 600.000000,0.002790;
>>>> 640.000000,0.003198;
>>>> 680.000000,0.003761;
>>>> 720.000000,0.004606;
>>>> 760.000000,0.006072;
>>>> 800.000000,0.010155;
>>>> 810.000000,0.011452;
>>>> 820.000000,0.013161;
>>>> 830.000000,0.015538;
>>>> 840.000000,0.019126;
>>>> 850.000000,0.025426;
>>>> 860.000000,0.043639;
>>>> 900.000000,0.975519;
>>>> 940.000000,0.999402;
>>>> 980.000000,0.999985;
>>>> 1020.000000,1.000000;
>>>> 1060.000000,1.000000;
>>>> 1100.000000,1.000000;
>>>> 1140.000000,1.000000;
>>>> 1180.000000,1.000000;
>>>>
>>>> Any idea? Adaptive stepping is great, but I'd still like to be able to fix the time step.
>>>>
>>>>
>>>>



More information about the petsc-users mailing list