xref: /petsc/src/mat/tests/ex54.c (revision 4b5966dda8e0e9d52f06e77e04a44e922c8134a2)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for parallel AIJ and BAIJ formats.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7d71ae5a4SJacob Faibussowitsch {
8*4b5966ddSBarry Smith   Mat          E, A, B, *submatA, *submatB;
9c4762a1bSJed Brown   PetscInt     bs = 1, m = 11, ov = 1, i, j, k, *rows, *cols, nd = 5, *idx, rstart, rend, sz, mm, nn, M, N, Mbs;
10c4762a1bSJed Brown   PetscMPIInt  size, rank;
11c4762a1bSJed Brown   PetscScalar *vals, rval;
12c4762a1bSJed Brown   IS          *is1, *is2;
13c4762a1bSJed Brown   PetscRandom  rdm;
14c4762a1bSJed Brown   Vec          xx, s1, s2;
15c4762a1bSJed Brown   PetscReal    s1norm, s2norm, rnorm, tol = 100 * PETSC_SMALL;
1682b5ce2aSStefano Zampini   PetscBool    flg, test_nd0 = PETSC_FALSE, emptynd;
17c4762a1bSJed Brown 
18327415f7SBarry Smith   PetscFunctionBeginUser;
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
22c4762a1bSJed Brown 
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL));
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &m, NULL));
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL));
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL));
279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_nd0", &test_nd0, NULL));
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   /* Create a AIJ matrix A */
309566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
319566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, m * bs, m * bs, PETSC_DECIDE, PETSC_DECIDE));
329566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATAIJ));
33*4b5966ddSBarry Smith   PetscCall(MatSetBlockSize(A, bs));
349566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(A, PETSC_DEFAULT, NULL));
359566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(A, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL));
369566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
379566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   /* Create a BAIJ matrix B */
409566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
419566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, m * bs, m * bs, PETSC_DECIDE, PETSC_DECIDE));
429566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, MATBAIJ));
439566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(B, bs, PETSC_DEFAULT, NULL));
449566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocation(B, bs, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL));
459566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
469566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
499566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
529566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
53c4762a1bSJed Brown   Mbs = M / bs;
54c4762a1bSJed Brown 
559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &rows));
569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &cols));
579566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs * bs, &vals));
589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M, &idx));
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Now set blocks of values */
61c4762a1bSJed Brown   for (i = 0; i < 40 * bs; i++) {
62*4b5966ddSBarry Smith     PetscInt nr = 1, nc = 1;
639566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
64c4762a1bSJed Brown     cols[0] = bs * (int)(PetscRealPart(rval) * Mbs);
659566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
66c4762a1bSJed Brown     rows[0] = rstart + bs * (int)(PetscRealPart(rval) * m);
67c4762a1bSJed Brown     for (j = 1; j < bs; j++) {
68*4b5966ddSBarry Smith       PetscCall(PetscRandomGetValue(rdm, &rval));
69*4b5966ddSBarry Smith       if (PetscRealPart(rval) > .5) rows[nr++] = rows[0] + j - 1;
70*4b5966ddSBarry Smith     }
71*4b5966ddSBarry Smith     for (j = 1; j < bs; j++) {
72*4b5966ddSBarry Smith       PetscCall(PetscRandomGetValue(rdm, &rval));
73*4b5966ddSBarry Smith       if (PetscRealPart(rval) > .5) cols[nc++] = cols[0] + j - 1;
74c4762a1bSJed Brown     }
75c4762a1bSJed Brown 
76*4b5966ddSBarry Smith     for (j = 0; j < nr * nc; j++) {
779566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rdm, &rval));
78c4762a1bSJed Brown       vals[j] = rval;
79c4762a1bSJed Brown     }
80*4b5966ddSBarry Smith     PetscCall(MatSetValues(A, nr, rows, nc, cols, vals, ADD_VALUES));
81*4b5966ddSBarry Smith     PetscCall(MatSetValues(B, nr, rows, nc, cols, vals, ADD_VALUES));
82c4762a1bSJed Brown   }
839566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
869566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
87c4762a1bSJed Brown 
88*4b5966ddSBarry Smith   /* Test MatConvert_MPIAIJ_MPIBAIJ handles incompletely filled blocks */
89*4b5966ddSBarry Smith   PetscCall(MatConvert(A, MATBAIJ, MAT_INITIAL_MATRIX, &E));
90*4b5966ddSBarry Smith   PetscCall(MatDestroy(&E));
91*4b5966ddSBarry Smith 
92c4762a1bSJed Brown   /* Test MatIncreaseOverlap() */
939566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd, &is1));
949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd, &is2));
95c4762a1bSJed Brown 
9682b5ce2aSStefano Zampini   emptynd = PETSC_FALSE;
9782b5ce2aSStefano Zampini   if (rank == 0 && test_nd0) emptynd = PETSC_TRUE; /* test case */
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   for (i = 0; i < nd; i++) {
1009566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
101c4762a1bSJed Brown     sz = (int)(PetscRealPart(rval) * m);
102c4762a1bSJed Brown     for (j = 0; j < sz; j++) {
1039566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rdm, &rval));
104c4762a1bSJed Brown       idx[j * bs] = bs * (int)(PetscRealPart(rval) * Mbs);
105c4762a1bSJed Brown       for (k = 1; k < bs; k++) idx[j * bs + k] = idx[j * bs] + k;
106c4762a1bSJed Brown     }
10782b5ce2aSStefano Zampini     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, emptynd ? 0 : sz * bs, idx, PETSC_COPY_VALUES, is1 + i));
10882b5ce2aSStefano Zampini     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, emptynd ? 0 : sz * bs, idx, PETSC_COPY_VALUES, is2 + i));
109c4762a1bSJed Brown   }
1109566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(A, nd, is1, ov));
1119566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(B, nd, is2, ov));
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   for (i = 0; i < nd; ++i) {
1149566063dSJacob Faibussowitsch     PetscCall(ISEqual(is1[i], is2[i], &flg));
115c4762a1bSJed Brown 
11648a46eb9SPierre Jolivet     if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "i=%" PetscInt_FMT ", flg=%d :bs=%" PetscInt_FMT " m=%" PetscInt_FMT " ov=%" PetscInt_FMT " nd=%" PetscInt_FMT " np=%d\n", i, flg, bs, m, ov, nd, size));
117c4762a1bSJed Brown   }
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   for (i = 0; i < nd; ++i) {
1209566063dSJacob Faibussowitsch     PetscCall(ISSort(is1[i]));
1219566063dSJacob Faibussowitsch     PetscCall(ISSort(is2[i]));
122c4762a1bSJed Brown   }
123c4762a1bSJed Brown 
1249566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_INITIAL_MATRIX, &submatB));
1259566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_INITIAL_MATRIX, &submatA));
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   /* Test MatMult() */
128c4762a1bSJed Brown   for (i = 0; i < nd; i++) {
1299566063dSJacob Faibussowitsch     PetscCall(MatGetSize(submatA[i], &mm, &nn));
1309566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx));
1319566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s1));
1329566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s2));
133c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
1349566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx, rdm));
1359566063dSJacob Faibussowitsch       PetscCall(MatMult(submatA[i], xx, s1));
1369566063dSJacob Faibussowitsch       PetscCall(MatMult(submatB[i], xx, s2));
1379566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1, NORM_2, &s1norm));
1389566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2, NORM_2, &s2norm));
139c4762a1bSJed Brown       rnorm = s2norm - s1norm;
14048a46eb9SPierre Jolivet       if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", rank, (double)s1norm, (double)s2norm));
141c4762a1bSJed Brown     }
1429566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
1439566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
1449566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
145c4762a1bSJed Brown   }
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */
1489566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_REUSE_MATRIX, &submatA));
1499566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_REUSE_MATRIX, &submatB));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* Test MatMult() */
152c4762a1bSJed Brown   for (i = 0; i < nd; i++) {
1539566063dSJacob Faibussowitsch     PetscCall(MatGetSize(submatA[i], &mm, &nn));
1549566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx));
1559566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s1));
1569566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s2));
157c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
1589566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx, rdm));
1599566063dSJacob Faibussowitsch       PetscCall(MatMult(submatA[i], xx, s1));
1609566063dSJacob Faibussowitsch       PetscCall(MatMult(submatB[i], xx, s2));
1619566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1, NORM_2, &s1norm));
1629566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2, NORM_2, &s2norm));
163c4762a1bSJed Brown       rnorm = s2norm - s1norm;
16448a46eb9SPierre Jolivet       if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", rank, (double)s1norm, (double)s2norm));
165c4762a1bSJed Brown     }
1669566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
1679566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
1689566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
169c4762a1bSJed Brown   }
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   /* Free allocated memory */
172c4762a1bSJed Brown   for (i = 0; i < nd; ++i) {
1739566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is1[i]));
1749566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is2[i]));
175c4762a1bSJed Brown   }
1769566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd, &submatA));
1779566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd, &submatB));
178c4762a1bSJed Brown 
1799566063dSJacob Faibussowitsch   PetscCall(PetscFree(is1));
1809566063dSJacob Faibussowitsch   PetscCall(PetscFree(is2));
1819566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
1829566063dSJacob Faibussowitsch   PetscCall(PetscFree(rows));
1839566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols));
1849566063dSJacob Faibussowitsch   PetscCall(PetscFree(vals));
1859566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1869566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1879566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
1889566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
189b122ec5aSJacob Faibussowitsch   return 0;
190c4762a1bSJed Brown }
191c4762a1bSJed Brown 
192c4762a1bSJed Brown /*TEST
193c4762a1bSJed Brown 
194c4762a1bSJed Brown    test:
195c4762a1bSJed Brown       nsize: {{1 3}}
19682b5ce2aSStefano Zampini       args: -mat_block_size {{1 3 4 6 8}} -ov {{1 3}} -mat_size {{11 13}} -nd 7
197c4762a1bSJed Brown       output_file: output/ex54.out
198c4762a1bSJed Brown 
199c4762a1bSJed Brown    test:
200c4762a1bSJed Brown       suffix: 2
201c4762a1bSJed Brown       args: -nd 2 -test_nd0
202c4762a1bSJed Brown       output_file: output/ex54.out
203c4762a1bSJed Brown 
204c4762a1bSJed Brown    test:
205c4762a1bSJed Brown       suffix: 3
206c4762a1bSJed Brown       nsize: 3
207c4762a1bSJed Brown       args: -nd 2 -test_nd0
208c4762a1bSJed Brown       output_file: output/ex54.out
209c4762a1bSJed Brown 
210c4762a1bSJed Brown TEST*/
211