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.\n\n";

 6:  #include petscmat.h

 10: int main(int argc,char **args)
 11: {
 12:   Mat         C,A;
 13:   int         i,j,m = 3,n = 2,rank,size,I,J,ierr,rstart,rend,nz,*idx;
 14:   PetscScalar v,*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,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&C);
 23:   MatSetFromOptions(C);
 24:   MatMPIAIJSetPreallocation(C,5,PETSC_NULL,5,PETSC_NULL);
 25:   MatSeqAIJSetPreallocation(C,5,PETSC_NULL);
 26:   for (i=0; i<m; i++) {
 27:     for (j=2*rank; j<2*rank+2; j++) {
 28:       v = -1.0;  I = j + n*i;
 29:       if (i>0)   {J = I - n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 30:       if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 31:       if (j>0)   {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 32:       if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 33:       v = 4.0; MatSetValues(C,1,&I,1,&I,&v,INSERT_VALUES);
 34:     }
 35:   }
 36:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 37:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
 38:   PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);
 39:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);
 40:   PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
 41:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);

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

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

 65:   MatDestroy(A);
 66:   MatDestroy(C);
 67:   PetscFinalize();
 68:   return 0;
 69: }