[petsc-users] Slepc: Computing first 4 smallest eigenvalues and eigenvectors of a large graph Laplacian
Barry Smith
bsmith at mcs.anl.gov
Sun Mar 26 12:30:06 CDT 2017
The flop rate for the MatMult and MatSolve are much lower than I would expect. Try using the AIJ matrix instead of SBAIJ.
> On Mar 26, 2017, at 6:25 AM, Bodhisatta Pramanik <bodhi91 at iastate.edu> wrote:
>
> *First set a new PETSC_ARCH say arch-opt and ./configure with the additional argument --with-debugging=0 then recompile PETSc with the new PETSC_ARCH and recompile your code then send the new -log_view output.*
> I have attached the new log summary with this file. The computation time is a lot faster now. It takes around 23 seconds to compute the first 3 eigenvectors of a 200K sized matrix.
>
> Thanks,
> Bodhi
>
> On Sun, Mar 26, 2017 at 3:14 AM, Bodhisatta Pramanik <bodhi91 at iastate.edu> wrote:
> *The error is in EPSComputeError(), not in EPSSolve(). Did you modify the state of the EPS object between EPSSolve() and EPSComputeError()?*
>
> I am not modifying the eps object. I have copied a part of the code where I am computing the eigenvectors.
>
> EPSCreate(PETSC_COMM_WORLD,&eps);
> EPSSetOperators(eps,m.Lpl,NULL);
> EPSSetProblemType(eps,EPS_HEP);
> EPSSetFromOptions(eps);
>
> EPSSolve(eps); //Solving for the eigenvalues
>
> EPSGetIterationNumber(eps,&its);
> PetscPrintf(PETSC_COMM_WORLD,"Number of iterations of the method: %D\n",its);
> EPSGetType(eps,&type);
> PetscPrintf(PETSC_COMM_WORLD,"Solution Method: %s\n\n",type);
> EPSGetTolerances(eps,&get_tol,&maxit);
> PetscPrintf(PETSC_COMM_WORLD,"Stopping condition: tol=%.4g, maxit=%D\n",(double)get_tol,maxit);
> EPSGetConverged(eps,&nconv);
> PetscPrintf(PETSC_COMM_WORLD, "Number of converged eigenpairs: %D\n\n",nconv);
> if(nconv>0)
> {
> PetscPrintf(PETSC_COMM_WORLD,
> " k ||Ax-kx||/||kx||\n"
> "-------------- ----------------\n");
> for(i=0;i<nconv;i++)
> {
> EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
> EPSComputeError(eps,i,EPS_ERROR_RELATIVE,&error);
> re = PetscRealPart(kr);
> im = PetscImaginaryPart(kr);
>
> re = kr;
> im = ki;
> if(im!=0.0) {
> PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi %12g\n",(double)re,(double)im,(double)error);
> }
> else {
> PetscPrintf(PETSC_COMM_WORLD," %12f %12g\n",(double)re,(double)error);
> }
> }
> PetscPrintf(PETSC_COMM_WORLD,"\n");
> }
>
> EPSDestroy(&eps);
> VecDestroy(&xr);
> VecDestroy(&xi);
> }
>
> This is what I pass through my command line:
> ./RunPart -eps_type jd -eps_nev 3 -st_ksp_type cg -st_ksp_rtol 0.001 -eps_tol 0.001 -st_pc_type bjacobi -eps_smallest_real
>
>
> *In graph partitioning, I would strongly recommend deflating the eigenvector corresponding to the zero eigenvalue, as is done in ex11.c:*
> I have tried doing this but for some reason the eigenvalues fail to converge.
>
>
>
> On Sun, Mar 26, 2017 at 2:41 AM, Jose E. Roman <jroman at dsic.upv.es> wrote:
>
> > El 26 mar 2017, a las 6:08, Bodhisatta Pramanik <bodhi91 at iastate.edu> escribió:
> >
> > Send all the output in the error message:
> >
> > [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> > [0]PETSC ERROR: Object is in wrong state
> > [0]PETSC ERROR: Must call EPSSolve() first: 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.7.5, Jan, 01, 2017
> > [0]PETSC ERROR: ./RunPart on a arch-linux2-c-debug named research-5.ece.iastate.edu by bodhi91 Sat Mar 25 22:33:56 2017
> > [0]PETSC ERROR: Configure options --with-cc=gcc --with-css=g++ -with-fc=gfortran --download-fblaslapack --download-mpich
> > [0]PETSC ERROR: #10468 EPSComputeError() line 643 in /tmp/Bodhi/slepc-3.7.3/src/eps/interface/epssolve.c
> > 0.000000 2.0795e-317
> >
>
> The error is in EPSComputeError(), not in EPSSolve(). Did you modify the state of the EPS object between EPSSolve() and EPSComputeError()?
>
> In graph partitioning, I would strongly recommend deflating the eigenvector corresponding to the zero eigenvalue, as is done in ex11.c:
> http://slepc.upv.es/documentation/current/src/eps/examples/tutorials/ex11.c.html
>
> Jose
>
>
>
>
> --
> Bodhisatta Pramanik,
> Graduate Student,
> Department of Electrical and Computer Engineering,
> 301 Durham,
> Iowa State University,
> Ames,Iowa 50011,
> bodhi91 at iastate.edu
> 515-735-6300
>
>
>
> --
> Bodhisatta Pramanik,
> Graduate Student,
> Department of Electrical and Computer Engineering,
> 301 Durham,
> Iowa State University,
> Ames,Iowa 50011,
> bodhi91 at iastate.edu
> 515-735-6300
> <New Log Summary>
More information about the petsc-users
mailing list