xref: /petsc/src/mat/tutorials/ex1.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Reads a PETSc matrix and vector from a file and reorders it.\n\
3c4762a1bSJed Brown   -f0 <input_file> : first file to load (small system)\n\
4c4762a1bSJed Brown   -f1 <input_file> : second file to load (larger system)\n\n";
5c4762a1bSJed Brown 
6c4762a1bSJed Brown /*T
7c4762a1bSJed Brown    Concepts: Mat^ordering a matrix - loading a binary matrix and vector;
8c4762a1bSJed Brown    Concepts: Mat^loading a binary matrix and vector;
9c4762a1bSJed Brown    Concepts: Vectors^loading a binary vector;
10c4762a1bSJed Brown    Concepts: PetscLog^preloading executable
11c4762a1bSJed Brown    Processors: 1
12c4762a1bSJed Brown T*/
13c4762a1bSJed Brown 
14c4762a1bSJed Brown /*
15c4762a1bSJed Brown   Include "petscmat.h" so that we can use matrices.
16c4762a1bSJed Brown   automatically includes:
17c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h    - vectors
18c4762a1bSJed Brown      petscmat.h    - matrices
19c4762a1bSJed Brown      petscis.h     - index sets            petscviewer.h - viewers
20c4762a1bSJed Brown */
21c4762a1bSJed Brown #include <petscmat.h>
22c4762a1bSJed Brown 
23c4762a1bSJed Brown int main(int argc,char **args)
24c4762a1bSJed Brown {
25c4762a1bSJed Brown   Mat             A;                      /* matrix */
26c4762a1bSJed Brown   PetscViewer     fd;                     /* viewer */
27c4762a1bSJed Brown   char            file[PETSC_MAX_PATH_LEN];           /* input file name */
28c4762a1bSJed Brown   IS              isrow,iscol;            /* row and column permutations */
29c4762a1bSJed Brown   MatOrderingType rtype = MATORDERINGRCM;
30c4762a1bSJed Brown   PetscBool       flg;
31c4762a1bSJed Brown 
32*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
33c4762a1bSJed Brown   /*
34c4762a1bSJed Brown      Determine files from which we read the two linear systems
35c4762a1bSJed Brown      (matrix and right-hand-side vector).
36c4762a1bSJed Brown   */
375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg));
3828b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   /* -----------------------------------------------------------
41c4762a1bSJed Brown                   Beginning of loop
42c4762a1bSJed Brown      ----------------------------------------------------------- */
43c4762a1bSJed Brown   /*
44c4762a1bSJed Brown      Loop through the reordering 2 times.
45c4762a1bSJed Brown       - The intention here is to preload and solve a small system;
46c4762a1bSJed Brown         then load another (larger) system and solve it as well.
47c4762a1bSJed Brown         This process preloads the instructions with the smaller
48c4762a1bSJed Brown         system so that more accurate performance monitoring (via
49c4762a1bSJed Brown         -log_view) can be done with the larger one (that actually
50c4762a1bSJed Brown         is the system of interest).
51c4762a1bSJed Brown   */
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* - - - - - - - - - - - New Stage - - - - - - - - - - - - -
54c4762a1bSJed Brown                          Load system i
55c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /*
58c4762a1bSJed Brown      Open binary file.  Note that we use FILE_MODE_READ to indicate
59c4762a1bSJed Brown      reading from this file.
60c4762a1bSJed Brown   */
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /*
64c4762a1bSJed Brown      Load the matrix; then destroy the viewer.
65c4762a1bSJed Brown   */
665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATSEQAIJ));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,fd));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
70c4762a1bSJed Brown 
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOrdering(A,rtype,&isrow,&iscol));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(isrow,PETSC_VIEWER_STDOUT_WORLD));
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /*
75c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
76c4762a1bSJed Brown      are no longer needed.
77c4762a1bSJed Brown   */
785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isrow));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iscol));
81c4762a1bSJed Brown 
82*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
83*b122ec5aSJacob Faibussowitsch   return 0;
84c4762a1bSJed Brown }
85c4762a1bSJed Brown 
86c4762a1bSJed Brown /*TEST
87c4762a1bSJed Brown 
88c4762a1bSJed Brown    test:
89dfd57a17SPierre Jolivet       requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES)
90c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/medium -viewer_binary_skip_info
91c4762a1bSJed Brown 
92c4762a1bSJed Brown TEST*/
93