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