xref: /petsc/src/mat/tests/ex122.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown static char help[] = "Test MatMatMult() for AIJ and Dense matrices.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc,char **argv)
6c4762a1bSJed Brown {
7c20d7725SJed Brown   Mat            A,B,C;
8c20d7725SJed Brown   PetscInt       M=10,N=5;
9c4762a1bSJed Brown   PetscErrorCode ierr;
10c4762a1bSJed Brown   PetscRandom    r;
11c4762a1bSJed Brown   PetscBool      equal=PETSC_FALSE;
12c4762a1bSJed Brown   PetscReal      fill = 1.0;
13c4762a1bSJed Brown   PetscInt       nza,am,an;
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
16*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
17*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL));
18*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL));
19c4762a1bSJed Brown 
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&r));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(r));
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /* create a aij matrix A */
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,M));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATAIJ));
27c4762a1bSJed Brown   nza  = (PetscInt)(.3*M); /* num of nozeros in each row of A */
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(A,nza,NULL));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(A,nza,NULL,nza,NULL));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(A,r));
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   /* create a dense matrix B */
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(A,&am,&an));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(B,PETSC_DECIDE,am,N,PETSC_DECIDE));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(B,MATDENSE));
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqDenseSetPreallocation(B,NULL));
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIDenseSetPreallocation(B,NULL));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(B,r));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&r));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Test MatMatMult() */
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C));
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C));
47*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual(B,A,C,10,&equal));
482c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != B*A");
49c4762a1bSJed Brown 
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
53c4762a1bSJed Brown   ierr = PetscFinalize();
54c4762a1bSJed Brown   return ierr;
55c4762a1bSJed Brown }
56c4762a1bSJed Brown 
57c4762a1bSJed Brown /*TEST
58c4762a1bSJed Brown 
59c4762a1bSJed Brown    test:
60c4762a1bSJed Brown       output_file: output/ex122.out
61c4762a1bSJed Brown 
62c4762a1bSJed Brown TEST*/
63