xref: /petsc/src/mat/tests/ex124.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1c4762a1bSJed Brown static char help[] = "Check the difference of the two matrices \n\
2c4762a1bSJed Brown Reads PETSc matrix A and B, then check B=A-B \n\
3c4762a1bSJed Brown Input parameters include\n\
4c4762a1bSJed Brown   -fA <input_file> -fB <input_file> \n\n";
5c4762a1bSJed Brown 
6c4762a1bSJed Brown #include <petscmat.h>
7c4762a1bSJed Brown 
8ad0e0ea3SStefano Zampini int main(int argc,char **args)
9c4762a1bSJed Brown {
10c4762a1bSJed Brown   Mat            A,B;
11c4762a1bSJed Brown   PetscViewer    fd;
12c4762a1bSJed Brown   char           file[2][PETSC_MAX_PATH_LEN];
13c4762a1bSJed Brown   PetscBool      flg;
14c4762a1bSJed Brown   PetscErrorCode ierr;
15c4762a1bSJed Brown   PetscMPIInt    size;
16c4762a1bSJed Brown   PetscInt       ma,na,mb,nb;
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
195f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
202c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   /* read the two matrices, A and B */
235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-fA",file[0],sizeof(file[0]),&flg));
24*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -fA options");
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-fB",file[1],sizeof(file[1]),&flg));
26*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -fP options");
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* Load matrices */
295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[0],FILE_MODE_READ,&fd));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,fd));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
33c4762a1bSJed Brown   printf("\n A:\n");
34c4762a1bSJed Brown   printf("----------------------\n");
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(A,&ma,&na));
37c4762a1bSJed Brown 
385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[1],FILE_MODE_READ,&fd));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(B,fd));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
42c4762a1bSJed Brown   printf("\n B:\n");
43c4762a1bSJed Brown   printf("----------------------\n");
445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(B,&mb,&nb));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* Compute B = -A + B */
482c71b3e2SJacob Faibussowitsch   PetscCheckFalse(ma != mb || na != nb,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"nonconforming matrix size");
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(B,-1.0,A,DIFFERENT_NONZERO_PATTERN));
50c4762a1bSJed Brown   printf("\n B - A:\n");
51c4762a1bSJed Brown   printf("----------------------\n");
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
53c4762a1bSJed Brown 
545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
56c4762a1bSJed Brown   ierr = PetscFinalize();
57c4762a1bSJed Brown   return ierr;
58c4762a1bSJed Brown }
59