xref: /petsc/src/mat/tests/ex56.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Test the use of MatSetValuesBlocked(), MatZeroRows() for rectangular MatBAIJ matrix, test MatSetValuesBlocked() for MatSBAIJ matrix (-test_mat_sbaij).";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7*d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   Mat         A;
9c4762a1bSJed Brown   PetscInt    bs = 3, m = 4, n = 6, i, j, val = 10, row[2], col[3], eval, rstart;
10c4762a1bSJed Brown   PetscMPIInt size, rank;
11c4762a1bSJed Brown   PetscScalar x[6][9], y[3][3], one = 1.0;
12c4762a1bSJed Brown   PetscBool   flg, testsbaij = PETSC_FALSE;
13c4762a1bSJed Brown 
14327415f7SBarry Smith   PetscFunctionBeginUser;
159566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
18c4762a1bSJed Brown 
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_mat_sbaij", &testsbaij));
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   if (testsbaij) {
229566063dSJacob Faibussowitsch     PetscCall(MatCreateSBAIJ(PETSC_COMM_WORLD, bs, m * bs, n * bs, PETSC_DECIDE, PETSC_DECIDE, 1, NULL, 1, NULL, &A));
23c4762a1bSJed Brown   } else {
249566063dSJacob Faibussowitsch     PetscCall(MatCreateBAIJ(PETSC_COMM_WORLD, bs, m * bs, n * bs, PETSC_DECIDE, PETSC_DECIDE, 1, NULL, 1, NULL, &A));
25c4762a1bSJed Brown   }
269566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
27c4762a1bSJed Brown   eval = 9;
28c4762a1bSJed Brown 
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-ass_extern", &flg));
30c4762a1bSJed Brown   if (flg && (size != 1)) rstart = m * ((rank + 1) % size);
31c4762a1bSJed Brown   else rstart = m * (rank);
32c4762a1bSJed Brown 
339371c9d4SSatish Balay   row[0] = rstart + 0;
349371c9d4SSatish Balay   row[1] = rstart + 2;
359371c9d4SSatish Balay   col[0] = rstart + 0;
369371c9d4SSatish Balay   col[1] = rstart + 1;
379371c9d4SSatish Balay   col[2] = rstart + 3;
38c4762a1bSJed Brown   for (i = 0; i < 6; i++) {
39c4762a1bSJed Brown     for (j = 0; j < 9; j++) x[i][j] = (PetscScalar)val++;
40c4762a1bSJed Brown   }
41c4762a1bSJed Brown 
429566063dSJacob Faibussowitsch   PetscCall(MatSetValuesBlocked(A, 2, row, 3, col, &x[0][0], INSERT_VALUES));
43c4762a1bSJed Brown 
449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /*
48c4762a1bSJed Brown   This option does not work for rectangular matrices
499566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
50c4762a1bSJed Brown   */
51c4762a1bSJed Brown 
529566063dSJacob Faibussowitsch   PetscCall(MatSetValuesBlocked(A, 2, row, 3, col, &x[0][0], INSERT_VALUES));
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   /* Do another MatSetValues to test the case when only one local block is specified */
55c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
56c4762a1bSJed Brown     for (j = 0; j < 3; j++) y[i][j] = (PetscScalar)(10 + i * eval + j);
57c4762a1bSJed Brown   }
589566063dSJacob Faibussowitsch   PetscCall(MatSetValuesBlocked(A, 1, row, 1, col, &y[0][0], INSERT_VALUES));
599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-zero_rows", &flg));
63c4762a1bSJed Brown   if (flg) {
64c4762a1bSJed Brown     col[0] = rstart * bs + 0;
65c4762a1bSJed Brown     col[1] = rstart * bs + 1;
66c4762a1bSJed Brown     col[2] = rstart * bs + 2;
679566063dSJacob Faibussowitsch     PetscCall(MatZeroRows(A, 3, col, one, 0, 0));
68c4762a1bSJed Brown   }
69c4762a1bSJed Brown 
709566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
739566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
74b122ec5aSJacob Faibussowitsch   return 0;
75c4762a1bSJed Brown }
76c4762a1bSJed Brown 
77c4762a1bSJed Brown /*TEST
78c4762a1bSJed Brown 
79c4762a1bSJed Brown    test:
808cc725e6SPierre Jolivet       filter: grep -v " MPI process"
81c4762a1bSJed Brown 
82c4762a1bSJed Brown    test:
83c4762a1bSJed Brown       suffix: 4
84c4762a1bSJed Brown       nsize: 3
85c4762a1bSJed Brown       args: -ass_extern
868cc725e6SPierre Jolivet       filter: grep -v " MPI process"
87c4762a1bSJed Brown 
88c4762a1bSJed Brown    test:
89c4762a1bSJed Brown       suffix: 5
90c4762a1bSJed Brown       nsize: 3
91c4762a1bSJed Brown       args: -ass_extern -zero_rows
928cc725e6SPierre Jolivet       filter: grep -v " MPI process"
93c4762a1bSJed Brown 
94c4762a1bSJed Brown TEST*/
95