Actual source code: ex21.c

  1: /*$Id: ex21.c,v 1.22 2001/08/07 03:03:07 balay Exp $*/

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

 6:  #include petscmat.h

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

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

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

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

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

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