[petsc-users] lsqr usage

Ildar Rakhmanov ildar.kk at gmail.com
Wed May 23 17:46:14 CDT 2012


Great, thank you.

Adding this

ierr = KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
  ierr = KSPSetType(ksp,KSPLSQR);//CHKERRQ(ierr);
  ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
  ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);

solved the problem

On Wed, May 23, 2012 at 6:29 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:

>
>  You are trying the default preconditioner ILU(0) which requires a square
> matrix.  Run with -pc_type none    or see the manual page
>
> KSPLSQR - This implements LSQR
>
>   Options Database Keys:
> +   -ksp_lsqr_set_standard_error  - Set Standard Error Estimates of
> Solution see KSPLSQRSetStandardErrorVec()
> .   -ksp_lsqr_monitor - Monitor residual norm and norm of residual of
> normal equations
> -   see KSPSolve()
>
>   Level: beginner
>
>   Notes:
>     This varient, when applied with no preconditioning is identical to the
> original algorithm in exact arithematic; however, in practice, with no
> preconditioning
>     due to inexact arithematic, it can converge differently. Hence when no
> preconditioner is used (PCType PCNONE) it automatically reverts to the
> original algorithm.
>
>     With the PETSc built-in preconditioners, such as ICC, one should call
> KSPSetOperators(ksp,A,A'*A,...) since the preconditioner needs to work
>     for the normal equations A'*A.
>
>     Supports only left preconditioning.
>
>   References:The original unpreconditioned algorithm can be found in Paige
> and Saunders, ACM Transactions on Mathematical Software, Vol 8, pp 43-71,
> 1982.
>     In exact arithmetic the LSQR method (with no preconditioning) is
> identical to the KSPCG algorithm applied to the normal equations.
>     The preconditioned varient was implemented by Bas van't Hof and is
> essentially a left preconditioning for the Normal Equations. It appears the
> implementation with preconditioner
>     track the true norm of the residual and uses that in the convergence
> test.
>
>   Developer Notes: How is this related to the KSPCGNE implementation? One
> difference is that KSPCGNE applies
>            the preconditioner transpose times the preconditioner,  so one
> does not need to pass A'*A as the third argument to KSPSetOperators().
>
>
>   For least squares problems without a zero to A*x = b, there are
> additional convergence tests for the residual of the normal equations,
> A'*(b - Ax), see KSPLSQRDefaultConverged()
>
> .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available
> types), KSP, KSPLSQRDefaultConverged()
>
>
>
>   Barry
>
> On May 23, 2012, at 4:58 PM, Ildar Rakhmanov wrote:
>
> > Here is full message:
> >
> > Matrix A after construction
> > row 0: (0, 1)  (1, 0)
> > row 1: (0, 0)  (1, 1)
> > row 2: (0, 1)  (1, 1)
> > u
> > 1
> > 1
> > b
> > 0
> > 0
> > 0
> > b
> > 1
> > 1
> > 2
> > [0]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> > [0]PETSC ERROR: Invalid argument!
> > [0]PETSC ERROR: Must be square matrix, rows 3 columns 2!
> > [0]PETSC ERROR:
> ------------------------------------------------------------------------
> > [0]PETSC ERROR: Petsc Release Version 3.1.0, Patch 8, Thu Mar 17
> 13:37:48 CDT 2011
> > [0]PETSC ERROR: See docs/changes/index.html for recent updates.
> > [0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> > [0]PETSC ERROR: See docs/index.html for manual pages.
> > [0]PETSC ERROR:
> ------------------------------------------------------------------------
> > [0]PETSC ERROR: ./lsqr on a linux-gnu named ildar by ildar Wed May 23
> 17:13:44 2012
> > [0]PETSC ERROR: Libraries linked from /home/ildar/pbsm/lib
> > [0]PETSC ERROR: Configure run at Sat Mar  3 12:01:22 2012
> > [0]PETSC ERROR: Configure options --prefix=/home/ildar/pbsm
> --download-hypre=1 --download-spooles=1 --download-plapack=1
> --download-spai=1 --download-blacs=1 --download-triangle=1
> --download-f-blas-lapack=1 --download-umfpack=1 --download-sowing=1
> --download-c2html=1 --download-superlu_dist=1 --download-parmetis=1
> --download-scalapack=1 --download-superlu=1 --download-mumps=1
> > [0]PETSC ERROR:
> ------------------------------------------------------------------------
> > [0]PETSC ERROR: MatGetOrdering() line 223 in src/mat/order/sorder.c
> > [0]PETSC ERROR: PCSetUp_ILU() line 194 in
> src/ksp/pc/impls/factor/ilu/ilu.c
> > [0]PETSC ERROR: PCSetUp() line 795 in src/ksp/pc/interface/precon.c
> > [0]PETSC ERROR: KSPSetUp() line 237 in src/ksp/ksp/interface/itfunc.c
> > [0]PETSC ERROR: KSPSolve() line 353 in src/ksp/ksp/interface/itfunc.c
> > [0]PETSC ERROR: main() line 109 in "unknowndirectory/"lsqr.c
> > application called MPI_Abort(MPI_COMM_WORLD, 62) - process 0
> >
> >
> >
> >
> > On Wed, May 23, 2012 at 5:29 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:
> >
> >  You cut and threw away the most important part of the error message.
> Where is the stack trace? Are we to guess where the error message was
> triggered? Please send the entire error message.
> >
> >   Barry
> >
> > On May 23, 2012, at 4:15 PM, Ildar Rakhmanov wrote:
> >
> > > Hi,
> > > I seems I am missing something. I want to use KSPLSQR.
> > > I set up lsqr.c (see attachment)
> > > it uses 3x2 matrix.
> > >
> > > The results is
> > >
> > > ildar at ildar:~/lsqr> ./lsqr
> > > Matrix A after construction
> > > row 0: (0, 1)  (1, 0)
> > > row 1: (0, 0)  (1, 1)
> > > row 2: (0, 1)  (1, 1)
> > > u
> > > 1
> > > 1
> > > b
> > > 0
> > > 0
> > > 0
> > > b
> > > 1
> > > 1
> > > 2
> > > [0]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> > > [0]PETSC ERROR: Invalid argument!
> > > [0]PETSC ERROR: Must be square matrix, rows 3 columns 2!
> > > [0]PETSC ERROR:
> ------------------------------------------------------------------------
> > >
> > >
> > >
> > > Could any one point out how to correctly use LSQR solver?
> > >
> > > <lsqr.c>
> >
> >
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20120523/2be144b4/attachment-0001.html>


More information about the petsc-users mailing list