about MatAssemblyBegin and MatAssemblyEnd
Yujie
recrusader at gmail.com
Mon Mar 3 10:59:09 CST 2008
Dear Barry:
Why do you call these two functions in your code? Because I have realized
the same method with yours to calculate MPIDENSE_MPIAIJ. However, I didn't
call MatAssemblyBegin and
MatAssemblyEnd after I call MatTranspose for matrix C.
In application, sometimes I don't know what do the functions I call in PETSc
do? So, I don't know whether I should
call these two functions although you give me these conditions.
Regarding this example, which condition should it belong to? thanks a
lot.
Regards,
Yujie
On 3/3/08, Barry Smith <bsmith at mcs.anl.gov> wrote:
>
>
> MatAssembly does
> 1) moves any values set on one process that belongs on another
> 2) finalizes the sparse matrix data structure (for example takes any
> holes out of the CSR format)
> doesn't need to do this for dense matrices
> 3) sets up whatever VecScatters that are need to do matrix vector
> products
> 4) marks the matrix as assembled and ready to be used
>
>
> Barry
>
>
> On Mar 2, 2008, at 11:47 PM, Yujie wrote:
>
> > Hi, everyone
> >
> > I just check the codes added by Barry for MPIDENSE_MPIAIJ. The codes
> > is as follows. I am wondering what context one may call
> > MatAssemblyBegin and MatAssemblyEnd?
> > From the webpage, I can find some information about two functions:
> > "Begins assembling the matrix. This routine should be called after
> > completing all calls to MatSetValues()."
> > I mean that I don't know how to judge whether I should call two
> > functions. Could you give me some advice? thanks a lot.
> >
> > 4723 */
> > 4724 PetscErrorCode MatMatMultNumeric_MPIDense_MPIAIJ(Mat A,Mat
> > B,Mat C)
> > 4725 {
> > 4726 PetscErrorCode ierr;
> > 4727 Mat_MPIDense *sub_c = (Mat_MPIDense*)C->data;
> > 4728 Mat At,Bt,Ct;
> > 4729
> > 4730 PetscFunctionBegin;
> > 4731 ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&At);CHKERRQ(ierr);
> > 4732 ierr = MatTranspose(B,MAT_INITIAL_MATRIX,&Bt);CHKERRQ(ierr);
> > 4733 ierr = MatMatMult(Bt,At,MAT_INITIAL_MATRIX,
> > 1.0,&Ct);CHKERRQ(ierr);
> > 4734 ierr = MatDestroy(At);CHKERRQ(ierr);
> > 4735 ierr = MatDestroy(Bt);CHKERRQ(ierr);
> > 4736 ierr = MatTranspose(Ct,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
> > 4737 ierr = MatDestroy(Ct);CHKERRQ(ierr);
> > 4738
> > 4739 C->assembled = PETSC_TRUE;
> > 4740 ierr = MatAssemblyBegin(sub_c-
> > >A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
> > 4741 ierr = MatAssemblyEnd(sub_c->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
> > 4742 if (!C->was_assembled) {
> > 4743 ierr = MatSetUpMultiply_MPIDense(C);CHKERRQ(ierr);
> > 4744 }
> > 4745 PetscFunctionReturn(0);
> > 4746 }
> >
> >
> >
> > Regards,
> > Yujie
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20080303/3126da01/attachment.htm>
More information about the petsc-users
mailing list