xref: /petsc/src/mat/tests/ex113.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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   PetscReal      fill=4.0;
14c4762a1bSJed Brown   char           file[128];
15c4762a1bSJed Brown   PetscBool      flg;
16c4762a1bSJed Brown 
17*327415f7SBarry Smith   PetscFunctionBeginUser;
189566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
19c4762a1bSJed Brown   /*  Load the matrix A */
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg));
2128b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate a file name for matrix A with the -f option.");
22c4762a1bSJed Brown 
239566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer));
249566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
259566063dSJacob Faibussowitsch   PetscCall(MatLoad(A,viewer));
269566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
27c4762a1bSJed Brown 
289566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A1));
299566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2));
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   /* dstMat = A*A1*A2 */
329566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A1,A2,MAT_INITIAL_MATRIX,fill,&Mtmp));
339566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A,Mtmp,MAT_INITIAL_MATRIX,fill,&dstMat));
349566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Mtmp));
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   /* dstMat += A1*A2 */
379566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A1,A2,MAT_INITIAL_MATRIX,fill,&Mtmp));
389566063dSJacob Faibussowitsch   PetscCall(MatAXPY(dstMat,1.0,Mtmp,SUBSET_NONZERO_PATTERN));
399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Mtmp));
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /* dstMat += A*A1 */
429566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A,A1,MAT_INITIAL_MATRIX,fill,&Mtmp));
439566063dSJacob Faibussowitsch   PetscCall(MatAXPY(dstMat, 1.0, Mtmp,SUBSET_NONZERO_PATTERN));
449566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Mtmp));
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   /* dstMat += A */
479566063dSJacob Faibussowitsch   PetscCall(MatAXPY(dstMat, 1.0, A,SUBSET_NONZERO_PATTERN));
48c4762a1bSJed Brown 
499566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
509566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A1));
519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A2));
529566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&dstMat));
539566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
54b122ec5aSJacob Faibussowitsch   return 0;
55c4762a1bSJed Brown }
56