1c4762a1bSJed Brown static char help[] = "Tests MatGetColumnVector()."; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 6d71ae5a4SJacob Faibussowitsch { 7c4762a1bSJed Brown Mat C; 8c4762a1bSJed Brown PetscInt i, j, m = 3, n = 2, Ii, J, col = 0; 9c4762a1bSJed Brown PetscMPIInt size, rank; 10c4762a1bSJed Brown PetscScalar v; 11c4762a1bSJed Brown Vec yy; 12c4762a1bSJed Brown 13327415f7SBarry Smith PetscFunctionBeginUser; 14*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 159566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-col", &col, NULL)); 16c4762a1bSJed Brown 179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 19c4762a1bSJed Brown n = 2 * size; 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* create the matrix for the five point stencil, YET AGAIN*/ 229566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 239566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n)); 249566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 259566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(C, 5, NULL)); 269566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(C, 5, NULL, 5, NULL)); 279566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 28c4762a1bSJed Brown 29c4762a1bSJed Brown for (i = 0; i < m; i++) { 30c4762a1bSJed Brown for (j = 2 * rank; j < 2 * rank + 2; j++) { 319371c9d4SSatish Balay v = -1.0; 329371c9d4SSatish Balay Ii = j + n * i; 339371c9d4SSatish Balay if (i > 0) { 349371c9d4SSatish Balay J = Ii - n; 359371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 369371c9d4SSatish Balay } 379371c9d4SSatish Balay if (i < m - 1) { 389371c9d4SSatish Balay J = Ii + n; 399371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 409371c9d4SSatish Balay } 419371c9d4SSatish Balay if (j > 0) { 429371c9d4SSatish Balay J = Ii - 1; 439371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 449371c9d4SSatish Balay } 459371c9d4SSatish Balay if (j < n - 1) { 469371c9d4SSatish Balay J = Ii + 1; 479371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 489371c9d4SSatish Balay } 499371c9d4SSatish Balay v = 4.0; 509371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 51c4762a1bSJed Brown } 52c4762a1bSJed Brown } 539566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 549566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 559566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 56c4762a1bSJed Brown 579566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C, NULL, &yy)); 589566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(yy)); 59c4762a1bSJed Brown 609566063dSJacob Faibussowitsch PetscCall(MatGetColumnVector(C, yy, col)); 61c4762a1bSJed Brown 629566063dSJacob Faibussowitsch PetscCall(VecView(yy, PETSC_VIEWER_STDOUT_WORLD)); 63c4762a1bSJed Brown 649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&yy)); 659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 669566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 67b122ec5aSJacob Faibussowitsch return 0; 68c4762a1bSJed Brown } 69c4762a1bSJed Brown 70c4762a1bSJed Brown /*TEST 71c4762a1bSJed Brown 72c4762a1bSJed Brown test: 73c4762a1bSJed Brown nsize: 3 74c4762a1bSJed Brown args: -col 7 75c4762a1bSJed Brown 76c4762a1bSJed Brown test: 77c4762a1bSJed Brown suffix: dense 78c4762a1bSJed Brown nsize: 3 79a5225ed3SStefano Zampini args: -col 7 -mat_type dense -vec_type {{mpi standard}} 80a5225ed3SStefano Zampini filter: grep -v type 81a5225ed3SStefano Zampini 82a5225ed3SStefano Zampini test: 83a5225ed3SStefano Zampini requires: cuda 84a5225ed3SStefano Zampini suffix: dense_cuda 85a5225ed3SStefano Zampini nsize: 3 86a5225ed3SStefano Zampini output_file: output/ex60_dense.out 87a5225ed3SStefano Zampini args: -col 7 -mat_type {{mpidense mpidensecuda}} -vec_type {{mpi standard cuda mpicuda}} 88a5225ed3SStefano Zampini filter: grep -v type 89c4762a1bSJed Brown 90c4762a1bSJed Brown TEST*/ 91