[petsc-users] Get LU decomposition of a rectangular matrix

Natacha BEREUX natacha.bereux at gmail.com
Tue Mar 13 07:04:32 CDT 2018


Hi Barry
Thanks for your answer,
I followed your suggestion (matrix type = MPIAIJ and superlu_dist on a
single processor)   but I still get the same problem.


Rectangular matrices seem to be forbidden in MatGetOrdering(), whatever
package is used for the LU decomposition

Here is the output :
 ./ex12f -f data/constraints.bin -pc_factor_mat_solver_package superlu_dist
-pc_type lu

   9   33    9   33  1.     0.3E+02 0.3E+02  0.     0.7E+04  1.      0.
[0]PETSC ERROR: --------------------- Error Message
--------------------------------------------------------------
[0]PETSC ERROR: Invalid argument
[0]PETSC ERROR: Must be square matrix, rows 33 columns 9
[0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for
trouble shooting.
[0]PETSC ERROR: Petsc Release Version 3.8.3, Dec, 09, 2017
[0]PETSC ERROR: ./ex12f on a
linux-debug-mumps-ml-hypre-superlu-superlu_dist named dsp0780444 by H03755
Tue Mar 13 08:57:22 2018
[0]PETSC ERROR: Configure options --with-mpi=1 --with-debugging=1
--with-mumps-lib="-L/home/H03755/dev/codeaster-prerequisites/v14/prerequisites//Mumps-511_consortium_aster/MPI/lib
-lzmumps -ldmumps -lmumps_common -lpord
-L/home/H03755/dev/codeaster-prerequisites/v14/prerequisites//Parmetis_aster-403_aster/lib
-lparmetis
-L/home/H03755/dev/codeaster-prerequisites/v14/prerequisites//Scotch_aster-604_aster6/MPI/lib
-lptscotch  -lptscotcherr  -lptscotcherrexit  -lptscotchparmetis  -lesmumps
-lscotch -lscotcherr -lscotcherrexit
-L/home/H03755/dev/codeaster-prerequisites/v14/prerequisites//Metis_aster-510_aster1/lib
-lmetis"
--with-mumps-include=/home/H03755/dev/codeaster-prerequisites/v14/prerequisites//Mumps-511_consortium_aster/MPI/include
--download-hypre=/home/H03755/Librairies/hypre-2.12.0.tar.gz
--download-ml=/home/H03755/Librairies/petsc-pkg-ml-e5040d11aa07.tar.gz
--with-openmp=0 --with-scalapack-lib="-lscalapack-openmpi -lblacs-openmpi
-lblacsF77init-openmpi -lblacsCinit-openmpi"
--with-blaslapack-lib="-llapack -lopenblas"
--PETSC_ARCH=linux-debug-mumps-ml-hypre-superlu-superlu_dist LIBS=-lgomp
--with-parmetis-dir=/home/H03755/dev/codeaster-prerequisites/v14/prerequisites//Parmetis_aster-403_aster
--with-metis-dir=/home/H03755/dev/codeaster-prerequisites/v14/prerequisites//Metis_aster-510_aster1
--prefix=/home/H03755/local/petsc/petsc-3.8.3
--download-superlu=/home/H03755/Librairies/superlu_dist-master.zip
--with-parmetis-dir=/home/H03755/dev/codeaster-prerequisites/v14/prerequisites//Parmetis_aster-403_aster
--with-metis-dir=/home/H03755/dev/codeaster-prerequisites/v14/prerequisites//Metis_aster-510_aster1
--download-superlu=/home/H03755/Librairies/superlu-a0819410c9eb779f9b296cdd95fbdfd96986ae10.tar.gz
[0]PETSC ERROR: #1 MatGetOrdering() line 243 in
/home/H03755/Librairies/petsc-3.8.3/src/mat/order/sorder.c
[0]PETSC ERROR: #2 MatGetOrdering() line 192 in
/home/H03755/Librairies/petsc-3.8.3/src/mat/order/sorder.c
[0]PETSC ERROR: #3 PCSetUp_LU() line 84 in
/home/H03755/Librairies/petsc-3.8.3/src/ksp/pc/impls/factor/lu/lu.c
[0]PETSC ERROR: #4 PCSetUp() line 924 in
/home/H03755/Librairies/petsc-3.8.3/src/ksp/pc/interface/precon.c
[0]PETSC ERROR: #5 User provided function() line 0 in User file


My code looks like :
call MatLoad(A,fd,ierr)
call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)
call KSPSetOperators(ksp,A,A,ierr)
call KSPGetPC(ksp,pc,ierr)
call PCSetFromOptions( pc, ierr)
call PCSetUp(pc, ierr)

Should I call directly a factorization routine ?  Instead of trying to
build a preconditionner ?


Best regards,
Natacha


On Mon, Mar 12, 2018 at 3:46 PM, Smith, Barry F. <bsmith at mcs.anl.gov> wrote:

>
>    Create the matrix on one process but make it a MPIAIJ matrix (not
> SeqAIJ or AIJ) then you should be able to do the factorization with
> SuperLU_Dist
>
>     Please let us know how it goes,
>
>     Barry
>
>
> > On Mar 12, 2018, at 4:35 AM, Natacha BEREUX <natacha.bereux at gmail.com>
> wrote:
> >
> > Dear Petsc user's
> >
> > I am trying to use SuperLU_dist to compute LU decomposition of a
> non-square matrix.
> >
> > My purpose is to build a basis of the nullspace of a matrix C . This
> matrix stores constraints (i.e. linear relationships between degrees of
> freedom): it is thus rectangular (size 9 x 33) and not square
> > I need to get the L U factors of C^T which is a tall skinny matrix.
> > It is possible when using directly SuperLu or SuperLU_dist
> > But when I use Petsc interface to SuperLU I get the following error :
> >
> > [0]PETSC ERROR: Invalid argument
> > [0]PETSC ERROR: Must be square matrix, rows 9 columns 33
> > [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html
> for trouble shooting.
> > [0]PETSC ERROR: Petsc Release Version 3.8.3, Dec, 09, 2017
> > [0]PETSC ERROR: #1 MatGetOrdering() line 243 in
> /home/H03755/Librairies/petsc-3.8.3/src/mat/order/sorder.c
> > [0]PETSC ERROR: #2 PCSetUp_LU() line 84 in /home/H03755/Librairies/petsc-
> 3.8.3/src/ksp/pc/impls/factor/lu/lu.c
> > [0]PETSC ERROR: #3 PCSetUp() line 924 in /home/H03755/Librairies/petsc-
> 3.8.3/src/ksp/pc/interface/precon.c
> > [0]PETSC ERROR: #4 KSPSetUp() line 381 in /home/H03755/Librairies/petsc-
> 3.8.3/src/ksp/ksp/interface/itfunc.c
> > [0]PETSC ERROR: #5 KSPSolve() line 612 in /home/H03755/Librairies/petsc-
> 3.8.3/src/ksp/ksp/interface/itfunc.c
> > [0]PETSC ERROR: #6 User provided function() line 0 in User file
> >
> > The use of non-square matrices is forbidden.
> >
> > Is there a way to specify that I want to use a rectangular matrix ? How
> can I proceed to call SuperLU through PETSc with a rectangular matrix ?
> >
> > Thanks a lot for your help
> > Best regards,
> >
> > Natacha
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20180313/fd55a686/attachment-0001.html>


More information about the petsc-users mailing list