[petsc-dev] problem with LSQR and CGNE

Raphaël Couturier raphael.couturier at univ-fcomte.fr
Thu Sep 3 13:21:55 CDT 2015


Hello all,

First of all, thank you Matthew and Barry for your answers in July.
Your solution using PCKSP is fine.
Now it works perfectly with the examples I have tried for ksp and snes. 
I made many tests with many processors.

I have implemented my own CGLS or LSQR method. In Petsc, I have seen 
that CGNE and LSQR are already implemented. However, I cannot use them 
because I have a problem of non square matrix.
In previous posts, I saw that the problem can come from the preconditioner.
In the attached code, I have my new solver (without preconditioner if I 
made no mistake) which is based on
  ierr = PCKSPGetKSP(pc,&sub_ksp);CHKERRQ(ierr);

and I need to use another ksp to minimize the residual (instead of my 
own CGLS). The matrix AS is of size (n,10) with n very large.
I have tried with that (line 262 in the code)
  {
         KSP ksp3;
         ierr = KSPCreate(PETSC_COMM_WORLD, &ksp3); CHKERRQ(ierr);
         ierr = KSPSetType(ksp3, KSPLSQR); CHKERRQ(ierr);
         ierr = KSPSetOperators(ksp3, AS, AS); CHKERRQ(ierr);
         ierr = KSPSetTolerances(ksp3, PETSC_DEFAULT, PETSC_DEFAULT, 
PETSC_DEFAULT, iter_max_minimization);
         ierr = PCSetType(ksp,PCNONE);
         ierr = KSPSolve(ksp3, b, Alpha); CHKERRQ(ierr);
       }

The problem on all the processors is:
[0]PETSC ERROR: No support for this operation for this object type
[0]PETSC ERROR: Only square matrices supported.
[0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html 
for trouble shooting.
[0]PETSC ERROR: Petsc Release Version 3.6.1, Jul, 22, 2015
[0]PETSC ERROR: ./ex15 on a arch-linux2-c-debug named extinction by 
couturie Thu Sep  3 13:22:03 2015
[0]PETSC ERROR: Configure options --download-openmpi --download-hypre 
--download-f2cblaslapack --with-fc=0
[0]PETSC ERROR: #2 MatGetDiagonalBlock_MPIDense() line 78 in 
/home/couturie/work/petsc-3.6.1/src/mat/impls/dense/mpi/mpidense.c
[0]PETSC ERROR: #3 MatGetDiagonalBlock() line 159 in 
/home/couturie/work/petsc-3.6.1/src/mat/interface/matrix.c
[0]PETSC ERROR: #4 PCSetUp_BJacobi() line 126 in 
/home/couturie/work/petsc-3.6.1/src/ksp/pc/impls/bjacobi/bjacobi.c
[0]PETSC ERROR: #5 PCSetUp() line 982 in 
/home/couturie/work/petsc-3.6.1/src/ksp/pc/interface/precon.c
[0]PETSC ERROR: #6 KSPSetUp() line 332 in 
/home/couturie/work/petsc-3.6.1/src/ksp/ksp/interface/itfunc.c
[0]PETSC ERROR: #7 KSPSolve() line 546 in 
/home/couturie/work/petsc-3.6.1/src/ksp/ksp/interface/itfunc.c
[0]PETSC ERROR: #8 KSPSolve_TSIRM() line 269 in 
/home/couturie/work/petsc-3.6.1/src/ksp/ksp/impls/gmres/tsirm.c
[0]PETSC ERROR: #9 KSPSolve() line 604 in 
/home/couturie/work/petsc-3.6.1/src/ksp/ksp/interface/itfunc.c
[0]PETSC ERROR: #10 main() line 200 in 
/home/couturie/work/petsc-3.6.1/src/ksp/ksp/examples/tutorials/ex15.c

Do you have an idea of the problem?
Is it a good idea to use the LSQR or CGLS method of Petsc?


Moreover I have another question. Do you think that the Petsc community 
would be interested to use this solver. If yes, I would be glad to 
improve it in order to integrate it in the source.

Thanks

Raphaël
-------------- next part --------------
A non-text attachment was scrubbed...
Name: tsirm.c
Type: text/x-csrc
Size: 9217 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20150903/57948fe2/attachment.bin>


More information about the petsc-dev mailing list