xref: /petsc/src/mat/tests/ex122.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1c4762a1bSJed Brown static char help[] = "Test MatMatMult() for AIJ and Dense matrices.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
7c20d7725SJed Brown   Mat         A, B, C;
8c20d7725SJed Brown   PetscInt    M = 10, N = 5;
9c4762a1bSJed Brown   PetscRandom r;
10c4762a1bSJed Brown   PetscBool   equal = PETSC_FALSE;
11c4762a1bSJed Brown   PetscReal   fill  = 1.0;
12c4762a1bSJed Brown   PetscInt    nza, am, an;
13c4762a1bSJed Brown 
14327415f7SBarry Smith   PetscFunctionBeginUser;
15c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL));
179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL));
189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-fill", &fill, NULL));
19c4762a1bSJed Brown 
209566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &r));
219566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(r));
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /* create a aij matrix A */
249566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
259566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, N, M));
269566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATAIJ));
27bbea24aaSStefano Zampini   nza = (PetscInt)(.3 * M); /* num of nonzeros in each row of A */
289566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(A, nza, NULL));
299566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(A, nza, NULL, nza, NULL));
309566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(A, r));
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   /* create a dense matrix B */
339566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &am, &an));
349566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
359566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, PETSC_DECIDE, am, N, PETSC_DECIDE));
369566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, MATDENSE));
379566063dSJacob Faibussowitsch   PetscCall(MatSeqDenseSetPreallocation(B, NULL));
389566063dSJacob Faibussowitsch   PetscCall(MatMPIDenseSetPreallocation(B, NULL));
399566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(B, r));
409566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&r));
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   /* Test MatMatMult() */
439566063dSJacob Faibussowitsch   PetscCall(MatMatMult(B, A, MAT_INITIAL_MATRIX, fill, &C));
449566063dSJacob Faibussowitsch   PetscCall(MatMatMult(B, A, MAT_REUSE_MATRIX, fill, &C));
459566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(B, A, C, 10, &equal));
4628b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "C != B*A");
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
499566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
509566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
519566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
52b122ec5aSJacob Faibussowitsch   return 0;
53c4762a1bSJed Brown }
54c4762a1bSJed Brown 
55c4762a1bSJed Brown /*TEST
56c4762a1bSJed Brown 
57c4762a1bSJed Brown    test:
58*3886731fSPierre Jolivet       output_file: output/empty.out
59c4762a1bSJed Brown 
60c4762a1bSJed Brown TEST*/
61