[petsc-users] Erro in solution

Barry Smith bsmith at mcs.anl.gov
Mon May 9 11:35:25 CDT 2011


   The code is very slow for one of two reasons

1) the matrix preallocation is not correct and it is spending a great deal of time in the MatSetValues() calls. To check this run the code with -info and grep the output for malloc. My guess is that you are not preallocating for the diagonal and hence the preallocation is not enough. See http://www.mcs.anl.gov/petsc/petsc-as/documentation/faq.html#efficient-assembly

2) the convergence of the solver is very slow. Run with -ksp_monitor and see how many iterations it is taking to converge. Is it thousands? If so, you will need to select a better preconditioner or use a direct solver -pc_type lu.

   Barry



On May 9, 2011, at 11:02 AM, Danesh Daroui wrote:

> 
> Hi again,
> 
> The issue is solved and it was due to a bug in my code. But now the
> problem is that PETSc is extremely slow when I try to solve even a small
> equation. I fill the matrix and right-hand-side and call the solver as
> below:
> 
> ierr=KSPCreate(PETSC_COMM_WORLD, &ksp);
> ierr=KSPSetOperators(ksp, Mp, Mp, DIFFERENT_NONZERO_PATTERN);
> ierr=KSPSetTolerances(ksp, 1.e-2/Msize, 1.e-50, PETSC_DEFAULT,
> 			PETSC_DEFAULT);
> ierr=KSPSetFromOptions(ksp);
> ierr=KSPSolve(ksp, bp, xp);
> 
> The Msize is the size of right-hand-side vector. I use C++ and am
> running PETSc in sequential mode so no MPI is used, but anyway, it
> shouldn't be that slow! Am I right or I am missing something?
> 
> Thanks,
> 
> Danesh
> 
> 
> 
> 
> On Sun, 2011-05-08 at 14:41 +0200, dan at ltu.se wrote:
>> Hi,
>> 
>> I checked my code and there is no zero on the diagonal. I also modified the code
>> so if there is any zero on the diagonal it will be stored anyway. After this
>> change the problem still exists! I am not sure if I have called PETSc functions
>> correctly or in correct sequence. But anyway, it returns the error:
>> 
>> [0]PETSC ERROR: Object is in wrong state!
>> 
>> 
>> that and the rest of the error message is exactly as I posted before. Any idea
>> to solve this problem?
>> 
>> Thanks,
>> 
>> D.
>> 
>> 
>> 
>> Quoting Barry Smith <bsmith at mcs.anl.gov>:
>> 
>>> 
>>>   Danesh,
>>> 
>>>    Some of the PETSc solvers (like the default ILU) require that you put
>>> entries on all the diagonal locations of the matrix, even if they are zero.
>>> So you should explicitly put a 0 on the diagonal locations that have 0 and
>>> the issue will go away.
>>> 
>>>    Barry
>>> 
>>> On May 6, 2011, at 1:15 PM, Danesh Daroui wrote:
>>> 
>>>> 
>>>> Dear PETSc users,
>>>> 
>>>> I have recently started to use PETSc in our code. As the first try, I
>>>> assemble the matrix row by row and then put each non-zero element in
>>>> each row into PETSc matrix. Then I start to solve the equation. I use
>>>> PETSc on a multi-core machine, so I will not use MPI in my program. I
>>>> compiled PETSc without MPI. The PETSc parts of my code are:
>>>> 
>>>> 
>>>>       Mat Mp;
>>>>       Vec xp, bp;
>>>>       KSP ksp;
>>>>       PetscErrorCode ierr;
>>>> 
>>>>       ierr=MatCreate(PETSC_COMM_WORLD, &Mp);
>>>>       ierr=MatSetSizes(Mp, Msize, Msize, Msize, Msize);
>>>>       ierr=MatSetType(Mp, MATSEQAIJ);
>>>>       ierr=MatSetUp(Mp);
>>>> 
>>>>       ierr=VecCreate(PETSC_COMM_WORLD, &xp);
>>>>       ierr=VecSetSizes(xp, Msize, Msize);
>>>>       ierr=VecSetFromOptions(xp);
>>>> 
>>>>       ierr=VecCreate(PETSC_COMM_WORLD, &bp);
>>>>       ierr=VecSetSizes(bp, Msize, Msize);
>>>>       ierr=VecSetFromOptions(bp);
>>>> 
>>>> 	.
>>>> 	.
>>>> 	.
>>>> 
>>>> 	// each non zero element of row "i" and column "k" is put in 
>>>> 	// PETSc matrix
>>>> 	if (M_row[k]!=0.0)
>>>> 		ierr=MatSetValues(Mp, 1, &i, 1, &k, &M_row[k],
>>>> 					ADD_VALUES);
>>>> 
>>>> 
>>>> 	ierr=MatAssemblyBegin(Mp, MAT_FINAL_ASSEMBLY);
>>>> 	ierr=MatAssemblyEnd(Mp, MAT_FINAL_ASSEMBLY);
>>>> 
>>>> 	.
>>>> 	.
>>>> 	.
>>>> 
>>>> 	// copy right-hand-side to PETSc vector
>>>> 	for (int i=0; i<Msize; ++i)
>>>> 		VecSetValue(bp, i, b[i], INSERT_VALUES);
>>>> 
>>>> 	VecAssemblyBegin(bp);
>>>> 	VecAssemblyEnd(bp);
>>>> 
>>>> 	ierr=KSPCreate(PETSC_COMM_WORLD, &ksp);
>>>> 	ierr=KSPSetOperators(ksp, Mp, Mp, DIFFERENT_NONZERO_PATTERN);
>>>> 	ierr=KSPSetTolerances(ksp, 1.e-2/Msize, 1.e-50, PETSC_DEFAULT, 
>>>> 				PETSC_DEFAULT);
>>>> 	ierr=KSPSetFromOptions(ksp);
>>>> 	ierr=KSPSolve(ksp, bp, xp);
>>>> 
>>>> 
>>>> 
>>>> 
>>>> 
>>>> The problem is that whenever I run the code I get following output:
>>>> 
>>>> 
>>>> 
>>>> 
>>>> [0]PETSC ERROR: --------------------- Error Message
>>>> ------------------------------------
>>>> [0]PETSC ERROR: Object is in wrong state!
>>>> [0]PETSC ERROR: Matrix is missing diagonal entry 309!
>>>> [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: Unknown Name on a linux-gnu named linux-sclt.site by dan
>>>> Fri May  6 19:34:19 2011
>>>> [0]PETSC ERROR: Libraries linked
>>>> from /home/danesh/petsc-3.1-p8/linux-gnu-intel/lib
>>>> [0]PETSC ERROR: Configure run at Fri May  6 15:50:26 2011
>>>> [0]PETSC ERROR: Configure options PETSC_ARCH=linux-gnu-intel
>>>> --with-scalar-type=complex --with-clanguage=c++ --with-mpi=0
>>>> --with-cc=icc --with-fc=gfortran --with-cxx=icpc
>>>> --with-blas-lapack-dir=/home/danesh/intel/mkl/lib/intel64/
>>>> [0]PETSC ERROR:
>>>> ------------------------------------------------------------------------
>>>> [0]PETSC ERROR: MatILUFactorSymbolic_SeqAIJ_ilu0() line 1627 in
>>>> src/mat/impls/aij/seq/aijfact.c
>>>> [0]PETSC ERROR: MatILUFactorSymbolic_SeqAIJ() line 1731 in
>>>> src/mat/impls/aij/seq/aijfact.c
>>>> [0]PETSC ERROR: MatILUFactorSymbolic() line 5464 in
>>>> src/mat/interface/matrix.c
>>>> [0]PETSC ERROR: PCSetUp_ILU() line 204 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
>>>> 
>>>> 
>>>> 
>>>> I am sure that the matrix is correct and the system is not singular
>>>> because the system can be solved using LAPACK's direct solver. I
>>>> probably format the PETSc matrix incorrectly or some PETSc parameters
>>>> are wrong. I would appreciate if you give me any clue about possible
>>>> problem(s).
>>>> 
>>>> Regards,
>>>> 
>>>> Danesh
>>>> 
>>>> 
>>> 
>>> 
>>> 
> 



More information about the petsc-users mailing list