1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests binary I/O of matrices and illustrates user-defined event logging.\n\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscmat.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown /* Note: Most applications would not read and write the same matrix within 7*c4762a1bSJed Brown the same program. This example is intended only to demonstrate 8*c4762a1bSJed Brown both input and output. */ 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown int main(int argc,char **args) 11*c4762a1bSJed Brown { 12*c4762a1bSJed Brown Mat C; 13*c4762a1bSJed Brown PetscScalar v; 14*c4762a1bSJed Brown PetscInt i,j,Ii,J,Istart,Iend,N,m = 4,n = 4; 15*c4762a1bSJed Brown PetscMPIInt rank,size; 16*c4762a1bSJed Brown PetscErrorCode ierr; 17*c4762a1bSJed Brown PetscViewer viewer; 18*c4762a1bSJed Brown #if defined(PETSC_USE_LOG) 19*c4762a1bSJed Brown PetscLogEvent MATRIX_GENERATE,MATRIX_READ; 20*c4762a1bSJed Brown #endif 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 23*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 27*c4762a1bSJed Brown N = m*n; 28*c4762a1bSJed Brown 29*c4762a1bSJed Brown /* PART 1: Generate matrix, then write it in binary format */ 30*c4762a1bSJed Brown 31*c4762a1bSJed Brown ierr = PetscLogEventRegister("Generate Matrix",0,&MATRIX_GENERATE);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = PetscLogEventBegin(MATRIX_GENERATE,0,0,0,0);CHKERRQ(ierr); 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown /* Generate matrix */ 35*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = MatSetFromOptions(C);CHKERRQ(ierr); 38*c4762a1bSJed Brown ierr = MatSetUp(C);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = MatGetOwnershipRange(C,&Istart,&Iend);CHKERRQ(ierr); 40*c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 41*c4762a1bSJed Brown v = -1.0; i = Ii/n; j = Ii - i*n; 42*c4762a1bSJed Brown if (i>0) {J = Ii - n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 43*c4762a1bSJed Brown if (i<m-1) {J = Ii + n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 44*c4762a1bSJed Brown if (j>0) {J = Ii - 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 45*c4762a1bSJed Brown if (j<n-1) {J = Ii + 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 46*c4762a1bSJed Brown v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,ADD_VALUES);CHKERRQ(ierr); 47*c4762a1bSJed Brown } 48*c4762a1bSJed Brown ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 51*c4762a1bSJed Brown 52*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"writing matrix in binary to matrix.dat ...\n");CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = MatView(C,viewer);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = PetscLogEventEnd(MATRIX_GENERATE,0,0,0,0);CHKERRQ(ierr); 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown /* PART 2: Read in matrix in binary format */ 60*c4762a1bSJed Brown 61*c4762a1bSJed Brown /* All processors wait until test matrix has been dumped */ 62*c4762a1bSJed Brown ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRQ(ierr); 63*c4762a1bSJed Brown 64*c4762a1bSJed Brown ierr = PetscLogEventRegister("Read Matrix",0,&MATRIX_READ);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = PetscLogEventBegin(MATRIX_READ,0,0,0,0);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"reading matrix in binary from matrix.dat ...\n");CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&viewer);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = MatLoad(C,viewer);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = PetscLogEventEnd(MATRIX_READ,0,0,0,0);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown /* Free data structures */ 75*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 76*c4762a1bSJed Brown 77*c4762a1bSJed Brown ierr = PetscFinalize(); 78*c4762a1bSJed Brown return ierr; 79*c4762a1bSJed Brown } 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown 82*c4762a1bSJed Brown 83*c4762a1bSJed Brown 84*c4762a1bSJed Brown /*TEST 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown test: 87*c4762a1bSJed Brown filter: grep -v "MPI processes" 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown TEST*/ 90