Actual source code: ex86.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  1: static char help[] = "Testing MatCreateMPIMatConcatenateSeqMat().\n\n";

  3: #include <petscmat.h>
  6: int main(int argc,char **argv)
  7: {
  9:   Mat            seqmat,mpimat;
 10:   PetscMPIInt    rank;
 11:   PetscScalar    value[3],*vals;
 12:   PetscInt       i,col[3],n=5,bs=1;
 13: 
 14:   PetscInitialize(&argc,&argv,(char*)0,help);
 15:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 16:   PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);

 18:   /* Create seqaij matrices of size (n+rank) by n */
 19:   MatCreate(PETSC_COMM_SELF,&seqmat);
 20:   MatSetSizes(seqmat,(n+rank)*bs,PETSC_DECIDE,PETSC_DECIDE,n*bs);
 21:   MatSetFromOptions(seqmat);
 22:   MatSeqAIJSetPreallocation(seqmat,3*bs,NULL);
 23:   MatSeqBAIJSetPreallocation(seqmat,bs,3,NULL);

 25:   if (bs == 1) {
 26:     value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
 27:     for (i=1; i<n-1; i++) {
 28:       col[0] = i-1; col[1] = i; col[2] = i+1;
 29:       MatSetValues(seqmat,1,&i,3,col,value,INSERT_VALUES);
 30:     }
 31:     i = n - 1; col[0] = n - 2; col[1] = n - 1;
 32:     MatSetValues(seqmat,1,&i,2,col,value,INSERT_VALUES);

 34:     i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
 35:     MatSetValues(seqmat,1,&i,2,col,value,INSERT_VALUES);
 36:   } else {
 37:     PetscInt *rows,*cols,j;
 38:     PetscMalloc3(bs*bs,&vals,bs,&rows,bs,&cols);
 39:     /* diagonal blocks */
 40:     for (i=0; i<bs*bs; i++) vals[i] = 2.0;
 41:     for (i=0; i<n*bs; i+=bs) {
 42:       for (j=0; j<bs; j++) {rows[j] = i+j; cols[j] = i+j;}
 43:       MatSetValues(seqmat,bs,rows,bs,cols,vals,INSERT_VALUES);
 44:     }
 45:     /* off-diagonal blocks */
 46:     for (i=0; i<bs*bs; i++) vals[i] = -1.0;
 47:     for (i=0; i<(n-1)*bs; i+=bs) {
 48:       for (j=0; j<bs; j++) {rows[j] = i+j; cols[j] = i+bs+j;}
 49:       MatSetValues(seqmat,bs,rows,bs,cols,vals,INSERT_VALUES);
 50:     }

 52:     PetscFree3(vals,rows,cols);
 53:   }
 54:   MatAssemblyBegin(seqmat,MAT_FINAL_ASSEMBLY);
 55:   MatAssemblyEnd(seqmat,MAT_FINAL_ASSEMBLY);
 56:   if (!rank) {
 57:     PetscPrintf(PETSC_COMM_SELF,"[%d] seqmat:\n",rank);
 58:     MatView(seqmat,PETSC_VIEWER_STDOUT_SELF);
 59:   }

 61:   MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,seqmat,PETSC_DECIDE,MAT_INITIAL_MATRIX,&mpimat);
 62:   MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,seqmat,PETSC_DECIDE,MAT_REUSE_MATRIX,&mpimat);
 63:   MatView(mpimat,PETSC_VIEWER_STDOUT_WORLD);

 65:   MatDestroy(&seqmat);
 66:   MatDestroy(&mpimat);
 67:   PetscFinalize();
 68:   return 0;
 69: }