1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests sequential and parallel MatMatMult() and MatAXPY(...,SUBSET_NONZERO_PATTERN) \n\ 3c4762a1bSJed Brown Input arguments are:\n\ 4c4762a1bSJed Brown -f <input_file> : file to load\n\n"; 5c4762a1bSJed Brown /* e.g., mpiexec -n 3 ./ex113 -f <file> */ 6c4762a1bSJed Brown 7c4762a1bSJed Brown #include <petscmat.h> 8c4762a1bSJed Brown 9c4762a1bSJed Brown int main(int argc,char **args) 10c4762a1bSJed Brown { 11c4762a1bSJed Brown Mat A,A1,A2,Mtmp,dstMat; 12c4762a1bSJed Brown PetscViewer viewer; 13c4762a1bSJed Brown PetscErrorCode ierr; 14c4762a1bSJed Brown PetscReal fill=4.0; 15c4762a1bSJed Brown char file[128]; 16c4762a1bSJed Brown PetscBool flg; 17c4762a1bSJed Brown 18c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 19c4762a1bSJed Brown /* Load the matrix A */ 20*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg)); 212c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate a file name for matrix A with the -f option."); 22c4762a1bSJed Brown 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer)); 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,viewer)); 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 27c4762a1bSJed Brown 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&A1)); 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* dstMat = A*A1*A2 */ 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A1,A2,MAT_INITIAL_MATRIX,fill,&Mtmp)); 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A,Mtmp,MAT_INITIAL_MATRIX,fill,&dstMat)); 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Mtmp)); 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* dstMat += A1*A2 */ 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A1,A2,MAT_INITIAL_MATRIX,fill,&Mtmp)); 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(dstMat,1.0,Mtmp,SUBSET_NONZERO_PATTERN)); 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Mtmp)); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* dstMat += A*A1 */ 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A,A1,MAT_INITIAL_MATRIX,fill,&Mtmp)); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(dstMat, 1.0, Mtmp,SUBSET_NONZERO_PATTERN)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Mtmp)); 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* dstMat += A */ 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(dstMat, 1.0, A,SUBSET_NONZERO_PATTERN)); 48c4762a1bSJed Brown 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A1)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A2)); 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&dstMat)); 53c4762a1bSJed Brown ierr = PetscFinalize(); 54c4762a1bSJed Brown return ierr; 55c4762a1bSJed Brown } 56