Actual source code: ex9.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  1: 
  2: static char help[] = "Tests MPI parallel matrix creation. Test MatCreateRedundantMatrix() \n\n";

  4: #include <petscmat.h>

  8: int main(int argc,char **args)
  9: {
 10:   Mat            C,Credundant;
 11:   MatInfo        info;
 12:   PetscMPIInt    rank,size,subsize;
 13:   PetscInt       i,j,m = 3,n = 2,low,high,iglobal;
 14:   PetscInt       Ii,J,ldim,nsubcomms;
 16:   PetscBool      flg_info,flg_mat;
 17:   PetscScalar    v,one = 1.0;
 18:   Vec            x,y;
 19:   MPI_Comm       subcomm;

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

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

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

 44:   /* Add extra elements (to illustrate variants of MatGetInfo) */
 45:   Ii   = n; J = n-2; v = 100.0;
 46:   MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);
 47:   Ii   = n-2; J = n; v = 100.0;
 48:   MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);

 50:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 51:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 53:   /* Form vectors */
 54:   MatCreateVecs(C,&x,&y);
 55:   VecGetLocalSize(x,&ldim);
 56:   VecGetOwnershipRange(x,&low,&high);
 57:   for (i=0; i<ldim; i++) {
 58:     iglobal = i + low;
 59:     v       = one*((PetscReal)i) + 100.0*rank;
 60:     VecSetValues(x,1,&iglobal,&v,INSERT_VALUES);
 61:   }
 62:   VecAssemblyBegin(x);
 63:   VecAssemblyEnd(x);

 65:   MatMult(C,x,y);

 67:   PetscOptionsHasName(NULL,NULL,"-view_info",&flg_info);
 68:   if (flg_info)  {
 69:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);
 70:     MatView(C,PETSC_VIEWER_STDOUT_WORLD);
 71: 
 72:     MatGetInfo(C,MAT_GLOBAL_SUM,&info);
 73:     PetscPrintf(PETSC_COMM_WORLD,"matrix information (global sums):\n\
 74:      nonzeros = %D, allocated nonzeros = %D\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated);
 75:     MatGetInfo (C,MAT_GLOBAL_MAX,&info);
 76:     PetscPrintf(PETSC_COMM_WORLD,"matrix information (global max):\n\
 77:      nonzeros = %D, allocated nonzeros = %D\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated);
 78:   }
 79: 
 80:   PetscOptionsHasName(NULL,NULL,"-view_mat",&flg_mat);
 81:   if (flg_mat) {
 82:     MatView(C,PETSC_VIEWER_STDOUT_WORLD);
 83:   }

 85:   /* Test MatCreateRedundantMatrix() */
 86:   nsubcomms = size;
 87:   PetscOptionsGetInt(NULL,NULL,"-nsubcomms",&nsubcomms,NULL);
 88:   MatCreateRedundantMatrix(C,nsubcomms,MPI_COMM_NULL,MAT_INITIAL_MATRIX,&Credundant);
 89:   MatCreateRedundantMatrix(C,nsubcomms,MPI_COMM_NULL,MAT_REUSE_MATRIX,&Credundant);

 91:   PetscObjectGetComm((PetscObject)Credundant,&subcomm);
 92:   MPI_Comm_size(subcomm,&subsize);
 93: 
 94:   if (subsize==2 && flg_mat) {
 95:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(subcomm),"\n[%d] Credundant:\n",rank);
 96:     MatView(Credundant,PETSC_VIEWER_STDOUT_(subcomm));
 97:   }
 98:   MatDestroy(&Credundant);
 99: 
100:   /* Test MatCreateRedundantMatrix() with user-provided subcomm */
101:   {
102:     PetscSubcomm psubcomm;

104:     PetscSubcommCreate(PETSC_COMM_WORLD,&psubcomm);
105:     PetscSubcommSetNumber(psubcomm,nsubcomms);
106:     PetscSubcommSetType(psubcomm,PETSC_SUBCOMM_CONTIGUOUS);
107:     /* enable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */
108:     PetscSubcommSetFromOptions(psubcomm);

110:     MatCreateRedundantMatrix(C,nsubcomms,PetscSubcommChild(psubcomm),MAT_INITIAL_MATRIX,&Credundant);
111:     MatCreateRedundantMatrix(C,nsubcomms,PetscSubcommChild(psubcomm),MAT_REUSE_MATRIX,&Credundant);

113:     PetscSubcommDestroy(&psubcomm);
114:     MatDestroy(&Credundant);
115:   }

117:   VecDestroy(&x);
118:   VecDestroy(&y);
119:   MatDestroy(&C);
120:   PetscFinalize();
121:   return 0;
122: }