xref: /petsc/src/dm/tests/ex35.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "MatLoad test for loading matrices that are created by DMCreateMatrix() and\n\
3c4762a1bSJed Brown                       stored in binary via MatView_MPI_DA.MatView_MPI_DA stores the matrix\n\
4c4762a1bSJed Brown                       in natural ordering. Hence MatLoad() has to read the matrix first in\n\
5c4762a1bSJed Brown                       natural ordering and then permute it back to the application ordering.This\n\
6c4762a1bSJed Brown                       example is used for testing the subroutine MatLoad_MPI_DA\n\n";
7c4762a1bSJed Brown 
8c4762a1bSJed Brown #include <petscdm.h>
9c4762a1bSJed Brown #include <petscdmda.h>
10c4762a1bSJed Brown 
11c4762a1bSJed Brown int main(int argc,char **argv)
12c4762a1bSJed Brown {
13c4762a1bSJed Brown   PetscInt       X = 10,Y = 8,Z=8;
14c4762a1bSJed Brown   DM             da;
15c4762a1bSJed Brown   PetscViewer    viewer;
16c4762a1bSJed Brown   Mat            A;
17c4762a1bSJed Brown 
18*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"temp.dat",FILE_MODE_WRITE,&viewer));
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   /* Read options */
225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-X",&X,NULL));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-Y",&Y,NULL));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-Z",&Z,NULL));
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   /* Create distributed array and get vectors */
275f80ce2aSJacob Faibussowitsch   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));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetMatType(da,MATMPIAIJ));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(da,&A));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(A,X));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,viewer));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
36c4762a1bSJed Brown 
375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"temp.dat",FILE_MODE_READ,&viewer));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(da,&A));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,viewer));
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /* Free memory */
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
45*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
46*b122ec5aSJacob Faibussowitsch   return 0;
47c4762a1bSJed Brown }
48