xref: /petsc/src/mat/tests/ex48.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1c4762a1bSJed Brown static char help[] = "Tests various routines in MatSeqBAIJ 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 {
7c4762a1bSJed Brown   Mat           A, B, C, D, Fact;
8c4762a1bSJed Brown   Vec           xx, s1, s2, yy;
9c4762a1bSJed Brown   PetscInt      m = 45, rows[2], cols[2], bs = 1, i, row, col, *idx, M;
10c4762a1bSJed Brown   PetscScalar   rval, vals1[4], vals2[4];
11c4762a1bSJed Brown   PetscRandom   rdm;
12c4762a1bSJed Brown   IS            is1, is2;
13c4762a1bSJed Brown   PetscReal     s1norm, s2norm, rnorm, tol = 1.e-4;
14c4762a1bSJed Brown   PetscBool     flg;
15c4762a1bSJed Brown   MatFactorInfo info;
16c4762a1bSJed Brown 
17327415f7SBarry Smith   PetscFunctionBeginUser;
18c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
19c4762a1bSJed Brown   /* Test MatSetValues() and MatGetValues() */
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL));
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &m, NULL));
22c4762a1bSJed Brown   M = m * bs;
239566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF, bs, M, M, 1, NULL, &A));
249566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
259566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, M, M, 15, NULL, &B));
269566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
279566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
289566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
299566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, M, &xx));
309566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(xx, &s1));
319566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(xx, &s2));
329566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(xx, &yy));
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   /* For each row add at least 15 elements */
35c4762a1bSJed Brown   for (row = 0; row < M; row++) {
36c4762a1bSJed Brown     for (i = 0; i < 25 * bs; i++) {
379566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rdm, &rval));
38c4762a1bSJed Brown       col = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M));
399566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &row, 1, &col, &rval, INSERT_VALUES));
409566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B, 1, &row, 1, &col, &rval, INSERT_VALUES));
41c4762a1bSJed Brown     }
42c4762a1bSJed Brown   }
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Now set blocks of values */
45c4762a1bSJed Brown   for (i = 0; i < 20 * bs; i++) {
469566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
47c4762a1bSJed Brown     cols[0]  = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M));
48c4762a1bSJed Brown     vals1[0] = rval;
499566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
50c4762a1bSJed Brown     cols[1]  = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M));
51c4762a1bSJed Brown     vals1[1] = rval;
529566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
53c4762a1bSJed Brown     rows[0]  = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M));
54c4762a1bSJed Brown     vals1[2] = rval;
559566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
56c4762a1bSJed Brown     rows[1]  = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M));
57c4762a1bSJed Brown     vals1[3] = rval;
589566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 2, rows, 2, cols, vals1, INSERT_VALUES));
599566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 2, rows, 2, cols, vals1, INSERT_VALUES));
60c4762a1bSJed Brown   }
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
639566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
649566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
659566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   /* Test MatNorm() */
689566063dSJacob Faibussowitsch   PetscCall(MatNorm(A, NORM_FROBENIUS, &s1norm));
699566063dSJacob Faibussowitsch   PetscCall(MatNorm(B, NORM_FROBENIUS, &s2norm));
70c4762a1bSJed Brown   rnorm = PetscAbsReal(s2norm - s1norm) / s2norm;
7148a46eb9SPierre Jolivet   if (rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatNorm_FROBENIUS()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs));
729566063dSJacob Faibussowitsch   PetscCall(MatNorm(A, NORM_INFINITY, &s1norm));
739566063dSJacob Faibussowitsch   PetscCall(MatNorm(B, NORM_INFINITY, &s2norm));
74c4762a1bSJed Brown   rnorm = PetscAbsReal(s2norm - s1norm) / s2norm;
7548a46eb9SPierre Jolivet   if (rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatNorm_INFINITY()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs));
769566063dSJacob Faibussowitsch   PetscCall(MatNorm(A, NORM_1, &s1norm));
779566063dSJacob Faibussowitsch   PetscCall(MatNorm(B, NORM_1, &s2norm));
78c4762a1bSJed Brown   rnorm = PetscAbsReal(s2norm - s1norm) / s2norm;
7948a46eb9SPierre Jolivet   if (rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatNorm_NORM_1()- NormA=%16.14e NormB=%16.14e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* MatShift() */
82c4762a1bSJed Brown   rval = 10 * s1norm;
839566063dSJacob Faibussowitsch   PetscCall(MatShift(A, rval));
849566063dSJacob Faibussowitsch   PetscCall(MatShift(B, rval));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* Test MatTranspose() */
879566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A, MAT_INPLACE_MATRIX, &A));
889566063dSJacob Faibussowitsch   PetscCall(MatTranspose(B, MAT_INPLACE_MATRIX, &B));
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   /* Now do MatGetValues()  */
91c4762a1bSJed Brown   for (i = 0; i < 30; i++) {
929566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
93c4762a1bSJed Brown     cols[0] = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M));
949566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
95c4762a1bSJed Brown     cols[1] = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M));
969566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
97c4762a1bSJed Brown     rows[0] = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M));
989566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rdm, &rval));
99c4762a1bSJed Brown     rows[1] = PetscMin(M - 1, (PetscInt)(PetscRealPart(rval) * M));
1009566063dSJacob Faibussowitsch     PetscCall(MatGetValues(A, 2, rows, 2, cols, vals1));
1019566063dSJacob Faibussowitsch     PetscCall(MatGetValues(B, 2, rows, 2, cols, vals2));
1029566063dSJacob Faibussowitsch     PetscCall(PetscArraycmp(vals1, vals2, 4, &flg));
10348a46eb9SPierre Jolivet     if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatGetValues bs = %" PetscInt_FMT "\n", bs));
104c4762a1bSJed Brown   }
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* Test MatMult(), MatMultAdd() */
107c4762a1bSJed Brown   for (i = 0; i < 40; i++) {
1089566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(xx, rdm));
1099566063dSJacob Faibussowitsch     PetscCall(VecSet(s2, 0.0));
1109566063dSJacob Faibussowitsch     PetscCall(MatMult(A, xx, s1));
1119566063dSJacob Faibussowitsch     PetscCall(MatMultAdd(A, xx, s2, s2));
1129566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1, NORM_2, &s1norm));
1139566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_2, &s2norm));
114c4762a1bSJed Brown     rnorm = s2norm - s1norm;
11548a46eb9SPierre Jolivet     if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatMult not equal to MatMultAdd Norm1=%e Norm2=%e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs));
116c4762a1bSJed Brown   }
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   /* Test MatMult() */
1199566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(A, B, 10, &flg));
12048a46eb9SPierre Jolivet   if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMult()\n"));
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   /* Test MatMultAdd() */
1239566063dSJacob Faibussowitsch   PetscCall(MatMultAddEqual(A, B, 10, &flg));
12448a46eb9SPierre Jolivet   if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMultAdd()\n"));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* Test MatMultTranspose() */
1279566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeEqual(A, B, 10, &flg));
12848a46eb9SPierre Jolivet   if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMultTranspose()\n"));
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* Test MatMultTransposeAdd() */
1319566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAddEqual(A, B, 10, &flg));
13248a46eb9SPierre Jolivet   if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMultTransposeAdd()\n"));
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   /* Test MatMatMult() */
1359566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, M, 40, NULL, &C));
1369566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(C, rdm));
137fb842aefSJose E. Roman   PetscCall(MatMatMult(A, C, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &D));
1389566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(A, C, D, 40, &flg));
13948a46eb9SPierre Jolivet   if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMatMult()\n"));
1409566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
1419566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, M, 40, NULL, &D));
142fb842aefSJose E. Roman   PetscCall(MatMatMult(A, C, MAT_REUSE_MATRIX, PETSC_DETERMINE, &D));
1439566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(A, C, D, 40, &flg));
14448a46eb9SPierre Jolivet   if (!flg) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMatMult()\n"));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* Do LUFactor() on both the matrices */
1479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M, &idx));
148c4762a1bSJed Brown   for (i = 0; i < M; i++) idx[i] = i;
1499566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, M, idx, PETSC_COPY_VALUES, &is1));
1509566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, M, idx, PETSC_COPY_VALUES, &is2));
1519566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
1529566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(is1));
1539566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(is2));
154c4762a1bSJed Brown 
1559566063dSJacob Faibussowitsch   PetscCall(MatFactorInfoInitialize(&info));
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   info.fill          = 2.0;
158c4762a1bSJed Brown   info.dtcol         = 0.0;
159c4762a1bSJed Brown   info.zeropivot     = 1.e-14;
160c4762a1bSJed Brown   info.pivotinblocks = 1.0;
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   if (bs < 4) {
1639566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(A, "petsc", MAT_FACTOR_LU, &Fact));
1649566063dSJacob Faibussowitsch     PetscCall(MatLUFactorSymbolic(Fact, A, is1, is2, &info));
1659566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(Fact, A, &info));
1669566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(yy, rdm));
1679566063dSJacob Faibussowitsch     PetscCall(MatForwardSolve(Fact, yy, xx));
1689566063dSJacob Faibussowitsch     PetscCall(MatBackwardSolve(Fact, xx, s1));
1699566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Fact));
1709566063dSJacob Faibussowitsch     PetscCall(VecScale(s1, -1.0));
1719566063dSJacob Faibussowitsch     PetscCall(MatMultAdd(A, s1, yy, yy));
1729566063dSJacob Faibussowitsch     PetscCall(VecNorm(yy, NORM_2, &rnorm));
17348a46eb9SPierre Jolivet     if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatForwardSolve/MatBackwardSolve - Norm1=%16.14e bs = %" PetscInt_FMT "\n", (double)rnorm, bs));
174c4762a1bSJed Brown   }
175c4762a1bSJed Brown 
1769566063dSJacob Faibussowitsch   PetscCall(MatLUFactor(B, is1, is2, &info));
1779566063dSJacob Faibussowitsch   PetscCall(MatLUFactor(A, is1, is2, &info));
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   /* Test MatSolveAdd() */
180c4762a1bSJed Brown   for (i = 0; i < 10; i++) {
1819566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(xx, rdm));
1829566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(yy, rdm));
1839566063dSJacob Faibussowitsch     PetscCall(MatSolveAdd(B, xx, yy, s2));
1849566063dSJacob Faibussowitsch     PetscCall(MatSolveAdd(A, xx, yy, s1));
1859566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1, NORM_2, &s1norm));
1869566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_2, &s2norm));
187c4762a1bSJed Brown     rnorm = s2norm - s1norm;
18848a46eb9SPierre Jolivet     if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatSolveAdd - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs));
189c4762a1bSJed Brown   }
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   /* Test MatSolveAdd() when x = A'b +x */
192c4762a1bSJed Brown   for (i = 0; i < 10; i++) {
1939566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(xx, rdm));
1949566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(s1, rdm));
1959566063dSJacob Faibussowitsch     PetscCall(VecCopy(s2, s1));
1969566063dSJacob Faibussowitsch     PetscCall(MatSolveAdd(B, xx, s2, s2));
1979566063dSJacob Faibussowitsch     PetscCall(MatSolveAdd(A, xx, s1, s1));
1989566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1, NORM_2, &s1norm));
1999566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_2, &s2norm));
200c4762a1bSJed Brown     rnorm = s2norm - s1norm;
20148a46eb9SPierre Jolivet     if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatSolveAdd(same) - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs));
202c4762a1bSJed Brown   }
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   /* Test MatSolve() */
205c4762a1bSJed Brown   for (i = 0; i < 10; i++) {
2069566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(xx, rdm));
2079566063dSJacob Faibussowitsch     PetscCall(MatSolve(B, xx, s2));
2089566063dSJacob Faibussowitsch     PetscCall(MatSolve(A, xx, s1));
2099566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1, NORM_2, &s1norm));
2109566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_2, &s2norm));
211c4762a1bSJed Brown     rnorm = s2norm - s1norm;
21248a46eb9SPierre Jolivet     if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatSolve - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs));
213c4762a1bSJed Brown   }
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   /* Test MatSolveTranspose() */
216c4762a1bSJed Brown   if (bs < 8) {
217c4762a1bSJed Brown     for (i = 0; i < 10; i++) {
2189566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xx, rdm));
2199566063dSJacob Faibussowitsch       PetscCall(MatSolveTranspose(B, xx, s2));
2209566063dSJacob Faibussowitsch       PetscCall(MatSolveTranspose(A, xx, s1));
2219566063dSJacob Faibussowitsch       PetscCall(VecNorm(s1, NORM_2, &s1norm));
2229566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2, NORM_2, &s2norm));
223c4762a1bSJed Brown       rnorm = s2norm - s1norm;
22448a46eb9SPierre Jolivet       if (rnorm < -tol || rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatSolveTranspose - Norm1=%16.14e Norm2=%16.14e bs = %" PetscInt_FMT "\n", (double)s1norm, (double)s2norm, bs));
225c4762a1bSJed Brown     }
226c4762a1bSJed Brown   }
227c4762a1bSJed Brown 
2289566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2299566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
2309566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
2319566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
2329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xx));
2339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s1));
2349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s2));
2359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&yy));
2369566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is1));
2379566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is2));
2389566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
2399566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
240b122ec5aSJacob Faibussowitsch   return 0;
241c4762a1bSJed Brown }
242c4762a1bSJed Brown 
243c4762a1bSJed Brown /*TEST
244c4762a1bSJed Brown 
245c4762a1bSJed Brown    test:
246c4762a1bSJed Brown       args: -mat_block_size {{1 2 3 4 5 6 7 8}}
247*3886731fSPierre Jolivet       output_file: output/empty.out
248c4762a1bSJed Brown 
249c4762a1bSJed Brown TEST*/
250