1 2 static char help[] = "MatLoad test for loading matrices that are created by DMCreateMatrix() and\n\ 3 stored in binary via MatView_MPI_DA.MatView_MPI_DA stores the matrix\n\ 4 in natural ordering. Hence MatLoad() has to read the matrix first in\n\ 5 natural ordering and then permute it back to the application ordering.This\n\ 6 example is used for testing the subroutine MatLoad_MPI_DA\n\n"; 7 8 #include <petscdm.h> 9 #include <petscdmda.h> 10 11 int main(int argc,char **argv) 12 { 13 PetscInt X = 10,Y = 8,Z=8; 14 PetscErrorCode ierr; 15 DM da; 16 PetscViewer viewer; 17 Mat A; 18 19 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 20 CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"temp.dat",FILE_MODE_WRITE,&viewer)); 21 22 /* Read options */ 23 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-X",&X,NULL)); 24 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-Y",&Y,NULL)); 25 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-Z",&Z,NULL)); 26 27 /* Create distributed array and get vectors */ 28 CHKERRQ(DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,X,Y,Z,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,NULL,&da)); 29 CHKERRQ(DMSetFromOptions(da)); 30 CHKERRQ(DMSetUp(da)); 31 CHKERRQ(DMSetMatType(da,MATMPIAIJ)); 32 CHKERRQ(DMCreateMatrix(da,&A)); 33 CHKERRQ(MatShift(A,X)); 34 CHKERRQ(MatView(A,viewer)); 35 CHKERRQ(MatDestroy(&A)); 36 CHKERRQ(PetscViewerDestroy(&viewer)); 37 38 CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"temp.dat",FILE_MODE_READ,&viewer)); 39 CHKERRQ(DMCreateMatrix(da,&A)); 40 CHKERRQ(MatLoad(A,viewer)); 41 42 /* Free memory */ 43 CHKERRQ(MatDestroy(&A)); 44 CHKERRQ(PetscViewerDestroy(&viewer)); 45 CHKERRQ(DMDestroy(&da)); 46 ierr = PetscFinalize(); 47 return ierr; 48 } 49