1 2 static char help[] = "Creates MatSeqBAIJ matrix of given BS for timing tests of MatMult().\n"; 3 4 #include <petscmat.h> 5 6 int main(int argc,char **args) 7 { 8 Mat A; 9 Vec x,y; 10 PetscErrorCode ierr; 11 PetscInt m=50000,bs=12,i,j,k,l,row,col,M; 12 PetscScalar rval,*vals; 13 PetscRandom rdm; 14 15 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 16 ierr = PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,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<25; 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 65 /*TEST 66 67 test: 68 requires: defined(PETSC_USING_64BIT_PTR) 69 args: -mat_block_size {{1 2 4 5 6 8 12 15}} 70 71 TEST*/ 72