Actual source code: ex9.c

  2: static char help[] = "Tests MPI parallel matrix creation.\n\n";

 4:  #include petscmat.h

  8: int main(int argc,char **args)
  9: {
 10:   Mat            C;
 11:   MatInfo        info;
 12:   PetscMPIInt    rank,size;
 13:   PetscInt       i,j,m = 3,n = 2,low,high,iglobal;
 14:   PetscInt       Ii,J,ldim;
 16:   PetscTruth     flg;
 17:   PetscScalar    v,one = 1.0;
 18:   Vec            u,b;
 19:   PetscInt       bs,ndiag,diag[7];  bs = 1,ndiag = 5;

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

 26:   MatCreate(PETSC_COMM_WORLD,&C);
 27:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
 28:   MatSetFromOptions(C);

 30:   diag[0] = n;
 31:   diag[1] = 1;
 32:   diag[2] = 0;
 33:   diag[3] = -1;
 34:   diag[4] = -n;
 35:   if (size>1) {ndiag = 7; diag[5] = 2; diag[6] = -2;}

 37:   /* Create the matrix for the five point stencil, YET AGAIN */
 38:   for (i=0; i<m; i++) {
 39:     for (j=2*rank; j<2*rank+2; j++) {
 40:       v = -1.0;  Ii = j + n*i;
 41:       if (i>0)   {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 42:       if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 43:       if (j>0)   {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 44:       if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 45:       v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
 46:     }
 47:   }

 49:   /* Add extra elements (to illustrate variants of MatGetInfo) */
 50:   Ii = n; J = n-2; v = 100.0;
 51:   MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);
 52:   Ii = n-2; J = n; v = 100.0;
 53:   MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);

 55:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 56:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 58:   /* Form vectors */
 59:   VecCreate(PETSC_COMM_WORLD,&u);
 60:   VecSetSizes(u,PETSC_DECIDE,m*n);
 61:   VecSetFromOptions(u);
 62:   VecDuplicate(u,&b);
 63:   VecGetLocalSize(u,&ldim);
 64:   VecGetOwnershipRange(u,&low,&high);
 65:   for (i=0; i<ldim; i++) {
 66:     iglobal = i + low;
 67:     v = one*((PetscReal)i) + 100.0*rank;
 68:     VecSetValues(u,1,&iglobal,&v,INSERT_VALUES);
 69:   }
 70:   VecAssemblyBegin(u);
 71:   VecAssemblyEnd(u);

 73:   MatMult(C,u,b);

 75:   VecView(u,PETSC_VIEWER_STDOUT_WORLD);
 76:   VecView(b,PETSC_VIEWER_STDOUT_WORLD);

 78:   VecDestroy(u);
 79:   VecDestroy(b);

 81:   PetscOptionsHasName(PETSC_NULL,"-view_info",&flg);
 82:   if (flg)  {PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);}
 83:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);

 85:   MatGetInfo(C,MAT_GLOBAL_SUM,&info);
 86:   PetscPrintf(PETSC_COMM_WORLD,"matrix information (global sums):\n\
 87:      nonzeros = %D, allocated nonzeros = %D\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated);
 88:   MatGetInfo (C,MAT_GLOBAL_MAX,&info);
 89:   PetscPrintf(PETSC_COMM_WORLD,"matrix information (global max):\n\
 90:      nonzeros = %D, allocated nonzeros = %D\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated);

 92:   MatDestroy(C);
 93:   PetscFinalize();
 94:   return 0;
 95: }