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>

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

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

 21:   /* create the matrix for the five point stencil, YET AGAIN*/
 22:   MatCreate(PETSC_COMM_WORLD,&C);
 23:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
 24:   MatSetFromOptions(C);
 25:   MatMPIAIJSetPreallocation(C,5,NULL,5,NULL);
 26:   MatSeqAIJSetPreallocation(C,5,NULL);
 27:   for (i=0; i<m; i++) {
 28:     for (j=2*rank; j<2*rank+2; j++) {
 29:       v = -1.0;  Ii = j + n*i;
 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:   }
 37:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 38:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
 39:   PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);
 40:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);
 41:   PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
 42:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);

 44:   MatGetOwnershipRange(C,&rstart,&rend);
 45:   for (i=rstart; i<rend; i++) {
 46:     MatGetRow(C,i,&nz,&idx,&values);
 47:     PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"[%d] get row %" PetscInt_FMT ": ",rank,i);
 48:     for (j=0; j<nz; j++) {
 49:       PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"%" PetscInt_FMT " %g  ",idx[j],(double)PetscRealPart(values[j]));
 50:     }
 51:     PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"\n");
 52:     MatRestoreRow(C,i,&nz,&idx,&values);
 53:   }
 54:   PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout);

 56:   MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A);
 57:   PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);
 58:   MatView(A,PETSC_VIEWER_STDOUT_WORLD);
 59:   PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
 60:   MatView(A,PETSC_VIEWER_STDOUT_WORLD);

 62:   MatDestroy(&A);
 63:   MatDestroy(&C);
 64:   PetscFinalize();
 65:   return 0;
 66: }

 68: /*TEST

 70:    test:
 71:       args: -mat_type seqaij

 73: TEST*/