xref: /petsc/src/mat/tests/ex60.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
1 
2 static char help[] = "Tests MatGetColumnVector().";
3 
4 #include <petscmat.h>
5 
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;
11   PetscErrorCode ierr;
12   PetscScalar    v;
13   Vec            yy;
14 
15   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
16   ierr = PetscOptionsGetInt(NULL,NULL,"-col",&col,NULL);CHKERRQ(ierr);
17 
18   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
19   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
20   n    = 2*size;
21 
22   /* create the matrix for the five point stencil, YET AGAIN*/
23   ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
24   ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr);
25   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
26   ierr = MatSeqAIJSetPreallocation(C,5,NULL);CHKERRQ(ierr);
27   ierr = MatMPIAIJSetPreallocation(C,5,NULL,5,NULL);CHKERRQ(ierr);
28   ierr = MatSetUp(C);CHKERRQ(ierr);
29 
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; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
34       if (i<m-1) {J = Ii + n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
35       if (j>0)   {J = Ii - 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
36       if (j<n-1) {J = Ii + 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
37       v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
38     }
39   }
40   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
41   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
42   ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
43 
44   ierr = VecCreate(PETSC_COMM_WORLD,&yy);CHKERRQ(ierr);
45   ierr = VecSetSizes(yy,PETSC_DECIDE,m*n);CHKERRQ(ierr);
46   ierr = VecSetFromOptions(yy);CHKERRQ(ierr);
47 
48   ierr = MatGetColumnVector(C,yy,col);CHKERRQ(ierr);
49 
50   ierr = VecView(yy,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
51 
52   ierr = VecDestroy(&yy);CHKERRQ(ierr);
53   ierr = MatDestroy(&C);CHKERRQ(ierr);
54   ierr = PetscFinalize();
55   return ierr;
56 }
57 
58 
59 /*TEST
60 
61    test:
62       nsize: 3
63       args: -col 7
64 
65    test:
66       suffix: dense
67       nsize: 3
68       args: -col 7 -mat_type dense
69 
70 TEST*/
71