[petsc-users] about MATSOR

Sun, Hui hus003 at ucsd.edu
Wed Feb 25 15:02:03 CST 2015


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? 

Another question, there are pc_sor_its and pc_sor_lits in MatSOR. What is local its? 

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
>


More information about the petsc-users mailing list