xref: /petsc/src/mat/tests/ex238.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
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 
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
76679dcc1SBarry Smith   Mat         A;
86679dcc1SBarry Smith   Vec         x, y;
9760b957fSBarry Smith   PetscInt    m = 50000, bs = 12, i, j, k, l, row, col, M, its = 25;
106679dcc1SBarry Smith   PetscScalar rval, *vals;
116679dcc1SBarry Smith   PetscRandom rdm;
126679dcc1SBarry Smith 
13327415f7SBarry Smith   PetscFunctionBeginUser;
14c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL));
169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-its", &its, NULL));
179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_size", &m, NULL));
186679dcc1SBarry Smith   M = m * bs;
199566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF, bs, M, M, 27, NULL, &A));
209566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
216679dcc1SBarry Smith 
229566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
239566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
249566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, M, &x));
259566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &y));
266679dcc1SBarry Smith 
27e269983cSBarry Smith   /* For each block row insert at most 27 blocks */
289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs * bs, &vals));
296679dcc1SBarry Smith   for (i = 0; i < m; i++) {
306679dcc1SBarry Smith     row = i;
316679dcc1SBarry Smith     for (j = 0; j < 27; j++) {
329566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rdm, &rval));
336679dcc1SBarry Smith       col = (PetscInt)(PetscRealPart(rval) * m);
346679dcc1SBarry Smith       for (k = 0; k < bs; k++) {
356679dcc1SBarry Smith         for (l = 0; l < bs; l++) {
369566063dSJacob Faibussowitsch           PetscCall(PetscRandomGetValue(rdm, &rval));
376679dcc1SBarry Smith           vals[k * bs + l] = rval;
386679dcc1SBarry Smith         }
396679dcc1SBarry Smith       }
409566063dSJacob Faibussowitsch       PetscCall(MatSetValuesBlocked(A, 1, &row, 1, &col, vals, INSERT_VALUES));
416679dcc1SBarry Smith     }
426679dcc1SBarry Smith   }
439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
459566063dSJacob Faibussowitsch   PetscCall(PetscFree(vals));
466679dcc1SBarry Smith 
476679dcc1SBarry Smith   /* Time MatMult(), MatMultAdd() */
48760b957fSBarry Smith   for (i = 0; i < its; i++) {
499566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x, rdm));
509566063dSJacob Faibussowitsch     PetscCall(MatMult(A, x, y));
519566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x, rdm));
529566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(y, rdm));
539566063dSJacob Faibussowitsch     PetscCall(MatMultAdd(A, x, y, y));
546679dcc1SBarry Smith   }
556679dcc1SBarry Smith 
569566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
589566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
599566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
609566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
61b122ec5aSJacob Faibussowitsch   return 0;
626679dcc1SBarry Smith }
636679dcc1SBarry Smith 
646679dcc1SBarry Smith /*TEST
656679dcc1SBarry Smith 
66cc384ab3SSatish Balay    testset:
67dfd57a17SPierre Jolivet      requires: defined(PETSC_USING_64BIT_PTR)
68*3886731fSPierre Jolivet      output_file: output/empty.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