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 13*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 14*9566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF,m,n,m,n,NULL,&Ad)); 15*9566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_SELF,p,q,p,q,NULL,&Bd)); 16*9566063dSJacob Faibussowitsch PetscCall(MatSetRandom(Ad,NULL)); 17*9566063dSJacob Faibussowitsch PetscCall(MatSetRandom(Bd,NULL)); 18*9566063dSJacob Faibussowitsch PetscCall(MatChop(Ad,0.2)); 19*9566063dSJacob Faibussowitsch PetscCall(MatChop(Bd,0.2)); 20*9566063dSJacob Faibussowitsch PetscCall(MatConvert(Ad,MATAIJ,MAT_INITIAL_MATRIX,&A)); 21*9566063dSJacob Faibussowitsch PetscCall(MatConvert(Bd,MATAIJ,MAT_INITIAL_MATRIX,&B)); 22*9566063dSJacob Faibussowitsch PetscCall(MatSeqAIJKron(A,B,MAT_INITIAL_MATRIX,&C)); 23*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A,NULL,"-A_view")); 24*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(B,NULL,"-B_view")); 25*9566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(C,NULL,"-C_view")); 26*9566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Bd,&Bv)); 27*9566063dSJacob Faibussowitsch PetscCall(MatCreateKAIJ(A,p,q,NULL,Bv,&K)); 28*9566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Bd,&Bv)); 29*9566063dSJacob Faibussowitsch PetscCall(MatMultEqual(C,K,10,&flg)); 3028b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"K*x != C*x"); 31*9566063dSJacob Faibussowitsch PetscCall(MatScale(A,1.3)); 32*9566063dSJacob Faibussowitsch PetscCall(MatScale(B,0.3)); 33*9566063dSJacob Faibussowitsch PetscCall(MatScale(Bd,0.3)); 34*9566063dSJacob Faibussowitsch PetscCall(MatSeqAIJKron(A,B,MAT_REUSE_MATRIX,&C)); 35*9566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(Bd,&Bv)); 36*9566063dSJacob Faibussowitsch PetscCall(MatKAIJSetT(K,p,q,Bv)); 37*9566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(Bd,&Bv)); 38*9566063dSJacob Faibussowitsch PetscCall(MatMultEqual(C,K,10,&flg)); 3928b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"K*x != C*x"); 40*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&K)); 41*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 42*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 43*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 44*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bd)); 45*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ad)); 46*9566063dSJacob 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 55ad7e164aSPierre Jolivet output_file: output/ex101.out 56ad7e164aSPierre Jolivet 57ad7e164aSPierre Jolivet TEST*/ 58