xref: /petsc/src/mat/tests/ex163.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatTransposeMatMult() on MatLoad() matrix \n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            A,C,Bdense,Cdense;
9c4762a1bSJed Brown   PetscViewer    fd;              /* viewer */
10c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN]; /* input file name */
11c4762a1bSJed Brown   PetscBool      flg,viewmats=PETSC_FALSE;
12c4762a1bSJed Brown   PetscMPIInt    rank,size;
13c4762a1bSJed Brown   PetscReal      fill=1.0;
14c4762a1bSJed Brown   PetscInt       m,n,i,j,BN=10,rstart,rend,*rows,*cols;
15c4762a1bSJed Brown   PetscScalar    *Barray,*Carray,rval,*array;
16c4762a1bSJed Brown   Vec            x,y;
17c4762a1bSJed Brown   PetscRandom    rand;
18c4762a1bSJed Brown 
19*327415f7SBarry Smith   PetscFunctionBeginUser;
209566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /* Determine file from which we read the matrix A */
249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg));
2528b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   /* Load matrix A */
289566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
299566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
309566063dSJacob Faibussowitsch   PetscCall(MatLoad(A,fd));
319566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   /* Print (for testing only) */
349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL, "-view_mats", &viewmats));
35c4762a1bSJed Brown   if (viewmats) {
36dd400576SPatrick Sanan     if (rank == 0) printf("A_aij:\n");
379566063dSJacob Faibussowitsch     PetscCall(MatView(A,0));
38c4762a1bSJed Brown   }
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   /* Test MatTransposeMatMult_aij_aij() */
419566063dSJacob Faibussowitsch   PetscCall(MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&C));
42c4762a1bSJed Brown   if (viewmats) {
43dd400576SPatrick Sanan     if (rank == 0) printf("\nC = A_aij^T * A_aij:\n");
449566063dSJacob Faibussowitsch     PetscCall(MatView(C,0));
45c4762a1bSJed Brown   }
469566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
479566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&m,&n));
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   /* create a dense matrix Bdense */
509566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&Bdense));
519566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Bdense,m,PETSC_DECIDE,PETSC_DECIDE,BN));
529566063dSJacob Faibussowitsch   PetscCall(MatSetType(Bdense,MATDENSE));
539566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Bdense));
549566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Bdense));
559566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Bdense,&rstart,&rend));
56c4762a1bSJed Brown 
579566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(m,&rows,BN,&cols,m*BN,&array));
58c4762a1bSJed Brown   for (i=0; i<m; i++) rows[i] = rstart + i;
599566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand));
609566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
61c4762a1bSJed Brown   for (j=0; j<BN; j++) {
62c4762a1bSJed Brown     cols[j] = j;
63c4762a1bSJed Brown     for (i=0; i<m; i++) {
649566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rand,&rval));
65c4762a1bSJed Brown       array[m*j+i] = rval;
66c4762a1bSJed Brown     }
67c4762a1bSJed Brown   }
689566063dSJacob Faibussowitsch   PetscCall(MatSetValues(Bdense,m,rows,BN,cols,array,INSERT_VALUES));
699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Bdense,MAT_FINAL_ASSEMBLY));
709566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Bdense,MAT_FINAL_ASSEMBLY));
719566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
729566063dSJacob Faibussowitsch   PetscCall(PetscFree3(rows,cols,array));
73c4762a1bSJed Brown   if (viewmats) {
74dd400576SPatrick Sanan     if (rank == 0) printf("\nBdense:\n");
759566063dSJacob Faibussowitsch     PetscCall(MatView(Bdense,0));
76c4762a1bSJed Brown   }
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /* Test MatTransposeMatMult_aij_dense() */
799566063dSJacob Faibussowitsch   PetscCall(MatTransposeMatMult(A,Bdense,MAT_INITIAL_MATRIX,fill,&C));
809566063dSJacob Faibussowitsch   PetscCall(MatTransposeMatMult(A,Bdense,MAT_REUSE_MATRIX,fill,&C));
81c4762a1bSJed Brown   if (viewmats) {
82dd400576SPatrick Sanan     if (rank == 0) printf("\nC=A^T*Bdense:\n");
839566063dSJacob Faibussowitsch     PetscCall(MatView(C,0));
84c4762a1bSJed Brown   }
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* Check accuracy */
879566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&Cdense));
889566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Cdense,n,PETSC_DECIDE,PETSC_DECIDE,BN));
899566063dSJacob Faibussowitsch   PetscCall(MatSetType(Cdense,MATDENSE));
909566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Cdense));
919566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Cdense));
929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Cdense,MAT_FINAL_ASSEMBLY));
939566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Cdense,MAT_FINAL_ASSEMBLY));
94c4762a1bSJed Brown 
959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
96c4762a1bSJed Brown   if (size == 1) {
979566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,m,NULL,&x));
989566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF,1,n,NULL,&y));
99c4762a1bSJed Brown   } else {
1009566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,m,PETSC_DECIDE,NULL,&x));
1019566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&y));
102c4762a1bSJed Brown   }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   /* Cdense[:,j] = A^T * Bdense[:,j] */
1059566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(Bdense,&Barray));
1069566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(Cdense,&Carray));
107c4762a1bSJed Brown   for (j=0; j<BN; j++) {
1089566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(x,Barray));
1099566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(y,Carray));
110c4762a1bSJed Brown 
1119566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A,x,y));
112c4762a1bSJed Brown 
1139566063dSJacob Faibussowitsch     PetscCall(VecResetArray(x));
1149566063dSJacob Faibussowitsch     PetscCall(VecResetArray(y));
115c4762a1bSJed Brown     Barray += m;
116c4762a1bSJed Brown     Carray += n;
117c4762a1bSJed Brown   }
1189566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(Bdense,&Barray));
1199566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(Cdense,&Carray));
120c4762a1bSJed Brown   if (viewmats) {
121dd400576SPatrick Sanan     if (rank == 0) printf("\nCdense:\n");
1229566063dSJacob Faibussowitsch     PetscCall(MatView(Cdense,0));
123c4762a1bSJed Brown   }
124c4762a1bSJed Brown 
1259566063dSJacob Faibussowitsch   PetscCall(MatEqual(C,Cdense,&flg));
126c4762a1bSJed Brown   if (!flg) {
127dd400576SPatrick Sanan     if (rank == 0) printf(" C != Cdense\n");
128c4762a1bSJed Brown   }
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* Free data structures */
1319566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1329566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1339566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Bdense));
1349566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Cdense));
1359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1369566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
1379566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
138b122ec5aSJacob Faibussowitsch   return 0;
139c4762a1bSJed Brown }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown /*TEST
142c4762a1bSJed Brown 
143c4762a1bSJed Brown    test:
144dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
145c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small
146c4762a1bSJed Brown       output_file: output/ex163.out
147c4762a1bSJed Brown 
148c4762a1bSJed Brown    test:
149c4762a1bSJed Brown       suffix: 2
150c4762a1bSJed Brown       nsize: 3
151dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
152c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small
153c4762a1bSJed Brown       output_file: output/ex163.out
154c4762a1bSJed Brown 
155c4762a1bSJed Brown TEST*/
156