xref: /petsc/src/mat/tests/ex248.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1ad7e164aSPierre Jolivet 
2ad7e164aSPierre Jolivet static char help[] = "Tests MatSeqAIJKron.\n\n";
3ad7e164aSPierre Jolivet 
4ad7e164aSPierre Jolivet #include <petscmat.h>
5ad7e164aSPierre Jolivet 
6ad7e164aSPierre Jolivet int main(int argc,char **argv)
7ad7e164aSPierre Jolivet {
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   PetscErrorCode    ierr;
13ad7e164aSPierre Jolivet 
14ad7e164aSPierre Jolivet   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
15*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateDense(PETSC_COMM_SELF,m,n,m,n,NULL,&Ad));
16*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateDense(PETSC_COMM_SELF,p,q,p,q,NULL,&Bd));
17*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(Ad,NULL));
18*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(Bd,NULL));
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatChop(Ad,0.2));
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatChop(Bd,0.2));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatConvert(Ad,MATAIJ,MAT_INITIAL_MATRIX,&A));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatConvert(Bd,MATAIJ,MAT_INITIAL_MATRIX,&B));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJKron(A,B,MAT_INITIAL_MATRIX,&C));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(A,NULL,"-A_view"));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(B,NULL,"-B_view"));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(C,NULL,"-C_view"));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseGetArrayRead(Bd,&Bv));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateKAIJ(A,p,q,NULL,Bv,&K));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseRestoreArrayRead(Bd,&Bv));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual(C,K,10,&flg));
312c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"K*x != C*x");
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatScale(A,1.3));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatScale(B,0.3));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatScale(Bd,0.3));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJKron(A,B,MAT_REUSE_MATRIX,&C));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseGetArrayRead(Bd,&Bv));
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatKAIJSetT(K,p,q,Bv));
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseRestoreArrayRead(Bd,&Bv));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual(C,K,10,&flg));
402c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"K*x != C*x");
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&K));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Bd));
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Ad));
47ad7e164aSPierre Jolivet   ierr = PetscFinalize();
48ad7e164aSPierre Jolivet   return ierr;
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