[petsc-users] TSSetExactFinalTime Behavior Without TSAdapt

Paul Urbanczyk gomer at stanford.edu
Tue Jun 7 14:49:33 CDT 2016


Hello,

See below for my initial code, before adding TSSetExactFinalTime.

That version yielded the following error:

[0]PETSC ERROR: --------------------- Error Message 
--------------------------------------------------------------
[0]PETSC ERROR: Object is in wrong state
[0]PETSC ERROR: You must call TSSetExactFinalTime() or use 
-ts_exact_final_time <stepover,interpolate,matchstep> before calling 
TSSolve()
[0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html 
for trouble shooting.
[0]PETSC ERROR: Petsc Development GIT revision: v3.7.1-449-g419beca  GIT 
Date: 2016-06-01 15:39:33 -0500
[0]PETSC ERROR: ./urbanSCFD on a arch-linux2-c-debug named prometheus by 
gomer Tue Jun  7 12:42:05 2016
[0]PETSC ERROR: Configure options
[0]PETSC ERROR: #1 TSSolve() line 3946 in 
/home/gomer/local/petsc/src/ts/interface/ts.c

I then added:

TSSetExactFinalTime(ts_solver,TS_EXACTFINALTIME_STEPOVER);

Running then yielded the following error:

[0]PETSC ERROR: 
------------------------------------------------------------------------
[0]PETSC ERROR: Caught signal number 11 SEGV: Segmentation Violation, 
probably memory access out of range
[0]PETSC ERROR: Try option -start_in_debugger or -on_error_attach_debugger
[0]PETSC ERROR: or see 
http://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind
[0]PETSC ERROR: or try http://valgrind.org on GNU/linux and Apple Mac OS 
X to find memory corruption errors
[0]PETSC ERROR: likely location of problem given in stack below
[0]PETSC ERROR: ---------------------  Stack Frames 
------------------------------------
[0]PETSC ERROR: Note: The EXACT line numbers in the stack are not available,
[0]PETSC ERROR:       INSTEAD the line number of the start of the function
[0]PETSC ERROR:       is given.
[0]PETSC ERROR: [0] TSAdaptChoose line 571 
/home/gomer/local/petsc/src/ts/adapt/interface/tsadapt.c
[0]PETSC ERROR: [0] TSStep_Euler line 20 
/home/gomer/local/petsc/src/ts/impls/explicit/euler/euler.c
[0]PETSC ERROR: [0] TSStep line 3707 
/home/gomer/local/petsc/src/ts/interface/ts.c
[0]PETSC ERROR: [0] TSSolve line 3928 
/home/gomer/local/petsc/src/ts/interface/ts.c
[0]PETSC ERROR: --------------------- Error Message 
--------------------------------------------------------------
[0]PETSC ERROR: Signal received
[0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html 
for trouble shooting.
[0]PETSC ERROR: Petsc Development GIT revision: v3.7.1-449-g419beca  GIT 
Date: 2016-06-01 15:39:33 -0500
[0]PETSC ERROR: ./urbanSCFD on a arch-linux2-c-debug named prometheus by 
gomer Tue Jun  7 12:44:17 2016
[0]PETSC ERROR: Configure options
[0]PETSC ERROR: #1 User provided function() line 0 in  unknown file
--------------------------------------------------------------------------

Finally, I added

TSAdapt ts_adapt;
TSGetAdapt(ts_solver,&ts_adapt);
TSAdaptSetType(ts_adapt,TSADAPTNONE);

After adding this, the program runs without error.

Hope that helps.

Please let me know if I was doing something wrong before to have caused 
this. It's not pressing for me, though, as my code now appears to be 
running ok.

-Paul

ORIGINAL VERSION WITHOUT TSSetExactFinalTime:

void cPETScSolverClass::SolveTSProblem()
{
     if(mpi_main.getRank() == 0)
     {
         std::cout << "*****STARTING SOLUTION*****" << std::endl;
     }

     // Set up the time-stepping solver object
     TSCreate(PETSC_COMM_WORLD,&ts_solver);

     // Set the problem type to be nonlinear
     TSSetProblemType(ts_solver,TS_NONLINEAR);

     // Point the solution to the solution vector (vector of conserved 
variables)
     TSSetSolution(ts_solver,solver_main.U_cons);

     // Set the time-stepping type
     TSSetType(ts_solver,TSEULER);

     // Set the initial time and the time step
     TSSetInitialTimeStep(ts_solver,0.0,0.001);

     // Set the total duration (number of steps and total time)
     TSSetDuration(ts_solver,10,0.01);

     // Point the solver to the monitor function
     TSMonitorSet(ts_solver,solver_main.TS_Monitor,this,NULL);

     // Set the distributed array context for the solver
     TSSetDM(ts_solver,main_da);

     // Get the nonlinear solver for the time-stepping solver
     TSGetSNES(ts_solver,&nl_solver);

     // Write out the initial condition before the solver is started
     WriteSolutionCSV("initial_condition.csv");

     // Run the solver!
     TSSolve(ts_solver,solver_main.U_cons);

     if(mpi_main.getRank() == 0)
     {
         std::cout << "-----COMPLETED SOLUTION-----" << std::endl;
     }

     // Write out the solution at the end
     WriteSolutionCSV("output_test.csv");

     // Destroy the time-stepping solver
     TSDestroy(&ts_solver);

}


NEW VERSION WITH TSSetExactFinalTime and TSAdaptSetType:

void cPETScSolverClass::SolveTSProblem()
{
     if(mpi_main.getRank() == 0)
     {
         std::cout << "*****STARTING SOLUTION*****" << std::endl;
     }

     // Set up the time-stepping solver object
     TSCreate(PETSC_COMM_WORLD,&ts_solver);

     // Set the problem type to be nonlinear
     TSSetProblemType(ts_solver,TS_NONLINEAR);

     // Point the solution to the solution vector (vector of conserved 
variables)
     TSSetSolution(ts_solver,solver_main.U_cons);

     // Set the time-stepping type
     TSSetType(ts_solver,TSEULER);

     // Set the initial time and the time step
     TSSetInitialTimeStep(ts_solver,0.0,0.001);

     // Set the total duration (number of steps and total time)
     TSSetDuration(ts_solver,1,0.001);

     // Set TS exact final time behavior
     TSSetExactFinalTime(ts_solver,TS_EXACTFINALTIME_STEPOVER);

     // Set TS adaptive time stepping behavior
     TSAdapt ts_adapt;
     TSGetAdapt(ts_solver,&ts_adapt);
     TSAdaptSetType(ts_adapt,TSADAPTNONE);

     // Point the solver to the function which evaluates the right hand side
TSSetRHSFunction(ts_solver,solver_main.RHS,solver_main.FormFunction,this);

     // Point the solver to the monitor function
     TSMonitorSet(ts_solver,solver_main.TS_Monitor,this,NULL);

     // Set the distributed array context for the solver
     TSSetDM(ts_solver,main_da);

     // Get the nonlinear solver for the time-stepping solver
     TSGetSNES(ts_solver,&nl_solver);

     // Write out the initial condition before the solver is started
     WriteSolutionCSV("initial_condition.csv");

     // Run the solver!
     TSSolve(ts_solver,solver_main.U_cons);

     if(mpi_main.getRank() == 0)
     {
         std::cout << "-----COMPLETED SOLUTION-----" << std::endl;
     }

     // Write out the solution at the end
     WriteSolutionCSV("output_test.csv");

     // Destroy the time-stepping solver
     TSDestroy(&ts_solver);

}


On 06/05/2016 06:44 PM, Barry Smith wrote:
>    Paul,
>
>     I don't think this was an expected change. Can you send code or walk us through the steps that led you to this?
>
>    Barry
>
>
>> On Jun 3, 2016, at 10:23 AM, Paul Urbanczyk <gomer at stanford.edu> wrote:
>>
>> Hello,
>>
>> I recently upgraded to the latest version of PETSc (3.7.1) from the git repository.
>>
>> Running my code gave me an error message about TSSolve being in the wrong state and needing to call TSSetExactFinalTime.
>>
>> I added this function call with TS_EXACTFINALTIME_STEPOVER, but when I re-ran my code, I got a seg fault and error messages pointing to TSAdaptChoose.
>>
>> After fiddling around for a while, I found that I now needed to explicitly set the time step adaptation through TSAdaptSetType.
>>
>> Can anyone confirm that this is so? If so, it might be helpful to note that in documentation or error messages somewhere.
>>
>> Thanks,
>>
>> Paul
>>
>>




More information about the petsc-users mailing list