xref: /petsc/src/mat/tests/ex60.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1 static char help[] = "Tests MatGetColumnVector().";
2 
3 #include <petscmat.h>
4 
main(int argc,char ** args)5 int main(int argc, char **args)
6 {
7   Mat         C;
8   PetscInt    i, j, m = 3, n = 2, Ii, J, col = 0;
9   PetscMPIInt size, rank;
10   PetscScalar v;
11   Vec         yy;
12 
13   PetscFunctionBeginUser;
14   PetscCall(PetscInitialize(&argc, &args, NULL, 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;
32       Ii = j + n * i;
33       if (i > 0) {
34         J = Ii - n;
35         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
36       }
37       if (i < m - 1) {
38         J = Ii + n;
39         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
40       }
41       if (j > 0) {
42         J = Ii - 1;
43         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
44       }
45       if (j < n - 1) {
46         J = Ii + 1;
47         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
48       }
49       v = 4.0;
50       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
51     }
52   }
53   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
54   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
55   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
56 
57   PetscCall(MatCreateVecs(C, NULL, &yy));
58   PetscCall(VecSetFromOptions(yy));
59 
60   PetscCall(MatGetColumnVector(C, yy, col));
61 
62   PetscCall(VecView(yy, PETSC_VIEWER_STDOUT_WORLD));
63 
64   PetscCall(VecDestroy(&yy));
65   PetscCall(MatDestroy(&C));
66   PetscCall(PetscFinalize());
67   return 0;
68 }
69 
70 /*TEST
71 
72    test:
73       nsize: 3
74       args: -col 7
75 
76    test:
77       suffix: dense
78       nsize: 3
79       args: -col 7 -mat_type dense -vec_type {{mpi standard}}
80       filter: grep -v type
81 
82    test:
83       requires: cuda
84       suffix: dense_cuda
85       nsize: 3
86       output_file: output/ex60_dense.out
87       args: -col 7 -mat_type {{mpidense mpidensecuda}} -vec_type {{mpi standard cuda mpicuda}}
88       filter: grep -v type
89 
90 TEST*/
91