[petsc-users] about MATSOR

Matthew Knepley knepley at gmail.com
Wed Feb 25 19:00:42 CST 2015


On Wed, Feb 25, 2015 at 6:36 PM, Sun, Hui <hus003 at ucsd.edu> wrote:

> By the way, how do you add up two matrices? Assume A and B are of same
> size, both sparse and of type Matmpi. I want to perform this operation:
> C=A+B. What should I do? I haven't found anything related in the user
> manual.
>

http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatAXPY.html

   Matt


> Best,
> Hui
>
>
> ________________________________________
> From: Sun, Hui
> Sent: Wednesday, February 25, 2015 2:50 PM
> To: Barry Smith
> Cc: petsc-users at mcs.anl.gov
> Subject: RE: [petsc-users] about MATSOR
>
> Thank you Barry. It is a step in forming my PC matrix for the schur
> complement. I think I will try something else.
>
> Best,
> Hui
>
>
> ________________________________________
> From: Barry Smith [bsmith at mcs.anl.gov]
> Sent: Wednesday, February 25, 2015 1:12 PM
> To: Sun, Hui
> Cc: petsc-users at mcs.anl.gov
> Subject: Re: [petsc-users] about MATSOR
>
> > On Feb 25, 2015, at 3:02 PM, Sun, Hui <hus003 at ucsd.edu> wrote:
> >
> > Thank you, Barry. In fact, my matrix B is sparse. Is it possible to drop
> the terms less than a certain threshold from C, so that C can be sparse?
>
>   You can do that but it becomes much more complicated and expensive,
> especially assembling back the global sparse matrix.   Why do you want to
> do this? Are you doing some algebraic multigrid method with smoothed
> agglomeration or something? If so, this is big projection and it would
> likely be better to build off some existing tool rather than writing
> something from scratch.
>
> >
> > Another question, there are pc_sor_its and pc_sor_lits in MatSOR. What
> is local its?
>
>    In parallel it does its "global iterations" with communication between
> each iteration and inside each global iteration it does lits local
> iterations (without communication between them). (Some people call this
> hybrid parallel SOR.
>
>    Sequentially the number of SOR iterations it does is its*lits since
> there is never any communication between processes since there is only one.
>
> Barry
>
> >
> > Best,
> > Hui
> >
> > ________________________________________
> > From: Barry Smith [bsmith at mcs.anl.gov]
> > Sent: Wednesday, February 25, 2015 11:48 AM
> > To: Sun, Hui
> > Cc: petsc-users at mcs.anl.gov
> > Subject: Re: [petsc-users] about MATSOR
> >
> >   Let me try to understand what you wish to do. You have a sequential
> matrix A and a  matrix B and you wish to compute the (dense) matrix C where
> each column of C is obtained by running three iterations of SOR (using the
> matrix A) with the corresponding column of B as the right hand side for the
> SOR; using an initial guess of zero to start the SOR?
> >
> >   If this is the case then create a seq dense matrix C of the same size
> as B (it will automatically be initialized with zeros).  Say B has n rows
> >
> >   VecCreateSeq(PETSC_COMM_SELF,n,&b);
> >   VecCreateSeqWithArray(PETSC_COMM_SELF,1,n,NULL,&x);
> >   PetscScalar *xx,*carray;
> >   MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&C);
> >   MatDenseGetArray(C,&carray);
> >   loop over columns
> >      ierr = MatGetColumnVector(B,b,col);CHKERRQ(ierr);
> >      ierr = VecPlaceArray(x,carray + n*col);    /* this makes the vec x
> point to the correct column of entries in the matrix C.
> >> err = MatSOR(A, b, 1, (SOR_FORWARD_SWEEP|SOR_ZERO_INITIAL_GUESS), 0, 3,
> 0, x);CHKERRQ(ierr);
> >
> >
> >   Barry
> >
> > Note you could also do this in parallel but it is kind of bogus be the
> SOR is run independently on each process so I don't recommend it as useful.
> >>
> >
> >> On Feb 25, 2015, at 12:14 PM, Sun, Hui <hus003 at ucsd.edu> wrote:
> >>
> >> I want to do 3 steps of gauss seidel from a Mat A to another Mat B. Is
> there a way to do this?
> >>
> >> I mean, what I can think of is to get the column vectors of B by
> >> ierr = MatGetColumnVector(B,v,col);CHKERRQ(ierr);
> >>
> >>
> >>
> >> and apply Gauss seidel from A to v:
> >> ierr = MatSOR(A, v, 1, (SOR_FORWARD_SWEEP|SOR_ZERO_INITIAL_GUESS), 0,
> 3, 0, x);CHKERRQ(ierr);
> >>
> >>
> >>
> >> But then I have to form another matrix C whose columns are composed by
> v, and I'm not sure how to do that.
> >>
> >>
> >>
> >> Best,
> >>
> >> Hui
> >>




-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20150225/9406a005/attachment-0001.html>


More information about the petsc-users mailing list