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