xref: /petsc/src/mat/tests/ex248.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1ad7e164aSPierre Jolivet static char help[] = "Tests MatSeqAIJKron.\n\n";
2ad7e164aSPierre Jolivet 
3ad7e164aSPierre Jolivet #include <petscmat.h>
4ad7e164aSPierre Jolivet 
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
7ad7e164aSPierre Jolivet   Mat                A, B, C, K, Ad, Bd;
8ad7e164aSPierre Jolivet   const PetscScalar *Bv;
9ad7e164aSPierre Jolivet   PetscInt           n = 10, m = 20, p = 7, q = 17;
10ad7e164aSPierre Jolivet   PetscBool          flg;
11ad7e164aSPierre Jolivet 
12327415f7SBarry Smith   PetscFunctionBeginUser;
13c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
149566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_SELF, m, n, m, n, NULL, &Ad));
159566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_SELF, p, q, p, q, NULL, &Bd));
169566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(Ad, NULL));
179566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(Bd, NULL));
182ce66baaSPierre Jolivet   PetscCall(MatFilter(Ad, 0.2, PETSC_FALSE, PETSC_FALSE));
192ce66baaSPierre Jolivet   PetscCall(MatFilter(Bd, 0.2, PETSC_FALSE, PETSC_FALSE));
209566063dSJacob Faibussowitsch   PetscCall(MatConvert(Ad, MATAIJ, MAT_INITIAL_MATRIX, &A));
219566063dSJacob Faibussowitsch   PetscCall(MatConvert(Bd, MATAIJ, MAT_INITIAL_MATRIX, &B));
229566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJKron(A, B, MAT_INITIAL_MATRIX, &C));
239566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(A, NULL, "-A_view"));
249566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(B, NULL, "-B_view"));
259566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(C, NULL, "-C_view"));
269566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Bd, &Bv));
279566063dSJacob Faibussowitsch   PetscCall(MatCreateKAIJ(A, p, q, NULL, Bv, &K));
289566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Bd, &Bv));
299566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(C, K, 10, &flg));
3028b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "K*x != C*x");
319566063dSJacob Faibussowitsch   PetscCall(MatScale(A, 1.3));
329566063dSJacob Faibussowitsch   PetscCall(MatScale(B, 0.3));
339566063dSJacob Faibussowitsch   PetscCall(MatScale(Bd, 0.3));
349566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJKron(A, B, MAT_REUSE_MATRIX, &C));
359566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(Bd, &Bv));
369566063dSJacob Faibussowitsch   PetscCall(MatKAIJSetT(K, p, q, Bv));
379566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(Bd, &Bv));
389566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(C, K, 10, &flg));
3928b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "K*x != C*x");
409566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&K));
419566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
429566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
439566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
449566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Bd));
459566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Ad));
469566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
47b122ec5aSJacob Faibussowitsch   return 0;
48ad7e164aSPierre Jolivet }
49ad7e164aSPierre Jolivet 
50ad7e164aSPierre Jolivet /*TEST
51ad7e164aSPierre Jolivet 
52ad7e164aSPierre Jolivet     test:
53ad7e164aSPierre Jolivet       suffix: 1
54ad7e164aSPierre Jolivet       nsize: 1
55*3886731fSPierre Jolivet       output_file: output/empty.out
56ad7e164aSPierre Jolivet 
57ad7e164aSPierre Jolivet TEST*/
58