[petsc-users] PETSC ERROR: KSPComputeEigenvalues()

Barry Smith bsmith at mcs.anl.gov
Thu Apr 24 06:32:46 CDT 2014


  You must allocate r and c with a length at least numberEigen before call the routine.

   Barry

BTW: I wouldn’t use a numberEigen of 1; even if you really care about only one eigenvalue use a numberEigen of say 10. 

BTW 2: for accurate computations of eigenvalues you should use the package SLEPc which has many high quality eigenvalue routines. The ones in PETSc are only for crude estimates to compute the conditioning of operators and preconditioned operators.

On Apr 24, 2014, at 5:26 AM, Oo <wumeng07maths at qq.com> wrote:

> 
> Hi,
> 
> I met a problem when I tried to use "KSPComputeEigenvalues()".
>  
> The following is the solving part of my code.
> When this programme run to the line 
> " KSPComputeEigenvalues(ksp, numberEigen, r, c,neig);"
> Then, PETSc Errors appear.
> 
> Do you have any suggestions?
> 
> ////////////////////////////////////////////
>    //Solving the Linear System Part======================================================
>     cout<<"solve the Linear system"<<endl;
>     t0=(double)clock();
>     KSP ksp;
>     PC pc;
>     KSPCreate(PETSC_COMM_WORLD, &ksp);
>    // KSPSetType(ksp, KSPCG);
>     KSPSetType(ksp, KSPGMRES);
>     KSPSetOperators(ksp, coeff_matrix, coeff_matrix, DIFFERENT_NONZERO_PATTERN);
> 
>     PetscReal rtol=(1.e-4)*pow(h,4.0);
>     PetscReal stol=(1.e-3)*pow(h,4.0);
> 
>     KSPSetTolerances(ksp, rtol, stol, PETSC_DEFAULT, PETSC_DEFAULT);
> 
>     KSPGetPC(ksp, &pc);
>     PCSetType(pc, PCGAMG);
>     KSPSetFromOptions(ksp);
> 
> 
>     //====KSPSetComputeEigenvalues()
>     KSPSetComputeEigenvalues(ksp, PETSC_TRUE);
>     //KSPSetComputeEigenvalues(ksp, PETSC_FALSE);
>     //=============================
> 
> 
>     KSPSetUp(ksp);
> 
>     KSPSolve(ksp, load_vect, m_solut);
> 
> 
>      cout<<"Before KSPComputeEignvalues()"<<endl;
>     //======KSPComputeEigenvalues()
>     numberEigen=1;//4*P_FEMMesh->P_Vertices->size();
>     KSPComputeEigenvalues(ksp,numberEigen,r,c,neig);
>     cout<<"After KSPComputeEignvalues()"<<endl;
> 
>     cout<<"EignValues:"<<endl;
>     for(unsigned int j=0; j<numberEigen; j++)
>     {
>         cout<<r[numberEigen]<<" + "<<c[numberEigen]<<"i"<<endl;
>     }
>     //========================
> 
> 
> 
> Here is Petsc errors.
> when the programme runs to 
> "KSPComputeEigenvalues(ksp,numberEigen,r,c,neig);"
> These errors appear.
> 
> PETSC ERRORs
> [0]PETSC ERROR: --------------------- Error Message ------------------------------------
> [0]PETSC ERROR: Null argument, when expecting valid pointer!
> [0]PETSC ERROR: Null Pointer: Parameter # 3!
> [0]PETSC ERROR: ------------------------------------------------------------------------
> [0]PETSC ERROR: Petsc Release Version 3.4.3, Oct, 15, 2013 
> [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: /Users/wumeng/MyWork/BuildMSplineTools/bin/msplinePDE_PFEM_2 on a arch-darwin-c-debug named vis032b.sophia.inria.fr by wumeng Thu Apr 24 11:50:16 2014
> [0]PETSC ERROR: Libraries linked from /Users/wumeng/MyWork/PETSc/arch-darwin-c-debug/lib
> [0]PETSC ERROR: Configure run at Wed Dec 18 15:31:31 2013
> [0]PETSC ERROR: Configure options --with-cc=gcc --with-fc=gfortran --download-f-blas-lapack --download-mpich
> [0]PETSC ERROR: ------------------------------------------------------------------------
> [0]PETSC ERROR: KSPComputeEigenvalues() line 118 in /Users/wumeng/MyWork/PETSc/src/ksp/ksp/interface/itfunc.c
> After KSPComputeEignvalues()
> EignValues:
> [0]PETSC ERROR: ------------------------------------------------------------------------
> [0]PETSC ERROR: Caught signal number 11 SEGV: Segmentation Violation, probably memory access out of range
> [0]PETSC ERROR: Try option -start_in_debugger or -on_error_attach_debugger
> [0]PETSC ERROR: or see http://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind[0]PETSC ERROR: or try http://valgrind.org on GNU/linux and Apple Mac OS X to find memory corruption errors
> [0]PETSC ERROR: likely location of problem given in stack below
> [0]PETSC ERROR: ---------------------  Stack Frames ------------------------------------
> [0]PETSC ERROR: Note: The EXACT line numbers in the stack are not available,
> [0]PETSC ERROR:       INSTEAD the line number of the start of the function
> [0]PETSC ERROR:       is given.
> [0]PETSC ERROR: [0] KSPComputeEigenvalues line 116 /Users/wumeng/MyWork/PETSc/src/ksp/ksp/interface/itfunc.c
> [0]PETSC ERROR: --------------------- Error Message ------------------------------------
> [0]PETSC ERROR: Signal received!
> [0]PETSC ERROR: ------------------------------------------------------------------------
> [0]PETSC ERROR: Petsc Release Version 3.4.3, Oct, 15, 2013 
> [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: /Users/wumeng/MyWork/BuildMSplineTools/bin/msplinePDE_PFEM_2 on a arch-darwin-c-debug named vis032b.sophia.inria.fr by wumeng Thu Apr 24 11:50:16 2014
> [0]PETSC ERROR: Libraries linked from /Users/wumeng/MyWork/PETSc/arch-darwin-c-debug/lib
> [0]PETSC ERROR: Configure run at Wed Dec 18 15:31:31 2013
> [0]PETSC ERROR: Configure options --with-cc=gcc --with-fc=gfortran --download-f-blas-lapack --download-mpich
> [0]PETSC ERROR: ------------------------------------------------------------------------
> [0]PETSC ERROR: User provided function() line 0 in unknown directory unknown file
> application called MPI_Abort(MPI_COMM_WORLD, 59) - process 0
> [unset]: aborting job:
> application called MPI_Abort(MPI_COMM_WORLD, 59) - process 0
> 
> 
> 
> Thanks,
> 
> M.



More information about the petsc-users mailing list