xref: /petsc/src/mat/tests/ex238.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
16679dcc1SBarry Smith static char help[] = "Creates MatSeqBAIJ matrix of given BS for timing tests of MatMult().\n";
26679dcc1SBarry Smith 
36679dcc1SBarry Smith #include <petscmat.h>
46679dcc1SBarry Smith 
56679dcc1SBarry Smith int main(int argc,char **args)
66679dcc1SBarry Smith {
76679dcc1SBarry Smith   Mat            A;
86679dcc1SBarry Smith   Vec            x,y;
96679dcc1SBarry Smith   PetscErrorCode ierr;
10760b957fSBarry Smith   PetscInt       m=50000,bs=12,i,j,k,l,row,col,M, its = 25;
116679dcc1SBarry Smith   PetscScalar    rval,*vals;
126679dcc1SBarry Smith   PetscRandom    rdm;
136679dcc1SBarry Smith 
146679dcc1SBarry Smith   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
156679dcc1SBarry Smith   ierr = PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL);CHKERRQ(ierr);
16760b957fSBarry Smith   ierr = PetscOptionsGetInt(NULL,NULL,"-its",&its,NULL);CHKERRQ(ierr);
176679dcc1SBarry Smith   ierr = PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL);CHKERRQ(ierr);
186679dcc1SBarry Smith   M    = m*bs;
196679dcc1SBarry Smith   ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,27,NULL,&A);CHKERRQ(ierr);
20e269983cSBarry Smith   ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
216679dcc1SBarry Smith 
226679dcc1SBarry Smith   ierr = PetscRandomCreate(PETSC_COMM_SELF,&rdm);CHKERRQ(ierr);
236679dcc1SBarry Smith   ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr);
246679dcc1SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,M,&x);CHKERRQ(ierr);
256679dcc1SBarry Smith   ierr = VecDuplicate(x,&y);CHKERRQ(ierr);
266679dcc1SBarry Smith 
27e269983cSBarry Smith   /* For each block row insert at most 27 blocks */
286679dcc1SBarry Smith   ierr = PetscMalloc1(bs*bs,&vals);CHKERRQ(ierr);
296679dcc1SBarry Smith   for (i=0; i<m; i++) {
306679dcc1SBarry Smith     row = i;
316679dcc1SBarry Smith     for (j=0; j<27; j++) {
326679dcc1SBarry Smith       ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr);
336679dcc1SBarry Smith       col  = (PetscInt)(PetscRealPart(rval)*m);
346679dcc1SBarry Smith       for (k=0; k<bs; k++) {
356679dcc1SBarry Smith         for (l=0; l<bs; l++) {
366679dcc1SBarry Smith           ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr);
376679dcc1SBarry Smith           vals[k*bs + l] = rval;
386679dcc1SBarry Smith         }
396679dcc1SBarry Smith       }
406679dcc1SBarry Smith       ierr = MatSetValuesBlocked(A,1,&row,1,&col,vals,INSERT_VALUES);CHKERRQ(ierr);
416679dcc1SBarry Smith     }
426679dcc1SBarry Smith   }
436679dcc1SBarry Smith   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
446679dcc1SBarry Smith   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
456679dcc1SBarry Smith   ierr = PetscFree(vals);CHKERRQ(ierr);
466679dcc1SBarry Smith 
476679dcc1SBarry Smith   /* Time MatMult(), MatMultAdd() */
48760b957fSBarry Smith   for (i=0; i<its; i++) {
496679dcc1SBarry Smith     ierr  = VecSetRandom(x,rdm);CHKERRQ(ierr);
506679dcc1SBarry Smith     ierr  = MatMult(A,x,y);CHKERRQ(ierr);
516679dcc1SBarry Smith     ierr  = VecSetRandom(x,rdm);CHKERRQ(ierr);
526679dcc1SBarry Smith     ierr  = VecSetRandom(y,rdm);CHKERRQ(ierr);
536679dcc1SBarry Smith     ierr  = MatMultAdd(A,x,y,y);CHKERRQ(ierr);
546679dcc1SBarry Smith   }
556679dcc1SBarry Smith 
566679dcc1SBarry Smith   ierr = MatDestroy(&A);CHKERRQ(ierr);
576679dcc1SBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
586679dcc1SBarry Smith   ierr = VecDestroy(&y);CHKERRQ(ierr);
596679dcc1SBarry Smith   ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr);
606679dcc1SBarry Smith   ierr = PetscFinalize();
616679dcc1SBarry Smith   return ierr;
626679dcc1SBarry Smith }
636679dcc1SBarry Smith 
646679dcc1SBarry Smith /*TEST
656679dcc1SBarry Smith 
66cc384ab3SSatish Balay    testset:
67*dfd57a17SPierre Jolivet      requires: defined(PETSC_USING_64BIT_PTR)
68cc384ab3SSatish Balay      output_file: output/ex238_1.out
69cc384ab3SSatish Balay      test:
70cc384ab3SSatish Balay        suffix: 1
71760b957fSBarry Smith        args: -mat_block_size 1 -mat_size 1000 -its 2
72cc384ab3SSatish Balay      test:
73cc384ab3SSatish Balay        suffix: 2
74760b957fSBarry Smith        args: -mat_block_size 2 -mat_size 1000 -its 2
75cc384ab3SSatish Balay      test:
76cc384ab3SSatish Balay        suffix: 4
77760b957fSBarry Smith        args: -mat_block_size 4 -mat_size 1000 -its 2
78cc384ab3SSatish Balay      test:
79cc384ab3SSatish Balay        suffix: 5
80760b957fSBarry Smith        args: -mat_block_size 5 -mat_size 1000 -its 2
81cc384ab3SSatish Balay      test:
82cc384ab3SSatish Balay        suffix: 6
83760b957fSBarry Smith        args: -mat_block_size 6 -mat_size 1000 -its 2
84cc384ab3SSatish Balay      test:
85cc384ab3SSatish Balay        suffix: 8
86760b957fSBarry Smith        args: -mat_block_size 8 -mat_size 1000 -its 2
87cc384ab3SSatish Balay      test:
88cc384ab3SSatish Balay        suffix: 12
89760b957fSBarry Smith        args: -mat_block_size 12 -mat_size 1000 -its 2
90cc384ab3SSatish Balay      test:
91cc384ab3SSatish Balay        suffix: 15
92760b957fSBarry Smith        args: -mat_block_size 15 -mat_size 1000 -its 2
936679dcc1SBarry Smith 
946679dcc1SBarry Smith TEST*/
95