xref: /petsc/src/mat/tests/ex163.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1c4762a1bSJed Brown static char help[] = "Tests MatTransposeMatMult() on MatLoad() matrix \n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   Mat          A, C, Bdense, Cdense;
8c4762a1bSJed Brown   PetscViewer  fd;                       /* viewer */
9c4762a1bSJed Brown   char         file[PETSC_MAX_PATH_LEN]; /* input file name */
10c4762a1bSJed Brown   PetscBool    flg, viewmats = PETSC_FALSE;
11c4762a1bSJed Brown   PetscMPIInt  rank, size;
12c4762a1bSJed Brown   PetscReal    fill = 1.0;
13c4762a1bSJed Brown   PetscInt     m, n, i, j, BN = 10, rstart, rend, *rows, *cols;
14c4762a1bSJed Brown   PetscScalar *Barray, *Carray, rval, *array;
15c4762a1bSJed Brown   Vec          x, y;
16c4762a1bSJed Brown   PetscRandom  rand;
17c4762a1bSJed Brown 
18327415f7SBarry Smith   PetscFunctionBeginUser;
19c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   /* Determine file from which we read the matrix A */
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &flg));
2428b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option");
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   /* Load matrix A */
279566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
289566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
299566063dSJacob Faibussowitsch   PetscCall(MatLoad(A, fd));
309566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   /* Print (for testing only) */
339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-view_mats", &viewmats));
34c4762a1bSJed Brown   if (viewmats) {
35dd400576SPatrick Sanan     if (rank == 0) printf("A_aij:\n");
369566063dSJacob Faibussowitsch     PetscCall(MatView(A, 0));
37c4762a1bSJed Brown   }
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   /* Test MatTransposeMatMult_aij_aij() */
409566063dSJacob Faibussowitsch   PetscCall(MatTransposeMatMult(A, A, MAT_INITIAL_MATRIX, fill, &C));
41c4762a1bSJed Brown   if (viewmats) {
42dd400576SPatrick Sanan     if (rank == 0) printf("\nC = A_aij^T * A_aij:\n");
439566063dSJacob Faibussowitsch     PetscCall(MatView(C, 0));
44c4762a1bSJed Brown   }
459566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
469566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
47c4762a1bSJed Brown 
48c4762a1bSJed Brown   /* create a dense matrix Bdense */
499566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &Bdense));
509566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Bdense, m, PETSC_DECIDE, PETSC_DECIDE, BN));
519566063dSJacob Faibussowitsch   PetscCall(MatSetType(Bdense, MATDENSE));
529566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Bdense));
539566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Bdense));
549566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(Bdense, &rstart, &rend));
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(m, &rows, BN, &cols, m * BN, &array));
57c4762a1bSJed Brown   for (i = 0; i < m; i++) rows[i] = rstart + i;
589566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
599566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
60c4762a1bSJed Brown   for (j = 0; j < BN; j++) {
61c4762a1bSJed Brown     cols[j] = j;
62c4762a1bSJed Brown     for (i = 0; i < m; i++) {
639566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rand, &rval));
64c4762a1bSJed Brown       array[m * j + i] = rval;
65c4762a1bSJed Brown     }
66c4762a1bSJed Brown   }
679566063dSJacob Faibussowitsch   PetscCall(MatSetValues(Bdense, m, rows, BN, cols, array, INSERT_VALUES));
689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Bdense, MAT_FINAL_ASSEMBLY));
699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Bdense, MAT_FINAL_ASSEMBLY));
709566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
719566063dSJacob Faibussowitsch   PetscCall(PetscFree3(rows, cols, array));
72c4762a1bSJed Brown   if (viewmats) {
73dd400576SPatrick Sanan     if (rank == 0) printf("\nBdense:\n");
749566063dSJacob Faibussowitsch     PetscCall(MatView(Bdense, 0));
75c4762a1bSJed Brown   }
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /* Test MatTransposeMatMult_aij_dense() */
789566063dSJacob Faibussowitsch   PetscCall(MatTransposeMatMult(A, Bdense, MAT_INITIAL_MATRIX, fill, &C));
799566063dSJacob Faibussowitsch   PetscCall(MatTransposeMatMult(A, Bdense, MAT_REUSE_MATRIX, fill, &C));
80c4762a1bSJed Brown   if (viewmats) {
81dd400576SPatrick Sanan     if (rank == 0) printf("\nC=A^T*Bdense:\n");
829566063dSJacob Faibussowitsch     PetscCall(MatView(C, 0));
83c4762a1bSJed Brown   }
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /* Check accuracy */
869566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &Cdense));
879566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Cdense, n, PETSC_DECIDE, PETSC_DECIDE, BN));
889566063dSJacob Faibussowitsch   PetscCall(MatSetType(Cdense, MATDENSE));
899566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Cdense));
909566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Cdense));
919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Cdense, MAT_FINAL_ASSEMBLY));
929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Cdense, MAT_FINAL_ASSEMBLY));
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
95c4762a1bSJed Brown   if (size == 1) {
969566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &x));
979566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, n, NULL, &y));
98c4762a1bSJed Brown   } else {
999566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, m, PETSC_DECIDE, NULL, &x));
1009566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, n, PETSC_DECIDE, NULL, &y));
101c4762a1bSJed Brown   }
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   /* Cdense[:,j] = A^T * Bdense[:,j] */
1049566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(Bdense, &Barray));
1059566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(Cdense, &Carray));
106c4762a1bSJed Brown   for (j = 0; j < BN; j++) {
1079566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(x, Barray));
1089566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(y, Carray));
109c4762a1bSJed Brown 
1109566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A, x, y));
111c4762a1bSJed Brown 
1129566063dSJacob Faibussowitsch     PetscCall(VecResetArray(x));
1139566063dSJacob Faibussowitsch     PetscCall(VecResetArray(y));
114c4762a1bSJed Brown     Barray += m;
115c4762a1bSJed Brown     Carray += n;
116c4762a1bSJed Brown   }
1179566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(Bdense, &Barray));
1189566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(Cdense, &Carray));
119c4762a1bSJed Brown   if (viewmats) {
120dd400576SPatrick Sanan     if (rank == 0) printf("\nCdense:\n");
1219566063dSJacob Faibussowitsch     PetscCall(MatView(Cdense, 0));
122c4762a1bSJed Brown   }
123c4762a1bSJed Brown 
1249566063dSJacob Faibussowitsch   PetscCall(MatEqual(C, Cdense, &flg));
125c4762a1bSJed Brown   if (!flg) {
126dd400576SPatrick Sanan     if (rank == 0) printf(" C != Cdense\n");
127c4762a1bSJed Brown   }
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* Free data structures */
1309566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1319566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1329566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Bdense));
1339566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Cdense));
1349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
1369566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
137b122ec5aSJacob Faibussowitsch   return 0;
138c4762a1bSJed Brown }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown /*TEST
141c4762a1bSJed Brown 
142c4762a1bSJed Brown    test:
143dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
144c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small
145*3886731fSPierre Jolivet       output_file: output/empty.out
146c4762a1bSJed Brown 
147c4762a1bSJed Brown    test:
148c4762a1bSJed Brown       suffix: 2
149c4762a1bSJed Brown       nsize: 3
150dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
151c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small
152*3886731fSPierre Jolivet       output_file: output/empty.out
153c4762a1bSJed Brown 
154c4762a1bSJed Brown TEST*/
155