[petsc-users] Pipelined solvers PIPECG etc fail

Sascha Schnepp mail at saschaschnepp.net
Thu Nov 20 08:49:58 CST 2014


Hi,

this is an error that occurs when calling ksp->converged and having 
ksp->normtype == KSP_NORM_NONE

In pipecg.c the error occurs here:
line99:  ierr       = (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); /* test for convergence */

It can be fixed (though there might be a nicer way) using this replacement code, wich explicitly calls KSPConvergendSkip for the respective norm:
  if (ksp->normtype == KSP_NORM_NONE) {
    ierr = KSPConvergedSkip (ksp,0,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr);
  } else {
    ierr = (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr);
  }


The issue could be related to this code in 
KSPCreate_PIPECG
  ierr = KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);CHKERRQ(ierr);
  ierr = KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);CHKERRQ(ierr);
  ierr = KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);CHKERRQ(ierr);
  ierr = KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,2);CHKERRQ(ierr);
where all norms have the same priority.

For the non-pipelined version of CG the unpreconditioned norm has a priority level of 3.

The modifications are independent. Choosing different default priorities will still make line 99 fail if type NONE is chosen explicitly.

Cheers,
Sascha


Begin forwarded message:

> On Thu, Nov 20, 2014 at 4:05 AM, H?kon Strandenes <haakon at hakostra.net>
> wrote:
> 
>> Hi,
>> 
>> First: thanks to the list and Barry for solving my previous problem on
>> loading vectors from HDF5 files. It works great!
>> 
>> Now a second problem: I want to try the pipelined KSP solvers, PIPECG in
>> particular. However, my program fails on the first KSPSolve() with the
>> message below. The same happens if I try PIPECR.
>> 
> 
> These are very new, so there could be bugs lurking. It would help our
> testing if this is
> reproducuble on a standard example. Does this happen if you try to solve
> KSP ex2
> or SNES ex5 with PIPECG?
> 
>   Thanks,
> 
>     Matt
> 
> 
>> I have tried this on three MPI implementations, on and two of them fully
>> support MPI 3 (I have read http://www.mcs.anl.gov/petsc/
>> documentation/faq.html#pipelined). I am using latest PETSc from Git repo.
>> 
>> I have tried to turn off the preconditioner without any effect.
>> 
>> The equation system is the result of a finite-volume discretization of a
>> Poisson-like equation if that helps anything. A normal CG solver works like
>> a charm, converging fast and effective.
>> 
>> Are there any tricks to get this working? Anything special I have to set
>> and/or define to get this solver running?
>> 
>> Thanks in advance.
>> 
>> Best regards,
>> H?kon Strandenes
>> 
>> 
>> Error message:
>> 
>> [0]PETSC ERROR: --------------------- Error Message
>> --------------------------------------------------------------
>> [0]PETSC ERROR: Object is in wrong state
>> [0]PETSC ERROR: Use KSPConvergedSkip() with KSPNormType of KSP_NORM_NONE
>> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html
>> for trouble shooting.
>> [0]PETSC ERROR: Petsc Release Version 3.5.2, unknown
>> [0]PETSC ERROR: madns on a arch-linux2-c-debug named hakostra-kontor by
>> hakostra Thu Nov 20 11:04:28 2014
>> [0]PETSC ERROR: Configure options --prefix=/opt/petsc/git-debug/install
>> --with-shared-libraries --download-hypre --with-64-bit-indices --with-hdf5
>> --with-hdf5-dir=/opt/HDF5/1.8.13/install --with-debugging
>> --with-debugger=gdb
>> [0]PETSC ERROR: #1 KSPConvergedDefault() line 710 in
>> /opt/petsc/git-debug/source/src/ksp/ksp/interface/iterativ.c
>> [0]PETSC ERROR: #2 KSPSolve_PIPECG() line 99 in
>> /opt/petsc/git-debug/source/src/ksp/ksp/impls/cg/pipecg/pipecg.c
>> [0]PETSC ERROR: #3 KSPSolve() line 459 in /opt/petsc/git-debug/source/
>> src/ksp/ksp/interface/itfunc.c
>> [0]PETSC ERROR: #4 findPressure() line 53 in /opt/MaDNS/source/src/FVM2.cpp
>> 
> 



More information about the petsc-users mailing list