xref: /petsc/src/mat/tests/ex113.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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 */
205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg));
21*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate a file name for matrix A with the -f option.");
22c4762a1bSJed Brown 
235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,viewer));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
27c4762a1bSJed Brown 
285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&A1));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&A2));
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   /* dstMat = A*A1*A2 */
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A1,A2,MAT_INITIAL_MATRIX,fill,&Mtmp));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,Mtmp,MAT_INITIAL_MATRIX,fill,&dstMat));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Mtmp));
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   /* dstMat += A1*A2 */
375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A1,A2,MAT_INITIAL_MATRIX,fill,&Mtmp));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(dstMat,1.0,Mtmp,SUBSET_NONZERO_PATTERN));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Mtmp));
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /* dstMat += A*A1 */
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,A1,MAT_INITIAL_MATRIX,fill,&Mtmp));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(dstMat, 1.0, Mtmp,SUBSET_NONZERO_PATTERN));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Mtmp));
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   /* dstMat += A */
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(dstMat, 1.0, A,SUBSET_NONZERO_PATTERN));
48c4762a1bSJed Brown 
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A1));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A2));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&dstMat));
53c4762a1bSJed Brown   ierr = PetscFinalize();
54c4762a1bSJed Brown   return ierr;
55c4762a1bSJed Brown }
56