1 static char help[] = "Creates MatSeqBAIJ matrix of given BS for timing tests of MatMult().\n"; 2 3 #include <petscmat.h> 4 5 int main(int argc,char **args) 6 { 7 Mat A; 8 Vec x,y; 9 PetscErrorCode ierr; 10 PetscInt m=50000,bs=12,i,j,k,l,row,col,M, its = 25; 11 PetscScalar rval,*vals; 12 PetscRandom rdm; 13 14 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 15 ierr = PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL);CHKERRQ(ierr); 16 ierr = PetscOptionsGetInt(NULL,NULL,"-its",&its,NULL);CHKERRQ(ierr); 17 ierr = PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL);CHKERRQ(ierr); 18 M = m*bs; 19 ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,27,NULL,&A);CHKERRQ(ierr); 20 ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 21 22 ierr = PetscRandomCreate(PETSC_COMM_SELF,&rdm);CHKERRQ(ierr); 23 ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); 24 ierr = VecCreateSeq(PETSC_COMM_SELF,M,&x);CHKERRQ(ierr); 25 ierr = VecDuplicate(x,&y);CHKERRQ(ierr); 26 27 /* For each block row insert at most 27 blocks */ 28 ierr = PetscMalloc1(bs*bs,&vals);CHKERRQ(ierr); 29 for (i=0; i<m; i++) { 30 row = i; 31 for (j=0; j<27; j++) { 32 ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 33 col = (PetscInt)(PetscRealPart(rval)*m); 34 for (k=0; k<bs; k++) { 35 for (l=0; l<bs; l++) { 36 ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 37 vals[k*bs + l] = rval; 38 } 39 } 40 ierr = MatSetValuesBlocked(A,1,&row,1,&col,vals,INSERT_VALUES);CHKERRQ(ierr); 41 } 42 } 43 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 44 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 45 ierr = PetscFree(vals);CHKERRQ(ierr); 46 47 /* Time MatMult(), MatMultAdd() */ 48 for (i=0; i<its; i++) { 49 ierr = VecSetRandom(x,rdm);CHKERRQ(ierr); 50 ierr = MatMult(A,x,y);CHKERRQ(ierr); 51 ierr = VecSetRandom(x,rdm);CHKERRQ(ierr); 52 ierr = VecSetRandom(y,rdm);CHKERRQ(ierr); 53 ierr = MatMultAdd(A,x,y,y);CHKERRQ(ierr); 54 } 55 56 ierr = MatDestroy(&A);CHKERRQ(ierr); 57 ierr = VecDestroy(&x);CHKERRQ(ierr); 58 ierr = VecDestroy(&y);CHKERRQ(ierr); 59 ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 60 ierr = PetscFinalize(); 61 return ierr; 62 } 63 64 /*TEST 65 66 testset: 67 requires: defined(PETSC_USING_64BIT_PTR) 68 output_file: output/ex238_1.out 69 test: 70 suffix: 1 71 args: -mat_block_size 1 -mat_size 1000 -its 2 72 test: 73 suffix: 2 74 args: -mat_block_size 2 -mat_size 1000 -its 2 75 test: 76 suffix: 4 77 args: -mat_block_size 4 -mat_size 1000 -its 2 78 test: 79 suffix: 5 80 args: -mat_block_size 5 -mat_size 1000 -its 2 81 test: 82 suffix: 6 83 args: -mat_block_size 6 -mat_size 1000 -its 2 84 test: 85 suffix: 8 86 args: -mat_block_size 8 -mat_size 1000 -its 2 87 test: 88 suffix: 12 89 args: -mat_block_size 12 -mat_size 1000 -its 2 90 test: 91 suffix: 15 92 args: -mat_block_size 15 -mat_size 1000 -its 2 93 94 TEST*/ 95