xref: /petsc/src/mat/tests/ex198.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1c4762a1bSJed Brown static char help[] = "Test MatMatMatMult\n\
2c4762a1bSJed Brown Reads PETSc matrix A B and C, then comput D=A*B*C \n\
3c4762a1bSJed Brown Input parameters include\n\
4c4762a1bSJed Brown   -fA <input_file> -fB <input_file> -fC <input_file> \n\n";
5c4762a1bSJed Brown 
6c4762a1bSJed Brown #include <petscmat.h>
7c4762a1bSJed Brown 
8c4762a1bSJed Brown int main(int argc,char **args)
9c4762a1bSJed Brown {
10c4762a1bSJed Brown   Mat            A,B,C,D,BC,ABC;
11c4762a1bSJed Brown   PetscViewer    fd;
12c4762a1bSJed Brown   char           file[3][PETSC_MAX_PATH_LEN];
13c4762a1bSJed Brown   PetscBool      flg;
14c4762a1bSJed Brown   PetscErrorCode ierr;
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
17c4762a1bSJed Brown   /* read matrices A, B and C */
185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-fA",file[0],sizeof(file[0]),&flg));
19*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_USER,"Must indicate binary file with the -fA options");
20c4762a1bSJed Brown 
215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-fB",file[1],sizeof(file[1]),&flg));
22*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_USER,"Must indicate binary file with the -fB options");
23c4762a1bSJed Brown 
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-fC",file[2],sizeof(file[2]),&flg));
25*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_USER,"Must indicate binary file with the -fC options");
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   /* Load matrices */
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[0],FILE_MODE_READ,&fd));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,fd));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
32c4762a1bSJed Brown 
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[1],FILE_MODE_READ,&fd));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(B,fd));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
37c4762a1bSJed Brown 
385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[2],FILE_MODE_READ,&fd));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(C,fd));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /* Test MatMatMult() */
445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(B,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&BC));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(A,BC,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&ABC));
46c4762a1bSJed Brown 
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMatMult(A,B,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&D));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMatMult(A,B,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&D));
495f80ce2aSJacob Faibussowitsch   /* CHKERRQ(MatView(D,PETSC_VIEWER_STDOUT_WORLD)); */
50c4762a1bSJed Brown 
515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatEqual(ABC,D,&flg));
52*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"ABC != D");
53c4762a1bSJed Brown 
545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&ABC));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&BC));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
60c4762a1bSJed Brown   ierr = PetscFinalize();
61c4762a1bSJed Brown   return ierr;
62c4762a1bSJed Brown }
63c4762a1bSJed Brown 
64c4762a1bSJed Brown /*TEST
65c4762a1bSJed Brown 
66c4762a1bSJed Brown    test:
67dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
68c4762a1bSJed Brown       args: -fA ${DATAFILESPATH}/matrices/matmatmatmult/A.bin -fB ${DATAFILESPATH}/matrices/matmatmatmult/B.bin -fC ${DATAFILESPATH}/matrices/matmatmatmult/C.bin
69c4762a1bSJed Brown       output_file: output/ex198.out
70c4762a1bSJed Brown 
71c4762a1bSJed Brown    test:
72c4762a1bSJed Brown       suffix: 2
73c4762a1bSJed Brown       nsize: 3
74dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
75c4762a1bSJed Brown       args: -fA ${DATAFILESPATH}/matrices/matmatmatmult/A.bin -fB ${DATAFILESPATH}/matrices/matmatmatmult/B.bin -fC ${DATAFILESPATH}/matrices/matmatmatmult/C.bin
76c4762a1bSJed Brown       output_file: output/ex198.out
77c4762a1bSJed Brown 
78c4762a1bSJed Brown TEST*/
79