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*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 145f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PETSC_COMM_SELF,m,n,m,n,NULL,&Ad)); 155f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PETSC_COMM_SELF,p,q,p,q,NULL,&Bd)); 165f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(Ad,NULL)); 175f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(Bd,NULL)); 185f80ce2aSJacob Faibussowitsch CHKERRQ(MatChop(Ad,0.2)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(MatChop(Bd,0.2)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(Ad,MATAIJ,MAT_INITIAL_MATRIX,&A)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(Bd,MATAIJ,MAT_INITIAL_MATRIX,&B)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJKron(A,B,MAT_INITIAL_MATRIX,&C)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(A,NULL,"-A_view")); 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(B,NULL,"-B_view")); 255f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(C,NULL,"-C_view")); 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayRead(Bd,&Bv)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateKAIJ(A,p,q,NULL,Bv,&K)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayRead(Bd,&Bv)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(C,K,10,&flg)); 3028b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"K*x != C*x"); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,1.3)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(B,0.3)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(Bd,0.3)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJKron(A,B,MAT_REUSE_MATRIX,&C)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayRead(Bd,&Bv)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatKAIJSetT(K,p,q,Bv)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayRead(Bd,&Bv)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(C,K,10,&flg)); 3928b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"K*x != C*x"); 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&K)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Bd)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Ad)); 46*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 47*b122ec5aSJacob 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