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; 18589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-fA",file,sizeof(file),&flg);CHKERRQ(ierr); 19*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Input fileA not specified"); 20c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 21c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 22c4762a1bSJed Brown ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 23c4762a1bSJed Brown ierr = MatLoad(A,viewer);CHKERRQ(ierr); 24c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 25c4762a1bSJed Brown 26589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-fB",file,sizeof(file),&flg);CHKERRQ(ierr); 27*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Input fileB not specified"); 28c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = MatSetType(B,MATDENSE);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = MatLoad(B,viewer);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 33c4762a1bSJed Brown 34c4762a1bSJed Brown ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&AT);CHKERRQ(ierr); 351e1ea65dSPierre Jolivet ierr = MatMatMult(AT,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C);CHKERRQ(ierr); 36c4762a1bSJed Brown 37c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-view_C",&flg);CHKERRQ(ierr); 38c4762a1bSJed Brown if (flg) { 39c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"C.dat",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = MatView(C,viewer);CHKERRQ(ierr); 42c4762a1bSJed Brown ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 43c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 44c4762a1bSJed Brown } 45c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = MatDestroy(&AT);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = PetscFinalize(); 50c4762a1bSJed Brown return ierr; 51c4762a1bSJed Brown } 52