Actual source code: ex21.c

  2: static char help[] = "Tests converting a parallel AIJ formatted matrix to the parallel Row format.\n\
  3:  This also tests MatGetRow() and MatRestoreRow() for the parallel case.\n\n";

 5:  #include petscmat.h

  9: int main(int argc,char **args)
 10: {
 11:   Mat               C,A;
 12:   PetscInt          i,j,m = 3,n = 2,I,J,rstart,rend,nz;
 13:   PetscMPIInt       rank,size;
 14:   PetscErrorCode    ierr;
 15:   const PetscInt    *idx;
 16:   PetscScalar       v;
 17:   const PetscScalar *values;

 19:   PetscInitialize(&argc,&args,(char *)0,help);
 20:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 21:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 22:   n = 2*size;

 24:   /* create the matrix for the five point stencil, YET AGAIN*/
 25:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&C);
 26:   MatSetFromOptions(C);
 27:   MatMPIAIJSetPreallocation(C,5,PETSC_NULL,5,PETSC_NULL);
 28:   MatSeqAIJSetPreallocation(C,5,PETSC_NULL);
 29:   for (i=0; i<m; i++) {
 30:     for (j=2*rank; j<2*rank+2; j++) {
 31:       v = -1.0;  I = j + n*i;
 32:       if (i>0)   {J = I - n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 33:       if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 34:       if (j>0)   {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 35:       if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 36:       v = 4.0; MatSetValues(C,1,&I,1,&I,&v,INSERT_VALUES);
 37:     }
 38:   }
 39:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 40:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
 41:   PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);
 42:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);
 43:   PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
 44:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);

 46:   MatGetOwnershipRange(C,&rstart,&rend);
 47:   for (i=rstart; i<rend; i++) {
 48:     MatGetRow(C,i,&nz,&idx,&values);
 49:     PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"[%d] get row %D: ",rank,i);
 50:     for (j=0; j<nz; j++) {
 51: #if defined(PETSC_USE_COMPLEX)
 52:       PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"%D %g  ",idx[j],PetscRealPart(values[j]));
 53: #else
 54:       PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"%D %g  ",idx[j],values[j]);
 55: #endif
 56:     }
 57:     PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"\n");
 58:     MatRestoreRow(C,i,&nz,&idx,&values);
 59:   }
 60:   PetscSynchronizedFlush(PETSC_COMM_WORLD);

 62:   MatConvert(C,MATSAME,&A);
 63:   PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);
 64:   MatView(A,PETSC_VIEWER_STDOUT_WORLD);
 65:   PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
 66:   MatView(A,PETSC_VIEWER_STDOUT_WORLD);

 68:   MatDestroy(A);
 69:   MatDestroy(C);
 70:   PetscFinalize();
 71:   return 0;
 72: }