xref: /petsc/src/mat/tests/ex124.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Check the difference of the two matrices \n\
2*c4762a1bSJed Brown Reads PETSc matrix A and B, then check B=A-B \n\
3*c4762a1bSJed Brown Input parameters include\n\
4*c4762a1bSJed Brown   -fA <input_file> -fB <input_file> \n\n";
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown #include <petscmat.h>
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown #undef WRITEFILE
9*c4762a1bSJed Brown PetscInt main(PetscInt argc,char **args)
10*c4762a1bSJed Brown {
11*c4762a1bSJed Brown   Mat            A,B;
12*c4762a1bSJed Brown   PetscViewer    fd;
13*c4762a1bSJed Brown   char           file[2][PETSC_MAX_PATH_LEN];
14*c4762a1bSJed Brown   PetscBool      flg;
15*c4762a1bSJed Brown   PetscErrorCode ierr;
16*c4762a1bSJed Brown   PetscMPIInt    size;
17*c4762a1bSJed Brown   PetscInt       ma,na,mb,nb;
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
20*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
21*c4762a1bSJed Brown   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown   /* read the two matrices, A and B */
24*c4762a1bSJed Brown   ierr = PetscOptionsGetString(NULL,NULL,"-fA",file[0],PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
25*c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -fA options");
26*c4762a1bSJed Brown   ierr = PetscOptionsGetString(NULL,NULL,"-fB",file[1],PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
27*c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -fP options");
28*c4762a1bSJed Brown 
29*c4762a1bSJed Brown   /* Load matrices */
30*c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[0],FILE_MODE_READ,&fd);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = MatLoad(A,fd);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
34*c4762a1bSJed Brown   printf("\n A:\n");
35*c4762a1bSJed Brown   printf("----------------------\n");
36*c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
37*c4762a1bSJed Brown   ierr = MatGetSize(A,&ma,&na);CHKERRQ(ierr);
38*c4762a1bSJed Brown 
39*c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[1],FILE_MODE_READ,&fd);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = MatLoad(B,fd);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
43*c4762a1bSJed Brown   printf("\n B:\n");
44*c4762a1bSJed Brown   printf("----------------------\n");
45*c4762a1bSJed Brown   ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = MatGetSize(B,&mb,&nb);CHKERRQ(ierr);
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown   /* Compute B = -A + B */
49*c4762a1bSJed Brown   if (ma != mb || na != nb) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"nonconforming matrix size");
50*c4762a1bSJed Brown   ierr = MatAXPY(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
51*c4762a1bSJed Brown   printf("\n B - A:\n");
52*c4762a1bSJed Brown   printf("----------------------\n");
53*c4762a1bSJed Brown   ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
56*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = PetscFinalize();
58*c4762a1bSJed Brown   return ierr;
59*c4762a1bSJed Brown }
60