xref: /petsc/src/mat/tests/ex54.c (revision 609caa7c8c030312b00807b4f015fd827bb80932) !
1c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for parallel AIJ and BAIJ formats.\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          E, A, B, Bt, *submatA, *submatB;
8c4762a1bSJed Brown   PetscInt     bs = 1, m = 11, ov = 1, i, j, k, *rows, *cols, nd = 5, *idx, rstart, rend, sz, mm, nn, M, N, Mbs;
9c4762a1bSJed Brown   PetscMPIInt  size, rank;
10c4762a1bSJed Brown   PetscScalar *vals, rval;
11c4762a1bSJed Brown   IS          *is1, *is2;
12c4762a1bSJed Brown   PetscRandom  rdm;
13c4762a1bSJed Brown   Vec          xx, s1, s2;
14c4762a1bSJed Brown   PetscReal    s1norm, s2norm, rnorm, tol = 100 * PETSC_SMALL;
1582b5ce2aSStefano Zampini   PetscBool    flg, test_nd0 = PETSC_FALSE, emptynd;
16c4762a1bSJed Brown 
17327415f7SBarry Smith   PetscFunctionBeginUser;
18c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
21c4762a1bSJed Brown 
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL));
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &m, NULL));
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL));
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL));
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_nd0", &test_nd0, NULL));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* Create a AIJ matrix A */
299566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
309566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, m * bs, m * bs, PETSC_DECIDE, PETSC_DECIDE));
319566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATAIJ));
324b5966ddSBarry Smith   PetscCall(MatSetBlockSize(A, bs));
339566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(A, PETSC_DEFAULT, NULL));
349566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(A, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL));
359566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
369566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   /* Create a BAIJ matrix B */
399566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
409566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, m * bs, m * bs, PETSC_DECIDE, PETSC_DECIDE));
419566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, MATBAIJ));
429566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(B, bs, PETSC_DEFAULT, NULL));
439566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocation(B, bs, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL));
449566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(B));
459566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
46c4762a1bSJed Brown 
479566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rdm));
489566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
519566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
52c4762a1bSJed Brown   Mbs = M / bs;
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &rows));
559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &cols));
569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs * bs, &vals));
579566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M, &idx));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* Now set blocks of values */
60c4762a1bSJed Brown   for (i = 0; i < 40 * bs; i++) {
614b5966ddSBarry Smith     PetscInt nr = 1, nc = 1;
629566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
63c4762a1bSJed Brown     cols[0] = bs * (int)(PetscRealPart(rval) * Mbs);
649566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
65c4762a1bSJed Brown     rows[0] = rstart + bs * (int)(PetscRealPart(rval) * m);
66c4762a1bSJed Brown     for (j = 1; j < bs; j++) {
674b5966ddSBarry Smith       PetscCall(PetscRandomGetValue(rdm, &rval));
684b5966ddSBarry Smith       if (PetscRealPart(rval) > .5) rows[nr++] = rows[0] + j - 1;
694b5966ddSBarry Smith     }
704b5966ddSBarry Smith     for (j = 1; j < bs; j++) {
714b5966ddSBarry Smith       PetscCall(PetscRandomGetValue(rdm, &rval));
724b5966ddSBarry Smith       if (PetscRealPart(rval) > .5) cols[nc++] = cols[0] + j - 1;
73c4762a1bSJed Brown     }
74c4762a1bSJed Brown 
754b5966ddSBarry Smith     for (j = 0; j < nr * nc; j++) {
769566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rdm, &rval));
77c4762a1bSJed Brown       vals[j] = rval;
78c4762a1bSJed Brown     }
794b5966ddSBarry Smith     PetscCall(MatSetValues(A, nr, rows, nc, cols, vals, ADD_VALUES));
804b5966ddSBarry Smith     PetscCall(MatSetValues(B, nr, rows, nc, cols, vals, ADD_VALUES));
81c4762a1bSJed Brown   }
829566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
839566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
86c4762a1bSJed Brown 
875a2b941aSBarry Smith   /* Test MatConvert_MPIAIJ_MPI(S)BAIJ handles incompletely filled blocks */
884b5966ddSBarry Smith   PetscCall(MatConvert(A, MATBAIJ, MAT_INITIAL_MATRIX, &E));
894b5966ddSBarry Smith   PetscCall(MatDestroy(&E));
905a2b941aSBarry Smith   PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Bt));
915a2b941aSBarry Smith   PetscCall(MatAXPY(Bt, 1.0, B, DIFFERENT_NONZERO_PATTERN));
925a2b941aSBarry Smith   PetscCall(MatSetOption(Bt, MAT_SYMMETRIC, PETSC_TRUE));
935a2b941aSBarry Smith   PetscCall(MatConvert(Bt, MATSBAIJ, MAT_INITIAL_MATRIX, &E));
945a2b941aSBarry Smith   PetscCall(MatDestroy(&E));
955a2b941aSBarry Smith   PetscCall(MatDestroy(&Bt));
964b5966ddSBarry Smith 
97c4762a1bSJed Brown   /* Test MatIncreaseOverlap() */
989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd, &is1));
999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd, &is2));
100c4762a1bSJed Brown 
10182b5ce2aSStefano Zampini   emptynd = PETSC_FALSE;
10282b5ce2aSStefano Zampini   if (rank == 0 && test_nd0) emptynd = PETSC_TRUE; /* test case */
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   for (i = 0; i < nd; i++) {
1059566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
106c4762a1bSJed Brown     sz = (int)(PetscRealPart(rval) * m);
107c4762a1bSJed Brown     for (j = 0; j < sz; j++) {
1089566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rdm, &rval));
109c4762a1bSJed Brown       idx[j * bs] = bs * (int)(PetscRealPart(rval) * Mbs);
110c4762a1bSJed Brown       for (k = 1; k < bs; k++) idx[j * bs + k] = idx[j * bs] + k;
111c4762a1bSJed Brown     }
11282b5ce2aSStefano Zampini     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, emptynd ? 0 : sz * bs, idx, PETSC_COPY_VALUES, is1 + i));
11382b5ce2aSStefano Zampini     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, emptynd ? 0 : sz * bs, idx, PETSC_COPY_VALUES, is2 + i));
114c4762a1bSJed Brown   }
1159566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(A, nd, is1, ov));
1169566063dSJacob Faibussowitsch   PetscCall(MatIncreaseOverlap(B, nd, is2, ov));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   for (i = 0; i < nd; ++i) {
1199566063dSJacob Faibussowitsch     PetscCall(ISEqual(is1[i], is2[i], &flg));
120c4762a1bSJed Brown 
12148a46eb9SPierre 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));
122c4762a1bSJed Brown   }
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   for (i = 0; i < nd; ++i) {
1259566063dSJacob Faibussowitsch     PetscCall(ISSort(is1[i]));
1269566063dSJacob Faibussowitsch     PetscCall(ISSort(is2[i]));
127c4762a1bSJed Brown   }
128c4762a1bSJed Brown 
1299566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_INITIAL_MATRIX, &submatB));
1309566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_INITIAL_MATRIX, &submatA));
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* Test MatMult() */
133c4762a1bSJed Brown   for (i = 0; i < nd; i++) {
1349566063dSJacob Faibussowitsch     PetscCall(MatGetSize(submatA[i], &mm, &nn));
1359566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx));
1369566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s1));
1379566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s2));
138c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
1399566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx, rdm));
1409566063dSJacob Faibussowitsch       PetscCall(MatMult(submatA[i], xx, s1));
1419566063dSJacob Faibussowitsch       PetscCall(MatMult(submatB[i], xx, s2));
1429566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1, NORM_2, &s1norm));
1439566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2, NORM_2, &s2norm));
144c4762a1bSJed Brown       rnorm = s2norm - s1norm;
14548a46eb9SPierre 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));
146c4762a1bSJed Brown     }
1479566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
1489566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
1499566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
150c4762a1bSJed Brown   }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */
1539566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_REUSE_MATRIX, &submatA));
1549566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B, nd, is2, is2, MAT_REUSE_MATRIX, &submatB));
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   /* Test MatMult() */
157c4762a1bSJed Brown   for (i = 0; i < nd; i++) {
1589566063dSJacob Faibussowitsch     PetscCall(MatGetSize(submatA[i], &mm, &nn));
1599566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, mm, &xx));
1609566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s1));
1619566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xx, &s2));
162c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
1639566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx, rdm));
1649566063dSJacob Faibussowitsch       PetscCall(MatMult(submatA[i], xx, s1));
1659566063dSJacob Faibussowitsch       PetscCall(MatMult(submatB[i], xx, s2));
1669566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1, NORM_2, &s1norm));
1679566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2, NORM_2, &s2norm));
168c4762a1bSJed Brown       rnorm = s2norm - s1norm;
16948a46eb9SPierre 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));
170c4762a1bSJed Brown     }
1719566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xx));
1729566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s1));
1739566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&s2));
174c4762a1bSJed Brown   }
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   /* Free allocated memory */
177c4762a1bSJed Brown   for (i = 0; i < nd; ++i) {
1789566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is1[i]));
1799566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is2[i]));
180c4762a1bSJed Brown   }
1819566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd, &submatA));
1829566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nd, &submatB));
183c4762a1bSJed Brown 
1849566063dSJacob Faibussowitsch   PetscCall(PetscFree(is1));
1859566063dSJacob Faibussowitsch   PetscCall(PetscFree(is2));
1869566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
1879566063dSJacob Faibussowitsch   PetscCall(PetscFree(rows));
1889566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols));
1899566063dSJacob Faibussowitsch   PetscCall(PetscFree(vals));
1909566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1929566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
1939566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
194b122ec5aSJacob Faibussowitsch   return 0;
195c4762a1bSJed Brown }
196c4762a1bSJed Brown 
197c4762a1bSJed Brown /*TEST
198c4762a1bSJed Brown 
199c4762a1bSJed Brown    test:
200c4762a1bSJed Brown       nsize: {{1 3}}
20182b5ce2aSStefano Zampini       args: -mat_block_size {{1 3 4 6 8}} -ov {{1 3}} -mat_size {{11 13}} -nd 7
202*3886731fSPierre Jolivet       output_file: output/empty.out
203c4762a1bSJed Brown 
204c4762a1bSJed Brown    test:
205c4762a1bSJed Brown       suffix: 2
206c4762a1bSJed Brown       args: -nd 2 -test_nd0
207*3886731fSPierre Jolivet       output_file: output/empty.out
208c4762a1bSJed Brown 
209c4762a1bSJed Brown    test:
210c4762a1bSJed Brown       suffix: 3
211c4762a1bSJed Brown       nsize: 3
212c4762a1bSJed Brown       args: -nd 2 -test_nd0
213*3886731fSPierre Jolivet       output_file: output/empty.out
214c4762a1bSJed Brown 
215c4762a1bSJed Brown TEST*/
216