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