Actual source code: ex86.c

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1: static char help[] = "Testing MatCreateMPIMatConcatenateSeqMat().\n\n";

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

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

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

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

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

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

 63:   MatDestroy(&seqmat);
 64:   MatDestroy(&mpimat);
 65:   PetscFinalize();
 66:   return ierr;
 67: }