Actual source code: ex25.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: static char help[] = "Tests MatTranspose()\n\n";

  4: #include <petscmat.h>

  8: int main(int argc,char **args)
  9: {
 10:   Mat            C,A;
 11:   PetscScalar    v;
 12:   PetscInt       i,j,m = 4,n = 4,Ii,J,Istart,Iend;
 13:   PetscMPIInt    rank,size;
 15:   PetscBool      equal=PETSC_FALSE;

 17:   PetscInitialize(&argc,&args,(char*)0,help);
 18:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 19:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 21:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 22:   n    = m;

 24:   MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,5,NULL,5,NULL,&C);

 26:   /* create the symmetric matrix for the five point stencil */
 27:   MatGetOwnershipRange(C,&Istart,&Iend);
 28:   for (Ii=Istart; Ii<Iend; Ii++) {
 29:     v = -1.0; i = Ii/n; j = Ii - i*n;
 30:     if (i>0)   {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 31:     if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 32:     if (j>0)   {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 33:     if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 34:     v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
 35:   }
 36:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 37:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 39:   MatTranspose(C,MAT_INITIAL_MATRIX,&A);

 41:   MatEqual(C,A,&equal);
 42:   if (!equal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"C != C^T");

 44:   MatDestroy(&C);
 45:   MatDestroy(&A);
 46:   PetscFinalize();
 47:   return 0;
 48: }