[petsc-users] about MatTranspose

Jed Brown jedbrown at mcs.anl.gov
Thu Sep 1 12:54:27 CDT 2011


On Thu, Sep 1, 2011 at 12:45, Likun Tan <likunt at andrew.cmu.edu> wrote:

> I still have some confusions. When computing B before MPI_Scan(), could i
> compute the values in parallel? After using MPI_Scan(), does that mean the
> columns of B will be gathered in one processor?
>

No, the scan is just computing the start column for a given process.


>
> Here are main steps i took based on your suggestion,
>
> PetscSplitOwnership(PETSC_COMM_WORLD, &n, &N);
>

Add these two lines here:

MPI_Scan(&n, &cstart, 1, MPIU_INT, MPI_SUM, PETSC_COMM_WORLD);
cstart -= n;


> MatCreateSeqDense(PETSC_COMM_SELF, M, n, PETSC_NULL, &B);
> MatCreateSeqDense(PETSC_COMM_SELF, M, M, PETSC_NULL, &A);
>

Good


> MatCreateSeqDense(PETSC_COMM_SELF, M, N, PETSC_NULL, &x);
>

Replace with

MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&x);


> for(j=0; j<N; j++)
>

change the loop bounds here:

for(j=0; j<n; j++)


> {
>   for(i=0; i<M; i++)
>   {
>

Good, now compute value as the value that goes in (i,cstart+j).

      MatSetValues(B, 1, &i, 1, &j, &value, INSERT_VALUES);
>   }
> }
> MatAssemblyBegin(...);
> MatAssemblyEnd(...)
>

This part is correct.


> MPI_Scan(&n, &cstart, 1, MPIU_INT, MPI_SUM, PETSC_COMM_WORLD);
> cstart -= n;
>

We did this already.


>
> MatConvert(...);
> MatCholeskyFactor(...);
> MatMatSolve(...);
>

Yes.


You can gather the matrix x onto all processes if you need the whole result
everywhere, but for performance reasons if you scale further, you should
avoid it if possible.
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20110901/ee8df973/attachment.htm>


More information about the petsc-users mailing list