xref: /petsc/src/mat/tests/ex31.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests binary I/O of matrices and illustrates user-defined event logging.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown /* Note:  Most applications would not read and write the same matrix within
7c4762a1bSJed Brown   the same program.  This example is intended only to demonstrate
8c4762a1bSJed Brown   both input and output. */
9c4762a1bSJed Brown 
10c4762a1bSJed Brown int main(int argc,char **args)
11c4762a1bSJed Brown {
12c4762a1bSJed Brown   Mat            C;
13c4762a1bSJed Brown   PetscScalar    v;
14c4762a1bSJed Brown   PetscInt       i,j,Ii,J,Istart,Iend,N,m = 4,n = 4;
15c4762a1bSJed Brown   PetscMPIInt    rank,size;
16c4762a1bSJed Brown   PetscViewer    viewer;
17c4762a1bSJed Brown #if defined(PETSC_USE_LOG)
18c4762a1bSJed Brown   PetscLogEvent MATRIX_GENERATE,MATRIX_READ;
19c4762a1bSJed Brown #endif
20c4762a1bSJed Brown 
21*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
225f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
235f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
26c4762a1bSJed Brown   N    = m*n;
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* PART 1:  Generate matrix, then write it in binary format */
29c4762a1bSJed Brown 
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventRegister("Generate Matrix",0,&MATRIX_GENERATE));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(MATRIX_GENERATE,0,0,0,0));
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   /* Generate matrix */
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(C,&Istart,&Iend));
39c4762a1bSJed Brown   for (Ii=Istart; Ii<Iend; Ii++) {
40c4762a1bSJed Brown     v = -1.0; i = Ii/n; j = Ii - i*n;
415f80ce2aSJacob Faibussowitsch     if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES));}
425f80ce2aSJacob Faibussowitsch     if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES));}
435f80ce2aSJacob Faibussowitsch     if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES));}
445f80ce2aSJacob Faibussowitsch     if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES));}
455f80ce2aSJacob Faibussowitsch     v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,ADD_VALUES));
46c4762a1bSJed Brown   }
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
50c4762a1bSJed Brown 
515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"writing matrix in binary to matrix.dat ...\n"));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_WRITE,&viewer));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,viewer));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(MATRIX_GENERATE,0,0,0,0));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /* PART 2:  Read in matrix in binary format */
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* All processors wait until test matrix has been dumped */
615f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD));
62c4762a1bSJed Brown 
635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventRegister("Read Matrix",0,&MATRIX_READ));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(MATRIX_READ,0,0,0,0));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"reading matrix in binary from matrix.dat ...\n"));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&viewer));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(C,viewer));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(MATRIX_READ,0,0,0,0));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* Free data structures */
745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
75c4762a1bSJed Brown 
76*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
77*b122ec5aSJacob Faibussowitsch   return 0;
78c4762a1bSJed Brown }
79c4762a1bSJed Brown 
80c4762a1bSJed Brown /*TEST
81c4762a1bSJed Brown 
82c4762a1bSJed Brown    test:
83c4762a1bSJed Brown       filter: grep -v "MPI processes"
84c4762a1bSJed Brown 
85c4762a1bSJed Brown TEST*/
86