1c4762a1bSJed Brown static char help[] = "Tests MatGetRowIJ for SeqAIJ, SeqBAIJ and SeqSBAIJ\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown
DumpCSR(Mat A,PetscInt shift,PetscBool symmetric,PetscBool compressed)5d71ae5a4SJacob Faibussowitsch PetscErrorCode DumpCSR(Mat A, PetscInt shift, PetscBool symmetric, PetscBool compressed)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown MatType type;
8c4762a1bSJed Brown PetscInt i, j, nr, bs = 1;
9c4762a1bSJed Brown const PetscInt *ia, *ja;
10c4762a1bSJed Brown PetscBool done, isseqbaij, isseqsbaij;
11c4762a1bSJed Brown
12c4762a1bSJed Brown PetscFunctionBeginUser;
139566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQBAIJ, &isseqbaij));
149566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQSBAIJ, &isseqsbaij));
1548a46eb9SPierre Jolivet if (isseqbaij || isseqsbaij) PetscCall(MatGetBlockSize(A, &bs));
169566063dSJacob Faibussowitsch PetscCall(MatGetType(A, &type));
179566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(A, shift, symmetric, compressed, &nr, &ia, &ja, &done));
189566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "===========================================================\n"));
199566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "CSR for %s: shift %" PetscInt_FMT " symmetric %" PetscInt_FMT " compressed %" PetscInt_FMT "\n", type, shift, (PetscInt)symmetric, (PetscInt)compressed));
20c4762a1bSJed Brown for (i = 0; i < nr; i++) {
219566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT ":", i + shift));
2248a46eb9SPierre Jolivet for (j = ia[i]; j < ia[i + 1]; j++) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscInt_FMT, ja[j - shift]));
239566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
24c4762a1bSJed Brown }
259566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(A, shift, symmetric, compressed, &nr, &ia, &ja, &done));
263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
27c4762a1bSJed Brown }
28c4762a1bSJed Brown
main(int argc,char ** args)29d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
30d71ae5a4SJacob Faibussowitsch {
31c4762a1bSJed Brown Mat A, B, C;
32c4762a1bSJed Brown PetscInt i, j, k, m = 3, n = 3, bs = 1;
33c4762a1bSJed Brown PetscMPIInt size;
34c4762a1bSJed Brown
35327415f7SBarry Smith PetscFunctionBeginUser;
36*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
38be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
419566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL));
42c4762a1bSJed Brown /* adjust sizes by block size */
43c4762a1bSJed Brown if (m % bs) m += bs - m % bs;
44c4762a1bSJed Brown if (n % bs) n += bs - n % bs;
45c4762a1bSJed Brown
469566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &A));
479566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, m * n, m * n, PETSC_DECIDE, PETSC_DECIDE));
489566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(A, bs));
499566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSEQAIJ));
509566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
519566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &B));
529566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, m * n, m * n, PETSC_DECIDE, PETSC_DECIDE));
539566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(B, bs));
549566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATSEQBAIJ));
559566063dSJacob Faibussowitsch PetscCall(MatSetUp(B));
569566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &C));
579566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, m * n, m * n, PETSC_DECIDE, PETSC_DECIDE));
589566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(C, bs));
599566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATSEQSBAIJ));
609566063dSJacob Faibussowitsch PetscCall(MatSetUp(C));
619566063dSJacob Faibussowitsch PetscCall(MatSetOption(C, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
62c4762a1bSJed Brown
63c4762a1bSJed Brown for (i = 0; i < m; i++) {
64c4762a1bSJed Brown for (j = 0; j < n; j++) {
65c4762a1bSJed Brown PetscScalar v = -1.0;
66c4762a1bSJed Brown PetscInt Ii = j + n * i, J;
67c4762a1bSJed Brown J = Ii - n;
68c4762a1bSJed Brown if (J >= 0) {
699566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
709566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &Ii, 1, &J, &v, INSERT_VALUES));
719566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
72c4762a1bSJed Brown }
73c4762a1bSJed Brown J = Ii + n;
74c4762a1bSJed Brown if (J < m * n) {
759566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
769566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &Ii, 1, &J, &v, INSERT_VALUES));
779566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
78c4762a1bSJed Brown }
79c4762a1bSJed Brown J = Ii - 1;
80c4762a1bSJed Brown if (J >= 0) {
819566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
829566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &Ii, 1, &J, &v, INSERT_VALUES));
839566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
84c4762a1bSJed Brown }
85c4762a1bSJed Brown J = Ii + 1;
86c4762a1bSJed Brown if (J < m * n) {
879566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
889566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &Ii, 1, &J, &v, INSERT_VALUES));
899566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
90c4762a1bSJed Brown }
91c4762a1bSJed Brown v = 4.0;
929566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
939566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
949566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
95c4762a1bSJed Brown }
96c4762a1bSJed Brown }
979566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
989566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
999566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1009566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
1019566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
1029566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
103c4762a1bSJed Brown
104c4762a1bSJed Brown /* test MatGetRowIJ for the three Mat types */
1059566063dSJacob Faibussowitsch PetscCall(MatView(A, NULL));
1069566063dSJacob Faibussowitsch PetscCall(MatView(B, NULL));
1079566063dSJacob Faibussowitsch PetscCall(MatView(C, NULL));
108c4762a1bSJed Brown for (i = 0; i < 2; i++) {
109c4762a1bSJed Brown PetscInt shift = i;
110c4762a1bSJed Brown for (j = 0; j < 2; j++) {
111c4762a1bSJed Brown PetscBool symmetric = ((j > 0) ? PETSC_FALSE : PETSC_TRUE);
112c4762a1bSJed Brown for (k = 0; k < 2; k++) {
113c4762a1bSJed Brown PetscBool compressed = ((k > 0) ? PETSC_FALSE : PETSC_TRUE);
1149566063dSJacob Faibussowitsch PetscCall(DumpCSR(A, shift, symmetric, compressed));
1159566063dSJacob Faibussowitsch PetscCall(DumpCSR(B, shift, symmetric, compressed));
1169566063dSJacob Faibussowitsch PetscCall(DumpCSR(C, shift, symmetric, compressed));
117c4762a1bSJed Brown }
118c4762a1bSJed Brown }
119c4762a1bSJed Brown }
1209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
1219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
1229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
1239566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
124b122ec5aSJacob Faibussowitsch return 0;
125c4762a1bSJed Brown }
126c4762a1bSJed Brown
127c4762a1bSJed Brown /*TEST
128c4762a1bSJed Brown
129c4762a1bSJed Brown test:
130c4762a1bSJed Brown
131c4762a1bSJed Brown test:
132c4762a1bSJed Brown suffix: 2
133c4762a1bSJed Brown args: -bs 2
134c4762a1bSJed Brown
135c4762a1bSJed Brown TEST*/
136