xref: /petsc/src/mat/tests/ex165.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests C=A^T*B via MatTranspose() and MatMatMult(). \n\
3c4762a1bSJed Brown                      Contributed by Alexander Grayver, Jan. 2012 \n\n";
4c4762a1bSJed Brown /* Example:
5c4762a1bSJed Brown   mpiexec -n <np> ./ex165 -fA A.dat -fB B.dat -view_C
6c4762a1bSJed Brown  */
7c4762a1bSJed Brown 
8c4762a1bSJed Brown #include <petscmat.h>
9*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
10*d71ae5a4SJacob Faibussowitsch {
11c4762a1bSJed Brown   Mat         A, AT, B, C;
12c4762a1bSJed Brown   PetscViewer viewer;
13c4762a1bSJed Brown   PetscBool   flg;
14c4762a1bSJed Brown   char        file[PETSC_MAX_PATH_LEN];
15c4762a1bSJed Brown 
16327415f7SBarry Smith   PetscFunctionBeginUser;
179566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-fA", file, sizeof(file), &flg));
1928b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Input fileA not specified");
209566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &viewer));
219566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
229566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATAIJ));
239566063dSJacob Faibussowitsch   PetscCall(MatLoad(A, viewer));
249566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
25c4762a1bSJed Brown 
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-fB", file, sizeof(file), &flg));
2728b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Input fileB not specified");
289566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &viewer));
299566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
309566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, MATDENSE));
319566063dSJacob Faibussowitsch   PetscCall(MatLoad(B, viewer));
329566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
33c4762a1bSJed Brown 
349566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &AT));
359566063dSJacob Faibussowitsch   PetscCall(MatMatMult(AT, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C));
36c4762a1bSJed Brown 
379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-view_C", &flg));
38c4762a1bSJed Brown   if (flg) {
399566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "C.dat", FILE_MODE_WRITE, &viewer));
409566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_NATIVE));
419566063dSJacob Faibussowitsch     PetscCall(MatView(C, viewer));
429566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
439566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
44c4762a1bSJed Brown   }
459566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
469566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
479566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&AT));
489566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
499566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
50b122ec5aSJacob Faibussowitsch   return 0;
51c4762a1bSJed Brown }
52