1c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for MatBAIJ format.\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 {
75a2b941aSBarry Smith Mat A, B, E, Bt, *submatA, *submatB;
8c4762a1bSJed Brown PetscInt bs = 1, m = 43, ov = 1, i, j, k, *rows, *cols, M, nd = 5, *idx, mm, nn, lsize;
9c4762a1bSJed Brown PetscScalar *vals, rval;
10c4762a1bSJed Brown IS *is1, *is2;
11c4762a1bSJed Brown PetscRandom rdm;
12c4762a1bSJed Brown Vec xx, s1, s2;
13c4762a1bSJed Brown PetscReal s1norm, s2norm, rnorm, tol = PETSC_SQRT_MACHINE_EPSILON;
14c4762a1bSJed Brown PetscBool flg;
15c4762a1bSJed Brown
16327415f7SBarry Smith PetscFunctionBeginUser;
17c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
189566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL));
199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &m, NULL));
209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL));
219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL));
22c4762a1bSJed Brown M = m * bs;
23c4762a1bSJed Brown
249566063dSJacob Faibussowitsch PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF, bs, M, M, 1, NULL, &A));
259566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
269566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, M, M, 15, NULL, &B));
274b5966ddSBarry Smith PetscCall(MatSetBlockSize(B, bs));
289566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
299566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
309566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm));
31c4762a1bSJed Brown
329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs, &rows));
339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs, &cols));
349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * bs, &vals));
359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M, &idx));
36c4762a1bSJed Brown
37c4762a1bSJed Brown /* Now set blocks of values */
38c4762a1bSJed Brown for (i = 0; i < 20 * bs; i++) {
394b5966ddSBarry Smith PetscInt nr = 1, nc = 1;
409566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval));
41c4762a1bSJed Brown cols[0] = bs * (int)(PetscRealPart(rval) * m);
429566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval));
43c4762a1bSJed Brown rows[0] = bs * (int)(PetscRealPart(rval) * m);
44c4762a1bSJed Brown for (j = 1; j < bs; j++) {
454b5966ddSBarry Smith PetscCall(PetscRandomGetValue(rdm, &rval));
464b5966ddSBarry Smith if (PetscRealPart(rval) > .5) rows[nr++] = rows[0] + j - 1;
474b5966ddSBarry Smith }
484b5966ddSBarry Smith for (j = 1; j < bs; j++) {
494b5966ddSBarry Smith PetscCall(PetscRandomGetValue(rdm, &rval));
504b5966ddSBarry Smith if (PetscRealPart(rval) > .5) cols[nc++] = cols[0] + j - 1;
51c4762a1bSJed Brown }
52c4762a1bSJed Brown
534b5966ddSBarry Smith for (j = 0; j < nr * nc; j++) {
549566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval));
55c4762a1bSJed Brown vals[j] = rval;
56c4762a1bSJed Brown }
574b5966ddSBarry Smith PetscCall(MatSetValues(A, nr, rows, nc, cols, vals, ADD_VALUES));
584b5966ddSBarry Smith PetscCall(MatSetValues(B, nr, rows, nc, cols, vals, ADD_VALUES));
59c4762a1bSJed Brown }
60c4762a1bSJed Brown
619566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
629566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
639566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
649566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
65c4762a1bSJed Brown
665a2b941aSBarry Smith /* Test MatConvert_SeqAIJ_Seq(S)BAIJ handles incompletely filled blocks */
674b5966ddSBarry Smith PetscCall(MatConvert(B, MATBAIJ, MAT_INITIAL_MATRIX, &E));
684b5966ddSBarry Smith PetscCall(MatDestroy(&E));
695a2b941aSBarry Smith PetscCall(MatTranspose(B, MAT_INITIAL_MATRIX, &Bt));
705a2b941aSBarry Smith PetscCall(MatAXPY(Bt, 1.0, B, DIFFERENT_NONZERO_PATTERN));
715a2b941aSBarry Smith PetscCall(MatSetOption(Bt, MAT_SYMMETRIC, PETSC_TRUE));
725a2b941aSBarry Smith PetscCall(MatConvert(Bt, MATSBAIJ, MAT_INITIAL_MATRIX, &E));
735a2b941aSBarry Smith PetscCall(MatDestroy(&E));
745a2b941aSBarry Smith PetscCall(MatDestroy(&Bt));
754b5966ddSBarry Smith
76c4762a1bSJed Brown /* Test MatIncreaseOverlap() */
779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is1));
789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd, &is2));
79c4762a1bSJed Brown
80c4762a1bSJed Brown for (i = 0; i < nd; i++) {
819566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval));
82c4762a1bSJed Brown lsize = (int)(PetscRealPart(rval) * m);
83c4762a1bSJed Brown for (j = 0; j < lsize; j++) {
849566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm, &rval));
85c4762a1bSJed Brown idx[j * bs] = bs * (int)(PetscRealPart(rval) * m);
86c4762a1bSJed Brown for (k = 1; k < bs; k++) idx[j * bs + k] = idx[j * bs] + k;
87c4762a1bSJed Brown }
889566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize * bs, idx, PETSC_COPY_VALUES, is1 + i));
899566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, lsize * bs, idx, PETSC_COPY_VALUES, is2 + i));
90c4762a1bSJed Brown }
919566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A, nd, is1, ov));
929566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(B, nd, is2, ov));
93c4762a1bSJed Brown
94c4762a1bSJed Brown for (i = 0; i < nd; ++i) {
959566063dSJacob Faibussowitsch PetscCall(ISEqual(is1[i], is2[i], &flg));
9628b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "i=%" PetscInt_FMT ", flg =%d", i, (int)flg);
97c4762a1bSJed Brown }
98c4762a1bSJed Brown
99c4762a1bSJed Brown for (i = 0; i < nd; ++i) {
1009566063dSJacob Faibussowitsch PetscCall(ISSort(is1[i]));
1019566063dSJacob Faibussowitsch PetscCall(ISSort(is2[i]));
102c4762a1bSJed Brown }
103c4762a1bSJed Brown
1049566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_INITIAL_MATRIX, &submatA));
1059566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_INITIAL_MATRIX, &submatB));
106c4762a1bSJed Brown
107c4762a1bSJed Brown /* Test MatMult() */
108c4762a1bSJed Brown for (i = 0; i < nd; i++) {
1099566063dSJacob Faibussowitsch PetscCall(MatGetSize(submatA[i], &mm, &nn));
1109566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx));
1119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &s1));
1129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &s2));
113c4762a1bSJed Brown for (j = 0; j < 3; j++) {
1149566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx, rdm));
1159566063dSJacob Faibussowitsch PetscCall(MatMult(submatA[i], xx, s1));
1169566063dSJacob Faibussowitsch PetscCall(MatMult(submatB[i], xx, s2));
1179566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_2, &s1norm));
1189566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_2, &s2norm));
119c4762a1bSJed Brown rnorm = s2norm - s1norm;
12048a46eb9SPierre Jolivet if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", (double)s1norm, (double)s2norm));
121c4762a1bSJed Brown }
1229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx));
1239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1));
1249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2));
125c4762a1bSJed Brown }
126c4762a1bSJed Brown /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */
1279566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_REUSE_MATRIX, &submatA));
1289566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_REUSE_MATRIX, &submatB));
129c4762a1bSJed Brown
130c4762a1bSJed Brown /* Test MatMult() */
131c4762a1bSJed Brown for (i = 0; i < nd; i++) {
1329566063dSJacob Faibussowitsch PetscCall(MatGetSize(submatA[i], &mm, &nn));
1339566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx));
1349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &s1));
1359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &s2));
136c4762a1bSJed Brown for (j = 0; j < 3; j++) {
1379566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx, rdm));
1389566063dSJacob Faibussowitsch PetscCall(MatMult(submatA[i], xx, s1));
1399566063dSJacob Faibussowitsch PetscCall(MatMult(submatB[i], xx, s2));
1409566063dSJacob Faibussowitsch PetscCall(VecNorm(s1, NORM_2, &s1norm));
1419566063dSJacob Faibussowitsch PetscCall(VecNorm(s2, NORM_2, &s2norm));
142c4762a1bSJed Brown rnorm = s2norm - s1norm;
14348a46eb9SPierre Jolivet if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n", (double)s1norm, (double)s2norm));
144c4762a1bSJed Brown }
1459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx));
1469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1));
1479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2));
148c4762a1bSJed Brown }
149c4762a1bSJed Brown
150c4762a1bSJed Brown /* Free allocated memory */
151c4762a1bSJed Brown for (i = 0; i < nd; ++i) {
1529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[i]));
1539566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2[i]));
154c4762a1bSJed Brown }
1559566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd, &submatA));
1569566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd, &submatB));
1579566063dSJacob Faibussowitsch PetscCall(PetscFree(is1));
1589566063dSJacob Faibussowitsch PetscCall(PetscFree(is2));
1599566063dSJacob Faibussowitsch PetscCall(PetscFree(idx));
1609566063dSJacob Faibussowitsch PetscCall(PetscFree(rows));
1619566063dSJacob Faibussowitsch PetscCall(PetscFree(cols));
1629566063dSJacob Faibussowitsch PetscCall(PetscFree(vals));
1639566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
1649566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
1659566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm));
1669566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
167b122ec5aSJacob Faibussowitsch return 0;
168c4762a1bSJed Brown }
169c4762a1bSJed Brown
170c4762a1bSJed Brown /*TEST
171c4762a1bSJed Brown
172c4762a1bSJed Brown test:
173c4762a1bSJed Brown args: -mat_block_size {{1 2 5 7 8}} -ov {{1 3}} -mat_size {{11 13}} -nd {{7}}
174*3886731fSPierre Jolivet output_file: output/empty.out
175c4762a1bSJed Brown
176c4762a1bSJed Brown TEST*/
177