[petsc-users] Inertia of Hermitian Matrix

Jose E. Roman jroman at dsic.upv.es
Fri Feb 23 05:02:41 CST 2018


Unfortunately MUMPS does not return inertia with complex matrices, it seems that it is not implemented. See the note "Usage with Complex Scalars" in section 3.4.5 of SLEPc users manual.

You could use multi-communicators with as many partitions as MPI processes, so that each process performs a factorization sequentially with PETSc's Cholesky (which returns inertia for complex Hermitian matrices). But this is useful only for matrices that are not too large.

Jose


> El 23 feb 2018, a las 5:39, Anthony Ruth <Anthony.J.Ruth.12 at nd.edu> escribió:
> 
> Hello,
> 
> I am trying to diagonalize a hermitian matrix using the Eigen Problem Solver in SLEPc, I run into errors on calls to MatGetInertia() with complex hermitian matrices that I did not see with real matrices. The complex and real versions were done with separate PETSC_ARCH. I do not know if the problem is with the set up of the matrix or more generally a problem calculating the inertia for a complex matrix. 
> The matrix is created by:
> 
> 	ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr);
> 	ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);CHKERRQ(ierr);
>         ierr = MatSetFromOptions(A);CHKERRQ(ierr);    
> 	ierr = MatSetUp(A);CHKERRQ(ierr);
> 	ierr = MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr);	
> 	ierr = MatGetOwnershipRange(A,&first_row,&last_row);CHKERRQ(ierr);
> 	ierr = MatSetValues(A,m,idxm,n,idxn,data,INSERT_VALUES);CHKERRQ(ierr);
> 	ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
> 	ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
> 
> For a hermitian matrix, all the eigenvalues are real, so I believe it is possible to calculate an inertia by looking at the signs of the diagonal entries. I believe if it was complex but not hermitian, the complex eigenvalues calculating inertia would be difficult.  Is there some problem with doing this through sparse iterative methods? Is there a certain place the matrix needs to be specified as hermitian besides upon assembly? 
> 
> Here is the error stack I see when running:
> 
> 
> Mat Object: 1 MPI processes
>   type: seqaij
> row 0: (0, 0.)  (1, 1. + 1. i) (2, 0.)  (3, 0.)  (4, 0.) 
> row 1: (0, 1. - 1. i) (1, 0.)  (2, 1. + 1. i) (3, 0.)  (4, 0.) 
> row 2: (0, 0.)  (1, 1. - 1. i) (2, 0.)  (3, 1. + 1. i) (4, 0.) 
> row 3: (0, 0.)  (1, 0.)  (2, 1. - 1. i) (3, 0.)  (4, 1. + 1. i)
> row 4: (0, 0.)  (1, 0.)  (2, 0.)  (3, 1. - 1. i) (4, 0.) 
> [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> [0]PETSC ERROR: No support for this operation for this object type
> [0]PETSC ERROR: Mat type mumps
> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
> [0]PETSC ERROR: Petsc Release Version 3.8.2, Nov, 09, 2017 
> [0]PETSC ERROR: Configure options --download-metis --download-mumps --download-parmetis --download-scalapack --with-scalar-type=complex
> [0]PETSC ERROR: #1 MatGetInertia() line 8416 in /home/anthony/DFTB+SIPs/petsc-3.8.2/src/mat/interface/matrix.c
> [0]PETSC ERROR: #2 EPSSliceGetInertia() line 333 in /home/anthony/DFTB+SIPs/slepc-3.8.1/src/eps/impls/krylov/krylovschur/ks-slice.c
> [0]PETSC ERROR: #3 EPSSetUp_KrylovSchur_Slice() line 459 in /home/anthony/DFTB+SIPs/slepc-3.8.1/src/eps/impls/krylov/krylovschur/ks-slice.c
> [0]PETSC ERROR: #4 EPSSetUp_KrylovSchur() line 146 in /home/anthony/DFTB+SIPs/slepc-3.8.1/src/eps/impls/krylov/krylovschur/krylovschur.c
> [0]PETSC ERROR: #5 EPSSetUp() line 165 in /home/anthony/DFTB+SIPs/slepc-3.8.1/src/eps/interface/epssetup.c
> [0]PETSC ERROR: #6 EPSSliceGetEPS() line 298 in /home/anthony/DFTB+SIPs/slepc-3.8.1/src/eps/impls/krylov/krylovschur/ks-slice.c
> [0]PETSC ERROR: #7 EPSSetUp_KrylovSchur_Slice() line 408 in /home/anthony/DFTB+SIPs/slepc-3.8.1/src/eps/impls/krylov/krylovschur/ks-slice.c
> [0]PETSC ERROR: #8 EPSSetUp_KrylovSchur() line 146 in /home/anthony/DFTB+SIPs/slepc-3.8.1/src/eps/impls/krylov/krylovschur/krylovschur.c
> [0]PETSC ERROR: #9 EPSSetUp() line 165 in /home/anthony/DFTB+SIPs/slepc-3.8.1/src/eps/interface/epssetup.c
> [0]PETSC ERROR: #10 SIPSolve() line 195 in /home/anthony/DFTB+SIPs/dftb-eig15/sips.c
> [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> [0]PETSC ERROR: Invalid argument
> [0]PETSC ERROR: Wrong type of object: Parameter # 1
> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
> [0]PETSC ERROR: Petsc Release Version 3.8.2, Nov, 09, 2017 
> [0]PETSC ERROR: Configure options --download-metis --download-mumps --download-parmetis --download-scalapack --with-scalar-type=complex
> [0]PETSC ERROR: #11 EPSGetConverged() line 257 in /home/anthony/DFTB+SIPs/slepc-3.8.1/src/eps/interface/epssolve.c
> [0]PETSC ERROR: #12 squareFromEPS() line 131 in /home/anthony/DFTB+SIPs/dftb-eig15/sips_square.c
> 
> 
> 
> regards,
> Anthony Ruth
> Condensed Matter Theory
> University of Notre Dame



More information about the petsc-users mailing list