xref: /petsc/src/mat/tests/ex113.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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 */
20589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg);CHKERRQ(ierr);
21*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate a file name for matrix A with the -f option.");
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
24c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
25c4762a1bSJed Brown   ierr = MatLoad(A,viewer);CHKERRQ(ierr);
26c4762a1bSJed Brown   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   ierr = MatDuplicate(A,MAT_COPY_VALUES,&A1);CHKERRQ(ierr);
29c4762a1bSJed Brown   ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr);
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   /* dstMat = A*A1*A2 */
32c4762a1bSJed Brown   ierr = MatMatMult(A1,A2,MAT_INITIAL_MATRIX,fill,&Mtmp);CHKERRQ(ierr);
33c4762a1bSJed Brown   ierr = MatMatMult(A,Mtmp,MAT_INITIAL_MATRIX,fill,&dstMat);CHKERRQ(ierr);
34c4762a1bSJed Brown   ierr = MatDestroy(&Mtmp);CHKERRQ(ierr);
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   /* dstMat += A1*A2 */
37c4762a1bSJed Brown   ierr = MatMatMult(A1,A2,MAT_INITIAL_MATRIX,fill,&Mtmp);CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = MatAXPY(dstMat,1.0,Mtmp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
39c4762a1bSJed Brown   ierr = MatDestroy(&Mtmp);CHKERRQ(ierr);
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /* dstMat += A*A1 */
42c4762a1bSJed Brown   ierr = MatMatMult(A,A1,MAT_INITIAL_MATRIX,fill,&Mtmp);CHKERRQ(ierr);
43c4762a1bSJed Brown   ierr = MatAXPY(dstMat, 1.0, Mtmp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
44c4762a1bSJed Brown   ierr = MatDestroy(&Mtmp);CHKERRQ(ierr);
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   /* dstMat += A */
47c4762a1bSJed Brown   ierr = MatAXPY(dstMat, 1.0, A,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
50c4762a1bSJed Brown   ierr = MatDestroy(&A1);CHKERRQ(ierr);
51c4762a1bSJed Brown   ierr = MatDestroy(&A2);CHKERRQ(ierr);
52c4762a1bSJed Brown   ierr = MatDestroy(&dstMat);CHKERRQ(ierr);
53c4762a1bSJed Brown   ierr = PetscFinalize();
54c4762a1bSJed Brown   return ierr;
55c4762a1bSJed Brown }
56c4762a1bSJed Brown 
57