[petsc-dev] TSComputeRHSJacobianConstant issue

Barry Smith bsmith at mcs.anl.gov
Sat Feb 27 17:21:37 CST 2016


  Jed and Emil,

  Yes, this is a terrible mess. The TSComputeRHSJacobianConstant() will only work for methods that do not change entries in the Jacobian matrix. And for some reason the runex3_4 and _5 which do not work are turned off in the makefile !  It looks like Rosw does work because it does not change the matrix? And yet it is still turned off? 

   I don't know what a good fix is. If we want to keep this capability then the TSComputeRHSJacobianConstant() interface has to keep a copy of the original user provided matrix and copy it into the matrices used by the method each time. This may lead to ugly special case code.

   Something has to be done about (simplest thing is just tossing the entire concept) because now there is a very broken API in TS that beginning users are likely to try and get upset with.

  Barry


> On Feb 26, 2016, at 11:04 AM, Dominic Meiser <dmeiser at txcorp.com> wrote:
> 
> Hi all,
> 
> I'm running into problems with TS with a time independent, linear
> rhs.  Following the user's manual, I use
> 
> ierr=TSSetProblemType(ts, TS_LINEAR);CHKERRQ(ierr);
> ierr=TSSetRHSFunction(ts, 0, TSComputeRHSFunctionLinear, 0);CHKERRQ(ierr);
> ierr=TSSetRHSJacobian(ts, mat, mat, TSComputeRHSJacobianConstant, 0);CHKERRQ(ierr);
> 
> But running the code I get:
> 
> [dmeiser at ivyamd src]$ ./sd -ts_monitor -snes_monitor -ksp_monitor_true_residual -snes_converged_reason -ts_type beuler -ts_max_steps 20 TS dt 0.01 time 0.
>    0 SNES Function norm 1.433485640883e+00 
>      0 KSP preconditioned resid norm 1.433485640883e-02 true resid norm 1.433485640883e+00 ||r(i)||/||b|| 1.000000000000e+00
>      1 KSP preconditioned resid norm 6.528150237447e-18 true resid norm 6.122666587510e-16 ||r(i)||/||b|| 4.271173992185e-16
>    1 SNES Function norm 1.419007435910e+02 
>  Nonlinear solve converged due to CONVERGED_ITS iterations 1
> 1 TS dt 0.01 time 0.01
>    0 SNES Function norm 1.433342292319e+02 
>  Nonlinear solve did not converge due to DIVERGED_LINEAR_SOLVE iterations 0
> [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> [0]PETSC ERROR:  
> [0]PETSC ERROR: TSStep has failed due to DIVERGED_NONLINEAR_SOLVE, increase -ts_max_snes_failures or make negative to attempt recovery
> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
> [0]PETSC ERROR: Petsc Development GIT revision: pre-tsfc-330-gda8116b  GIT Date: 2016-02-23 15:17:17 -0600
> [0]PETSC ERROR: ./sd on a cplx named ivyamd.txcorp.com by dmeiser Fri Feb 26 11:46:29 2016
> [0]PETSC ERROR: Configure options --with-scalar-type=complex --with-x=0 --with-debugging=1 --with-ssl=0 --with-fortran-kernels=0 --with-pthread=0 --download-fftw=1 --with-mpi=1 --download-mpich=1 --with-mpiuni-fortran-binding=0 --with-fortran-interfaces
> [0]PETSC ERROR: #1 TSStep() line 3346 in /scr_ivyamd/dmeiser/SD/petsc/src/ts/interface/ts.c
> [0]PETSC ERROR: #2 TSSolve() line 3544 in /scr_ivyamd/dmeiser/SD/petsc/src/ts/interface/ts.c
> [0]PETSC ERROR: #3 main() line 63 in /scr_ivyamd/dmeiser/SD/src/sd.c
> [0]PETSC ERROR: PETSc Option Table entries:
> [0]PETSC ERROR: -ksp_monitor_true_residual
> [0]PETSC ERROR: -snes_converged_reason
> [0]PETSC ERROR: -snes_monitor
> [0]PETSC ERROR: -ts_max_steps 2
> [0]PETSC ERROR: -ts_monitor
> [0]PETSC ERROR: -ts_type beuler
> [0]PETSC ERROR: ----------------End of Error Message -------send entire error message to petsc-maint at mcs.anl.gov----------
> application called MPI_Abort(MPI_COMM_WORLD, 91) - process 0
> [unset]: aborting job:
> application called MPI_Abort(MPI_COMM_WORLD, 91) - process 0
> 
> 
> It seems that the linear system is being solved just fine.  But
> apparently it's the wrong linear problem because the SNES
> function norm doesn't decrease.
> 
> I think the problem might be related to the following behaviour
> of `scr/ts/examples/tutorials/ex3.c`:
> 
> 
> [dmeiser at ivyamd tutorials]$ ./ex3 -nox -ts_type beuler -ts_max_steps 2 -ksp_converged_reason -snes_monitor
> Solving a linear TS problem on 1 processor
> Timestep   0: step size = 0.000143637, time = 0., 2-norm error = 1.01507e-15, max norm error = 3.10862e-15
>    0 SNES Function norm 2.018485956680e+03 
>    Linear solve converged due to CONVERGED_RTOL iterations 1
>    1 SNES Function norm 1.207665524289e+05 
> Timestep   1: step size = 0.000143637, time = 0.000143637, 2-norm error = 0.000599907, max norm error = 0.000863771
>    0 SNES Function norm 1.195422385486e+05 
>    Linear solve did not converge due to DIVERGED_PCSETUP_FAILED iterations 0
>                   PCSETUP_FAILED due to FACTOR_NUMERIC_ZEROPIVOT 
> [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> [0]PETSC ERROR:  
> [0]PETSC ERROR: TSStep has failed due to DIVERGED_NONLINEAR_SOLVE, increase -ts_max_snes_failures or make negative to attempt recovery
> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
> [0]PETSC ERROR: Petsc Development GIT revision: pre-tsfc-330-gda8116b  GIT Date: 2016-02-23 15:17:17 -0600
> [0]PETSC ERROR: ./ex3 on a cplx named ivyamd.txcorp.com by dmeiser Fri Feb 26 11:36:28 2016
> [0]PETSC ERROR: Configure options --with-scalar-type=complex --with-x=0 --with-debugging=1 --with-ssl=0 --with-fortran-kernels=0 --with-pthread=0 --download-fftw=1 --with-mpi=1 --download-mpich=1 --with-mpiuni-fortran-binding=0 --with-fortran-interfaces
> [0]PETSC ERROR: #1 TSStep() line 3346 in /scr_ivyamd/dmeiser/SD/petsc/src/ts/interface/ts.c
> [0]PETSC ERROR: #2 TSSolve() line 3544 in /scr_ivyamd/dmeiser/SD/petsc/src/ts/interface/ts.c
> [0]PETSC ERROR: #3 main() line 216 in /scr_ivyamd/dmeiser/SD/petsc/src/ts/examples/tutorials/ex3.c
> [0]PETSC ERROR: PETSc Option Table entries:
> [0]PETSC ERROR: -ksp_converged_reason
> [0]PETSC ERROR: -nox
> [0]PETSC ERROR: -snes_monitor
> [0]PETSC ERROR: -ts_max_steps 2
> [0]PETSC ERROR: -ts_type beuler
> [0]PETSC ERROR: ----------------End of Error Message -------send entire error message to petsc-maint at mcs.anl.gov----------
> application called MPI_Abort(MPI_COMM_WORLD, 91) - process 0
> [unset]: aborting job:
> application called MPI_Abort(MPI_COMM_WORLD, 91) - process 0
> 
> 
> I'm aware of Jed's commit e1244c6 which I thought was supposed to
> enable constant linear rhs without IFunction and IJacobian.  Note
> that e1244c6 also added tests that apparently produced very
> different output (see e.g.
> src/ts/examples/tutorials/output/ex3_4.out).
> 
> Is constant rhs without IFunction and IJacobian no longer
> supported or did it just breat at some point?  Do I need to call
> TSRHSJacobianSetReuse in addition to the lines of code above?  Is
> it better to treat such systems using TSComputeIFunctionLinear
> and TSComputeIJacobianConstant?
> 
> Thanks,
> Dominic
> 
> -- 
> Dominic Meiser
> Tech-X Corporation - 5621 Arapahoe Avenue - Boulder, CO 80303




More information about the petsc-dev mailing list