[petsc-users] Speed of MatGetRowSum
Derek Gaston
friedmud at gmail.com
Fri Apr 27 13:11:51 CDT 2018
Thanks for your response! I suppose I should have just checked the code
out myself :-)
Makes sense - I won't worry about it then unless it shows up in profiling
(which I doubt). I just wanted to make sure it wasn't an N^2 operation
(which it's not since my matrices are very sparse).
Thanks again,
Derek
On Fri, Apr 27, 2018 at 12:00 PM Smith, Barry F. <bsmith at mcs.anl.gov> wrote:
>
> Here is the current code: ;-)
>
> PetscErrorCode MatGetRowSum(Mat mat, Vec v)
> {
> Vec ones;
> PetscErrorCode ierr;
>
> PetscFunctionBegin;
> PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
> PetscValidType(mat,1);
> PetscValidHeaderSpecific(v,VEC_CLASSID,2);
> if (!mat->assembled)
> SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for
> unassembled matrix");
> MatCheckPreallocated(mat,1);
> ierr = MatCreateVecs(mat,&ones,NULL);CHKERRQ(ierr);
> ierr = VecSet(ones,1.);CHKERRQ(ierr);
> ierr = MatMult(mat,ones,v);CHKERRQ(ierr);
> ierr = VecDestroy(&ones);CHKERRQ(ierr);
> PetscFunctionReturn(0);
> }
>
> An optimized for AIJ one would avoid creating the matrix and avoid the
> multiples by 1 in the summation and avoid the communication.
>
> I wouldn't worry about its performance; it won't be the bottle neck
> generally. Writing a custom routine for AIJ (seq and MPI) would be easy,
> just structure it like the MatMult_xxx() routine.
>
> Barry
>
>
> > On Apr 27, 2018, at 11:58 AM, Derek Gaston <friedmud at gmail.com> wrote:
> >
> > I was wondering about the comment on MatGetRowSum that says:
> >
> > "This code is slow since it is not currently specialized for different
> formats"
> >
> > How slow are we talking about here? Is it actually going to loop over
> i,j from 0 to N on a sparse matrix?
> >
> > I need to get a vector that represents the "lumped" mass matrix: but I
> would like to get it fairly quickly. I would just compute the lumped
> matrix directly and then use MatGetDiag - but I actually want both
> (unlumped/consistent and lumped) and MatGetRowSum seemed like the best
> direction. Willing to take advice though!
> >
> > Derek
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20180427/ebdb69c2/attachment.html>
More information about the petsc-users
mailing list