xref: /petsc/src/mat/tests/ex91.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for sequential MatSBAIJ format. Derived from ex51.c\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   Mat          A, Atrans, sA, *submatA, *submatsA;
8c4762a1bSJed Brown   PetscInt     bs = 1, m = 43, ov = 1, i, j, k, *rows, *cols, M, nd = 5, *idx, mm, nn;
9c4762a1bSJed Brown   PetscMPIInt  size;
10c4762a1bSJed Brown   PetscScalar *vals, rval, one = 1.0;
11c4762a1bSJed Brown   IS          *is1, *is2;
12c4762a1bSJed Brown   PetscRandom  rand;
13c4762a1bSJed Brown   Vec          xx, s1, s2;
14c4762a1bSJed Brown   PetscReal    s1norm, s2norm, rnorm, tol = 10 * PETSC_SMALL;
15c4762a1bSJed Brown   PetscBool    flg;
16c4762a1bSJed Brown 
17327415f7SBarry Smith   PetscFunctionBeginUser;
18c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL));
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &m, NULL));
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL));
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   /* create a SeqBAIJ matrix A */
25c4762a1bSJed Brown   M = m * bs;
269566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF, bs, M, M, 1, NULL, &A));
279566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
289566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand));
299566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
30c4762a1bSJed Brown 
319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &rows));
329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &cols));
339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs * bs, &vals));
349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M, &idx));
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   /* Now set blocks of random values */
37c4762a1bSJed Brown   /* first, set diagonal blocks as zero */
38c4762a1bSJed Brown   for (j = 0; j < bs * bs; j++) vals[j] = 0.0;
39c4762a1bSJed Brown   for (i = 0; i < m; i++) {
409371c9d4SSatish Balay     cols[0] = i * bs;
419371c9d4SSatish Balay     rows[0] = i * bs;
42c4762a1bSJed Brown     for (j = 1; j < bs; j++) {
43c4762a1bSJed Brown       rows[j] = rows[j - 1] + 1;
44c4762a1bSJed Brown       cols[j] = cols[j - 1] + 1;
45c4762a1bSJed Brown     }
469566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, bs, rows, bs, cols, vals, ADD_VALUES));
47c4762a1bSJed Brown   }
48c4762a1bSJed Brown   /* second, add random blocks */
49c4762a1bSJed Brown   for (i = 0; i < 20 * bs; i++) {
509566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rand, &rval));
51c4762a1bSJed Brown     cols[0] = bs * (int)(PetscRealPart(rval) * m);
529566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rand, &rval));
53c4762a1bSJed Brown     rows[0] = bs * (int)(PetscRealPart(rval) * m);
54c4762a1bSJed Brown     for (j = 1; j < bs; j++) {
55c4762a1bSJed Brown       rows[j] = rows[j - 1] + 1;
56c4762a1bSJed Brown       cols[j] = cols[j - 1] + 1;
57c4762a1bSJed Brown     }
58c4762a1bSJed Brown 
59c4762a1bSJed Brown     for (j = 0; j < bs * bs; j++) {
609566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rand, &rval));
61c4762a1bSJed Brown       vals[j] = rval;
62c4762a1bSJed Brown     }
639566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, bs, rows, bs, cols, vals, ADD_VALUES));
64c4762a1bSJed Brown   }
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /* make A a symmetric matrix: A <- A^T + A */
709566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Atrans));
719566063dSJacob Faibussowitsch   PetscCall(MatAXPY(A, one, Atrans, DIFFERENT_NONZERO_PATTERN));
729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Atrans));
739566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Atrans));
749566063dSJacob Faibussowitsch   PetscCall(MatEqual(A, Atrans, &flg));
7528b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "A+A^T is non-symmetric");
769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Atrans));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /* create a SeqSBAIJ matrix sA (= A) */
799566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE));
809566063dSJacob Faibussowitsch   PetscCall(MatConvert(A, MATSEQSBAIJ, MAT_INITIAL_MATRIX, &sA));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /* Test sA==A through MatMult() */
83c4762a1bSJed Brown   for (i = 0; i < nd; i++) {
849566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A, &mm, &nn));
859566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx));
869566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s1));
879566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s2));
88c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
899566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx, rand));
909566063dSJacob Faibussowitsch       PetscCall(MatMult(A, xx, s1));
919566063dSJacob Faibussowitsch       PetscCall(MatMult(sA, xx, s2));
929566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1, NORM_2, &s1norm));
939566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2, NORM_2, &s2norm));
94c4762a1bSJed Brown       rnorm = s2norm - s1norm;
9548a46eb9SPierre Jolivet       if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", (double)s1norm, (double)s2norm));
96c4762a1bSJed Brown     }
979566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
989566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
999566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
100c4762a1bSJed Brown   }
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /* Test MatIncreaseOverlap() */
1039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd, &is1));
1049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd, &is2));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   for (i = 0; i < nd; i++) {
1079566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rand, &rval));
108c4762a1bSJed Brown     size = (int)(PetscRealPart(rval) * m);
109c4762a1bSJed Brown     for (j = 0; j < size; j++) {
1109566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rand, &rval));
111c4762a1bSJed Brown       idx[j * bs] = bs * (int)(PetscRealPart(rval) * m);
112c4762a1bSJed Brown       for (k = 1; k < bs; k++) idx[j * bs + k] = idx[j * bs] + k;
113c4762a1bSJed Brown     }
1149566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size * bs, idx, PETSC_COPY_VALUES, is1 + i));
1159566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, size * bs, idx, PETSC_COPY_VALUES, is2 + i));
116c4762a1bSJed Brown   }
117c4762a1bSJed Brown   /* for debugging */
118c4762a1bSJed Brown   /*
1199566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_SELF));
1209566063dSJacob Faibussowitsch   PetscCall(MatView(sA,PETSC_VIEWER_STDOUT_SELF));
121c4762a1bSJed Brown   */
122c4762a1bSJed Brown 
1239566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(A, nd, is1, ov));
1249566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(sA, nd, is2, ov));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   for (i = 0; i < nd; ++i) {
1279566063dSJacob Faibussowitsch     PetscCall(ISSort(is1[i]));
1289566063dSJacob Faibussowitsch     PetscCall(ISSort(is2[i]));
129c4762a1bSJed Brown   }
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   for (i = 0; i < nd; ++i) {
1329566063dSJacob Faibussowitsch     PetscCall(ISEqual(is1[i], is2[i], &flg));
13328b400f6SJacob Faibussowitsch     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "i=%" PetscInt_FMT ", is1 != is2", i);
134c4762a1bSJed Brown   }
135c4762a1bSJed Brown 
1369566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_INITIAL_MATRIX, &submatA));
1379566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(sA, nd, is2, is2, MAT_INITIAL_MATRIX, &submatsA));
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   /* Test MatMult() */
140c4762a1bSJed Brown   for (i = 0; i < nd; i++) {
1419566063dSJacob Faibussowitsch     PetscCall(MatGetSize(submatA[i], &mm, &nn));
1429566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx));
1439566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s1));
1449566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s2));
145c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
1469566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx, rand));
1479566063dSJacob Faibussowitsch       PetscCall(MatMult(submatA[i], xx, s1));
1489566063dSJacob Faibussowitsch       PetscCall(MatMult(submatsA[i], xx, s2));
1499566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1, NORM_2, &s1norm));
1509566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2, NORM_2, &s2norm));
151c4762a1bSJed Brown       rnorm = s2norm - s1norm;
15248a46eb9SPierre Jolivet       if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", (double)s1norm, (double)s2norm));
153c4762a1bSJed Brown     }
1549566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
1559566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
1569566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
157c4762a1bSJed Brown   }
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */
1609566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_REUSE_MATRIX, &submatA));
1619566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(sA, nd, is2, is2, MAT_REUSE_MATRIX, &submatsA));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* Test MatMult() */
164c4762a1bSJed Brown   for (i = 0; i < nd; i++) {
1659566063dSJacob Faibussowitsch     PetscCall(MatGetSize(submatA[i], &mm, &nn));
1669566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx));
1679566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s1));
1689566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s2));
169c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
1709566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx, rand));
1719566063dSJacob Faibussowitsch       PetscCall(MatMult(submatA[i], xx, s1));
1729566063dSJacob Faibussowitsch       PetscCall(MatMult(submatsA[i], xx, s2));
1739566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1, NORM_2, &s1norm));
1749566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2, NORM_2, &s2norm));
175c4762a1bSJed Brown       rnorm = s2norm - s1norm;
17648a46eb9SPierre Jolivet       if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", (double)s1norm, (double)s2norm));
177c4762a1bSJed Brown     }
1789566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
1799566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
1809566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
181c4762a1bSJed Brown   }
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   /* Free allocated memory */
184c4762a1bSJed Brown   for (i = 0; i < nd; ++i) {
1859566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is1[i]));
1869566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is2[i]));
187c4762a1bSJed Brown   }
1889566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd, &submatA));
1899566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd, &submatsA));
190c4762a1bSJed Brown 
1919566063dSJacob Faibussowitsch   PetscCall(PetscFree(is1));
1929566063dSJacob Faibussowitsch   PetscCall(PetscFree(is2));
1939566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
1949566063dSJacob Faibussowitsch   PetscCall(PetscFree(rows));
1959566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols));
1969566063dSJacob Faibussowitsch   PetscCall(PetscFree(vals));
1979566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1989566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sA));
1999566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
2009566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
201b122ec5aSJacob Faibussowitsch   return 0;
202c4762a1bSJed Brown }
203c4762a1bSJed Brown 
204c4762a1bSJed Brown /*TEST
205c4762a1bSJed Brown 
206c4762a1bSJed Brown    test:
207c4762a1bSJed Brown       args: -ov 2
208*3886731fSPierre Jolivet       output_file: output/empty.out
209c4762a1bSJed Brown 
210c4762a1bSJed Brown TEST*/
211