xref: /petsc/src/mat/tests/ex60.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatGetColumnVector().";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7*d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   Mat         C;
9c4762a1bSJed Brown   PetscInt    i, j, m = 3, n = 2, Ii, J, col = 0;
10c4762a1bSJed Brown   PetscMPIInt size, rank;
11c4762a1bSJed Brown   PetscScalar v;
12c4762a1bSJed Brown   Vec         yy;
13c4762a1bSJed Brown 
14327415f7SBarry Smith   PetscFunctionBeginUser;
159566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-col", &col, NULL));
17c4762a1bSJed Brown 
189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
20c4762a1bSJed Brown   n = 2 * size;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN*/
239566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
249566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
259566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
269566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(C, 5, NULL));
279566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(C, 5, NULL, 5, NULL));
289566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   for (i = 0; i < m; i++) {
31c4762a1bSJed Brown     for (j = 2 * rank; j < 2 * rank + 2; j++) {
329371c9d4SSatish Balay       v  = -1.0;
339371c9d4SSatish Balay       Ii = j + n * i;
349371c9d4SSatish Balay       if (i > 0) {
359371c9d4SSatish Balay         J = Ii - n;
369371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
379371c9d4SSatish Balay       }
389371c9d4SSatish Balay       if (i < m - 1) {
399371c9d4SSatish Balay         J = Ii + n;
409371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
419371c9d4SSatish Balay       }
429371c9d4SSatish Balay       if (j > 0) {
439371c9d4SSatish Balay         J = Ii - 1;
449371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
459371c9d4SSatish Balay       }
469371c9d4SSatish Balay       if (j < n - 1) {
479371c9d4SSatish Balay         J = Ii + 1;
489371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
499371c9d4SSatish Balay       }
509371c9d4SSatish Balay       v = 4.0;
519371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
52c4762a1bSJed Brown     }
53c4762a1bSJed Brown   }
549566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
559566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
569566063dSJacob Faibussowitsch   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
57c4762a1bSJed Brown 
589566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(C, NULL, &yy));
599566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(yy));
60c4762a1bSJed Brown 
619566063dSJacob Faibussowitsch   PetscCall(MatGetColumnVector(C, yy, col));
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(VecView(yy, PETSC_VIEWER_STDOUT_WORLD));
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&yy));
669566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
679566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
68b122ec5aSJacob Faibussowitsch   return 0;
69c4762a1bSJed Brown }
70c4762a1bSJed Brown 
71c4762a1bSJed Brown /*TEST
72c4762a1bSJed Brown 
73c4762a1bSJed Brown    test:
74c4762a1bSJed Brown       nsize: 3
75c4762a1bSJed Brown       args: -col 7
76c4762a1bSJed Brown 
77c4762a1bSJed Brown    test:
78c4762a1bSJed Brown       suffix: dense
79c4762a1bSJed Brown       nsize: 3
80a5225ed3SStefano Zampini       args: -col 7 -mat_type dense -vec_type {{mpi standard}}
81a5225ed3SStefano Zampini       filter: grep -v type
82a5225ed3SStefano Zampini 
83a5225ed3SStefano Zampini    test:
84a5225ed3SStefano Zampini       requires: cuda
85a5225ed3SStefano Zampini       suffix: dense_cuda
86a5225ed3SStefano Zampini       nsize: 3
87a5225ed3SStefano Zampini       output_file: output/ex60_dense.out
88a5225ed3SStefano Zampini       args: -col 7 -mat_type {{mpidense mpidensecuda}} -vec_type {{mpi standard cuda mpicuda}}
89a5225ed3SStefano Zampini       filter: grep -v type
90c4762a1bSJed Brown 
91c4762a1bSJed Brown TEST*/
92