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