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; 15ad7e164aSPierre Jolivet ierr = MatCreateDense(PETSC_COMM_SELF,m,n,m,n,NULL,&Ad);CHKERRQ(ierr); 16ad7e164aSPierre Jolivet ierr = MatCreateDense(PETSC_COMM_SELF,p,q,p,q,NULL,&Bd);CHKERRQ(ierr); 17ad7e164aSPierre Jolivet ierr = MatSetRandom(Ad,NULL);CHKERRQ(ierr); 18ad7e164aSPierre Jolivet ierr = MatSetRandom(Bd,NULL);CHKERRQ(ierr); 19ad7e164aSPierre Jolivet ierr = MatChop(Ad,0.2);CHKERRQ(ierr); 20ad7e164aSPierre Jolivet ierr = MatChop(Bd,0.2);CHKERRQ(ierr); 21ad7e164aSPierre Jolivet ierr = MatConvert(Ad,MATAIJ,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr); 22ad7e164aSPierre Jolivet ierr = MatConvert(Bd,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 23ad7e164aSPierre Jolivet ierr = MatSeqAIJKron(A,B,MAT_INITIAL_MATRIX,&C);CHKERRQ(ierr); 24ad7e164aSPierre Jolivet ierr = MatViewFromOptions(A,NULL,"-A_view");CHKERRQ(ierr); 25ad7e164aSPierre Jolivet ierr = MatViewFromOptions(B,NULL,"-B_view");CHKERRQ(ierr); 26ad7e164aSPierre Jolivet ierr = MatViewFromOptions(C,NULL,"-C_view");CHKERRQ(ierr); 27ad7e164aSPierre Jolivet ierr = MatDenseGetArrayRead(Bd,&Bv);CHKERRQ(ierr); 28ad7e164aSPierre Jolivet ierr = MatCreateKAIJ(A,p,q,NULL,Bv,&K);CHKERRQ(ierr); 29ad7e164aSPierre Jolivet ierr = MatDenseRestoreArrayRead(Bd,&Bv);CHKERRQ(ierr); 30ad7e164aSPierre Jolivet ierr = MatMultEqual(C,K,10,&flg);CHKERRQ(ierr); 31*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"K*x != C*x"); 32ad7e164aSPierre Jolivet ierr = MatScale(A,1.3);CHKERRQ(ierr); 33ad7e164aSPierre Jolivet ierr = MatScale(B,0.3);CHKERRQ(ierr); 34ad7e164aSPierre Jolivet ierr = MatScale(Bd,0.3);CHKERRQ(ierr); 35ad7e164aSPierre Jolivet ierr = MatSeqAIJKron(A,B,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 36ad7e164aSPierre Jolivet ierr = MatDenseGetArrayRead(Bd,&Bv);CHKERRQ(ierr); 37ad7e164aSPierre Jolivet ierr = MatKAIJSetT(K,p,q,Bv);CHKERRQ(ierr); 38ad7e164aSPierre Jolivet ierr = MatDenseRestoreArrayRead(Bd,&Bv);CHKERRQ(ierr); 39ad7e164aSPierre Jolivet ierr = MatMultEqual(C,K,10,&flg);CHKERRQ(ierr); 40*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"K*x != C*x"); 41ad7e164aSPierre Jolivet ierr = MatDestroy(&K);CHKERRQ(ierr); 42ad7e164aSPierre Jolivet ierr = MatDestroy(&C);CHKERRQ(ierr); 43ad7e164aSPierre Jolivet ierr = MatDestroy(&B);CHKERRQ(ierr); 44ad7e164aSPierre Jolivet ierr = MatDestroy(&A);CHKERRQ(ierr); 45ad7e164aSPierre Jolivet ierr = MatDestroy(&Bd);CHKERRQ(ierr); 46ad7e164aSPierre Jolivet ierr = MatDestroy(&Ad);CHKERRQ(ierr); 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