xref: /petsc/src/mat/tests/ex238.c (revision 6679dcc1765ed79a3d9c42bd96f478789854e2bd)
1*6679dcc1SBarry Smith 
2*6679dcc1SBarry Smith static char help[] = "Creates MatSeqBAIJ matrix of given BS for timing tests of MatMult().\n";
3*6679dcc1SBarry Smith 
4*6679dcc1SBarry Smith #include <petscmat.h>
5*6679dcc1SBarry Smith 
6*6679dcc1SBarry Smith int main(int argc,char **args)
7*6679dcc1SBarry Smith {
8*6679dcc1SBarry Smith   Mat            A;
9*6679dcc1SBarry Smith   Vec            x,y;
10*6679dcc1SBarry Smith   PetscErrorCode ierr;
11*6679dcc1SBarry Smith   PetscInt       m=50000,bs=12,i,j,k,l,row,col,M;
12*6679dcc1SBarry Smith   PetscScalar    rval,*vals;
13*6679dcc1SBarry Smith   PetscRandom    rdm;
14*6679dcc1SBarry Smith 
15*6679dcc1SBarry Smith   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
16*6679dcc1SBarry Smith   ierr = PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL);CHKERRQ(ierr);
17*6679dcc1SBarry Smith   ierr = PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL);CHKERRQ(ierr);
18*6679dcc1SBarry Smith   M    = m*bs;
19*6679dcc1SBarry Smith   ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,27,NULL,&A);CHKERRQ(ierr);
20*6679dcc1SBarry Smith   ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
21*6679dcc1SBarry Smith 
22*6679dcc1SBarry Smith   ierr = PetscRandomCreate(PETSC_COMM_SELF,&rdm);CHKERRQ(ierr);
23*6679dcc1SBarry Smith   ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr);
24*6679dcc1SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,M,&x);CHKERRQ(ierr);
25*6679dcc1SBarry Smith   ierr = VecDuplicate(x,&y);CHKERRQ(ierr);
26*6679dcc1SBarry Smith 
27*6679dcc1SBarry Smith   /* For each block row insert atleast 27 elements */
28*6679dcc1SBarry Smith   ierr = PetscMalloc1(bs*bs,&vals);CHKERRQ(ierr);
29*6679dcc1SBarry Smith   for (i=0; i<m; i++) {
30*6679dcc1SBarry Smith     row = i;
31*6679dcc1SBarry Smith     for (j=0; j<27; j++) {
32*6679dcc1SBarry Smith       ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr);
33*6679dcc1SBarry Smith       col  = (PetscInt)(PetscRealPart(rval)*m);
34*6679dcc1SBarry Smith       for (k=0; k<bs; k++) {
35*6679dcc1SBarry Smith         for (l=0; l<bs; l++) {
36*6679dcc1SBarry Smith           ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr);
37*6679dcc1SBarry Smith           vals[k*bs + l] = rval;
38*6679dcc1SBarry Smith         }
39*6679dcc1SBarry Smith       }
40*6679dcc1SBarry Smith       ierr = MatSetValuesBlocked(A,1,&row,1,&col,vals,INSERT_VALUES);CHKERRQ(ierr);
41*6679dcc1SBarry Smith     }
42*6679dcc1SBarry Smith   }
43*6679dcc1SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
44*6679dcc1SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
45*6679dcc1SBarry Smith   ierr = PetscFree(vals);CHKERRQ(ierr);
46*6679dcc1SBarry Smith 
47*6679dcc1SBarry Smith   /* Time MatMult(), MatMultAdd() */
48*6679dcc1SBarry Smith   for (i=0; i<25; i++) {
49*6679dcc1SBarry Smith     ierr  = VecSetRandom(x,rdm);CHKERRQ(ierr);
50*6679dcc1SBarry Smith     ierr  = MatMult(A,x,y);CHKERRQ(ierr);
51*6679dcc1SBarry Smith     ierr  = VecSetRandom(x,rdm);CHKERRQ(ierr);
52*6679dcc1SBarry Smith     ierr  = VecSetRandom(y,rdm);CHKERRQ(ierr);
53*6679dcc1SBarry Smith     ierr  = MatMultAdd(A,x,y,y);CHKERRQ(ierr);
54*6679dcc1SBarry Smith   }
55*6679dcc1SBarry Smith 
56*6679dcc1SBarry Smith   ierr = MatDestroy(&A);CHKERRQ(ierr);
57*6679dcc1SBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
58*6679dcc1SBarry Smith   ierr = VecDestroy(&y);CHKERRQ(ierr);
59*6679dcc1SBarry Smith   ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr);
60*6679dcc1SBarry Smith   ierr = PetscFinalize();
61*6679dcc1SBarry Smith   return ierr;
62*6679dcc1SBarry Smith }
63*6679dcc1SBarry Smith 
64*6679dcc1SBarry Smith 
65*6679dcc1SBarry Smith /*TEST
66*6679dcc1SBarry Smith 
67*6679dcc1SBarry Smith    test:
68*6679dcc1SBarry Smith       args: -mat_block_size {{1 2 4 5 6 8 12 15}}
69*6679dcc1SBarry Smith 
70*6679dcc1SBarry Smith TEST*/
71