xref: /petsc/src/mat/tests/ex165.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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>
9c4762a1bSJed Brown int main(int argc,char **args)
10c4762a1bSJed Brown {
11c4762a1bSJed Brown   PetscErrorCode ierr;
12c4762a1bSJed Brown   Mat            A,AT,B,C;
13c4762a1bSJed Brown   PetscViewer    viewer;
14c4762a1bSJed Brown   PetscBool      flg;
15c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN];
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-fA",file,sizeof(file),&flg));
19*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Input fileA not specified");
205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATAIJ));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,viewer));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
25c4762a1bSJed Brown 
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-fB",file,sizeof(file),&flg));
27*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Input fileB not specified");
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(B,MATDENSE));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(B,viewer));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
33c4762a1bSJed Brown 
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&AT));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(AT,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C));
36c4762a1bSJed Brown 
375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-view_C",&flg));
38c4762a1bSJed Brown   if (flg) {
395f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"C.dat",FILE_MODE_WRITE,&viewer));
405f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE));
415f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(C,viewer));
425f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPopFormat(viewer));
435f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&viewer));
44c4762a1bSJed Brown   }
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&AT));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
49c4762a1bSJed Brown   ierr = PetscFinalize();
50c4762a1bSJed Brown   return ierr;
51c4762a1bSJed Brown }
52