xref: /petsc/src/mat/tests/ex52.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests various routines in MatMPIBAIJ format.\n";
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    m = 2, bs = 1, M, row, col, start, end, i, j, k;
10c4762a1bSJed Brown   PetscMPIInt rank, size;
11c4762a1bSJed Brown   PetscScalar data = 100;
12c4762a1bSJed Brown   PetscBool   flg;
13c4762a1bSJed Brown 
14327415f7SBarry Smith   PetscFunctionBeginUser;
159566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   /* Test MatSetValues() and MatGetValues() */
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL));
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &m, NULL));
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   M = m * bs * size;
249566063dSJacob Faibussowitsch   PetscCall(MatCreateBAIJ(PETSC_COMM_WORLD, bs, PETSC_DECIDE, PETSC_DECIDE, M, M, PETSC_DECIDE, NULL, PETSC_DECIDE, NULL, &A));
25c4762a1bSJed Brown 
269566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &start, &end));
279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-column_oriented", &flg));
281baa6e33SBarry Smith   if (flg) PetscCall(MatSetOption(A, MAT_ROW_ORIENTED, PETSC_FALSE));
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   /* inproc assembly */
31c4762a1bSJed Brown   for (row = start; row < end; row++) {
3248a46eb9SPierre Jolivet     for (col = start; col < end; col++, data += 1) PetscCall(MatSetValues(A, 1, &row, 1, &col, &data, INSERT_VALUES));
33c4762a1bSJed Brown   }
349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
359566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* offproc assembly */
38c4762a1bSJed Brown   data = 5.0;
39c4762a1bSJed Brown   row  = (M + start - 1) % M;
4048a46eb9SPierre Jolivet   for (col = 0; col < M; col++) PetscCall(MatSetValues(A, 1, &row, 1, &col, &data, ADD_VALUES));
419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Test MatSetValuesBlocked() */
459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_setvaluesblocked", &flg));
46c4762a1bSJed Brown   if (flg) {
47c4762a1bSJed Brown     PetscScalar *bval;
48c4762a1bSJed Brown     row /= bs;
49c4762a1bSJed Brown     col = start / bs;
509566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bs * bs, &bval));
51c4762a1bSJed Brown     k = 1;
52a5b23f4aSJose E. Roman     /* row oriented - default */
53c4762a1bSJed Brown     for (i = 0; i < bs; i++) {
54c4762a1bSJed Brown       for (j = 0; j < bs; j++) {
559371c9d4SSatish Balay         bval[i * bs + j] = (PetscScalar)k;
569371c9d4SSatish Balay         k++;
57c4762a1bSJed Brown       }
58c4762a1bSJed Brown     }
599566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A, 1, &row, 1, &col, bval, INSERT_VALUES));
609566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
619566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
629566063dSJacob Faibussowitsch     PetscCall(PetscFree(bval));
63c4762a1bSJed Brown   }
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
669566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
679566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
68b122ec5aSJacob Faibussowitsch   return 0;
69c4762a1bSJed Brown }
70c4762a1bSJed Brown 
71c4762a1bSJed Brown /*TEST
72c4762a1bSJed Brown 
73c4762a1bSJed Brown    test:
74c4762a1bSJed Brown       suffix: 1
75c4762a1bSJed Brown       nsize: 3
76c4762a1bSJed Brown       args: -mat_block_size 2 -test_setvaluesblocked
77c4762a1bSJed Brown 
78c4762a1bSJed Brown    test:
79c4762a1bSJed Brown       suffix: 2
80c4762a1bSJed Brown       nsize: 3
81c4762a1bSJed Brown       args: -mat_block_size 2 -test_setvaluesblocked -column_oriented
82c4762a1bSJed Brown 
83c4762a1bSJed Brown    test:
84c4762a1bSJed Brown       suffix: 3
85c4762a1bSJed Brown       nsize: 3
86c4762a1bSJed Brown       args: -mat_block_size 1 -test_setvaluesblocked
87c4762a1bSJed Brown 
88c4762a1bSJed Brown    test:
89c4762a1bSJed Brown       suffix: 4
90c4762a1bSJed Brown       nsize: 3
91c4762a1bSJed Brown       args: -mat_block_size 1 -test_setvaluesblocked -column_oriented
92c4762a1bSJed Brown 
93c4762a1bSJed Brown TEST*/
94