Actual source code: ex60.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2: static char help[] = "Tests MatGetColumnVector().";

  4:  #include <petscmat.h>

  6: int main(int argc,char **args)
  7: {
  8:   Mat            C;
  9:   PetscInt       i,j,m = 3,n = 2,Ii,J,col = 0;
 10:   PetscMPIInt    size,rank;
 12:   PetscScalar    v;
 13:   Vec            yy;

 15:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 16:   PetscOptionsGetInt(NULL,NULL,"-col",&col,NULL);

 18:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 19:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 20:   n    = 2*size;

 22:   /* create the matrix for the five point stencil, YET AGAIN*/
 23:   MatCreate(PETSC_COMM_WORLD,&C);
 24:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
 25:   MatSetFromOptions(C);
 26:   MatSeqAIJSetPreallocation(C,5,NULL);
 27:   MatMPIAIJSetPreallocation(C,5,NULL,5,NULL);
 28:   MatSetUp(C);

 30:   for (i=0; i<m; i++) {
 31:     for (j=2*rank; j<2*rank+2; j++) {
 32:       v = -1.0;  Ii = j + n*i;
 33:       if (i>0)   {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 34:       if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 35:       if (j>0)   {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 36:       if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 37:       v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
 38:     }
 39:   }
 40:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 41:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
 42:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);

 44:   VecCreate(PETSC_COMM_WORLD,&yy);
 45:   VecSetSizes(yy,PETSC_DECIDE,m*n);
 46:   VecSetFromOptions(yy);

 48:   MatGetColumnVector(C,yy,col);

 50:   VecView(yy,PETSC_VIEWER_STDOUT_WORLD);

 52:   VecDestroy(&yy);
 53:   MatDestroy(&C);
 54:   PetscFinalize();
 55:   return ierr;
 56: }


 59: /*TEST

 61:    test:
 62:       nsize: 3
 63:       args: -col 7

 65:    test:
 66:       suffix: dense
 67:       nsize: 3
 68:       args: -col 7 -mat_type dense

 70: TEST*/