xref: /petsc/src/mat/tests/ex248.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1 
2 static char help[] = "Tests MatSeqAIJKron.\n\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **argv)
7 {
8   Mat               A,B,C,K,Ad,Bd;
9   const PetscScalar *Bv;
10   PetscInt          n = 10, m = 20, p = 7, q = 17;
11   PetscBool         flg;
12 
13   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
14   CHKERRQ(MatCreateDense(PETSC_COMM_SELF,m,n,m,n,NULL,&Ad));
15   CHKERRQ(MatCreateDense(PETSC_COMM_SELF,p,q,p,q,NULL,&Bd));
16   CHKERRQ(MatSetRandom(Ad,NULL));
17   CHKERRQ(MatSetRandom(Bd,NULL));
18   CHKERRQ(MatChop(Ad,0.2));
19   CHKERRQ(MatChop(Bd,0.2));
20   CHKERRQ(MatConvert(Ad,MATAIJ,MAT_INITIAL_MATRIX,&A));
21   CHKERRQ(MatConvert(Bd,MATAIJ,MAT_INITIAL_MATRIX,&B));
22   CHKERRQ(MatSeqAIJKron(A,B,MAT_INITIAL_MATRIX,&C));
23   CHKERRQ(MatViewFromOptions(A,NULL,"-A_view"));
24   CHKERRQ(MatViewFromOptions(B,NULL,"-B_view"));
25   CHKERRQ(MatViewFromOptions(C,NULL,"-C_view"));
26   CHKERRQ(MatDenseGetArrayRead(Bd,&Bv));
27   CHKERRQ(MatCreateKAIJ(A,p,q,NULL,Bv,&K));
28   CHKERRQ(MatDenseRestoreArrayRead(Bd,&Bv));
29   CHKERRQ(MatMultEqual(C,K,10,&flg));
30   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"K*x != C*x");
31   CHKERRQ(MatScale(A,1.3));
32   CHKERRQ(MatScale(B,0.3));
33   CHKERRQ(MatScale(Bd,0.3));
34   CHKERRQ(MatSeqAIJKron(A,B,MAT_REUSE_MATRIX,&C));
35   CHKERRQ(MatDenseGetArrayRead(Bd,&Bv));
36   CHKERRQ(MatKAIJSetT(K,p,q,Bv));
37   CHKERRQ(MatDenseRestoreArrayRead(Bd,&Bv));
38   CHKERRQ(MatMultEqual(C,K,10,&flg));
39   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"K*x != C*x");
40   CHKERRQ(MatDestroy(&K));
41   CHKERRQ(MatDestroy(&C));
42   CHKERRQ(MatDestroy(&B));
43   CHKERRQ(MatDestroy(&A));
44   CHKERRQ(MatDestroy(&Bd));
45   CHKERRQ(MatDestroy(&Ad));
46   CHKERRQ(PetscFinalize());
47   return 0;
48 }
49 
50 /*TEST
51 
52     test:
53       suffix: 1
54       nsize: 1
55       output_file: output/ex101.out
56 
57 TEST*/
58