xref: /petsc/src/mat/tests/ex60.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
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   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-col",&col,NULL));
17 
18   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
19   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
20   n    = 2*size;
21 
22   /* create the matrix for the five point stencil, YET AGAIN*/
23   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
24   CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
25   CHKERRQ(MatSetFromOptions(C));
26   CHKERRQ(MatSeqAIJSetPreallocation(C,5,NULL));
27   CHKERRQ(MatMPIAIJSetPreallocation(C,5,NULL,5,NULL));
28   CHKERRQ(MatSetUp(C));
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; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
34       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
35       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
36       if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
37       v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
38     }
39   }
40   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
41   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
42   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
43 
44   CHKERRQ(MatCreateVecs(C,NULL,&yy));
45   CHKERRQ(VecSetFromOptions(yy));
46 
47   CHKERRQ(MatGetColumnVector(C,yy,col));
48 
49   CHKERRQ(VecView(yy,PETSC_VIEWER_STDOUT_WORLD));
50 
51   CHKERRQ(VecDestroy(&yy));
52   CHKERRQ(MatDestroy(&C));
53   ierr = PetscFinalize();
54   return ierr;
55 }
56 
57 /*TEST
58 
59    test:
60       nsize: 3
61       args: -col 7
62 
63    test:
64       suffix: dense
65       nsize: 3
66       args: -col 7 -mat_type dense -vec_type {{mpi standard}}
67       filter: grep -v type
68 
69    test:
70       requires: cuda
71       suffix: dense_cuda
72       nsize: 3
73       output_file: output/ex60_dense.out
74       args: -col 7 -mat_type {{mpidense mpidensecuda}} -vec_type {{mpi standard cuda mpicuda}}
75       filter: grep -v type
76 
77 TEST*/
78