xref: /petsc/src/mat/tests/ex206.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Reads binary matrix - twice\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown int main(int argc,char **args)
5c4762a1bSJed Brown {
6c4762a1bSJed Brown   Mat               A;
7c4762a1bSJed Brown   PetscViewer       fd;                        /* viewer */
8c4762a1bSJed Brown   char              file[PETSC_MAX_PATH_LEN];  /* input file name */
9c4762a1bSJed Brown   PetscBool         flg;
10c4762a1bSJed Brown 
11c4762a1bSJed Brown   PetscInitialize(&argc,&args,(char*)0,help);
12c4762a1bSJed Brown 
135f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg));
1428b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
15c4762a1bSJed Brown 
165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "First MatLoad! \n"));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,fd));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
23c4762a1bSJed Brown 
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f2",file,sizeof(file),&flg));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Second MatLoad! \n"));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,fd));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&fd));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
30c4762a1bSJed Brown 
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
32*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
33*b122ec5aSJacob Faibussowitsch   return 0;
34c4762a1bSJed Brown }
35c4762a1bSJed Brown 
36c4762a1bSJed Brown /*TEST
37c4762a1bSJed Brown 
38c4762a1bSJed Brown    test:
39c4762a1bSJed Brown       suffix: aij_1
40dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
41c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_type aij -mat_block_size 1
42c4762a1bSJed Brown       filter: grep -v Mat_
43c4762a1bSJed Brown 
44c4762a1bSJed Brown    test:
45c4762a1bSJed Brown       suffix: aij_2
46c4762a1bSJed Brown       nsize: 2
47dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
48c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_type aij -mat_block_size 1
49c4762a1bSJed Brown       filter: grep -v Mat_
50c4762a1bSJed Brown 
51c4762a1bSJed Brown    test:
52c4762a1bSJed Brown       suffix: aij_2_d
53c4762a1bSJed Brown       nsize: 2
54dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
55c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -f2 ${DATAFILESPATH}/matrices/smallbs2 -mat_type aij -mat_block_size 1
56c4762a1bSJed Brown       filter: grep -v Mat_
57c4762a1bSJed Brown 
58c4762a1bSJed Brown    test:
59c4762a1bSJed Brown       suffix: baij_1_2
60dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
61c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_type baij -mat_block_size 2
62c4762a1bSJed Brown       filter: grep -v Mat_
63c4762a1bSJed Brown 
64c4762a1bSJed Brown    test:
65c4762a1bSJed Brown       suffix: baij_2_1_d
66c4762a1bSJed Brown       nsize: 2
67dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
68c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -f2 ${DATAFILESPATH}/matrices/smallbs2 -mat_type baij -mat_block_size 1
69c4762a1bSJed Brown       filter: grep -v Mat_
70c4762a1bSJed Brown 
71c4762a1bSJed Brown    test:
72c4762a1bSJed Brown       suffix: baij_2_2
73c4762a1bSJed Brown       nsize: 2
74dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
75c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_type baij -mat_block_size 2
76c4762a1bSJed Brown       filter: grep -v Mat_
77c4762a1bSJed Brown 
78c4762a1bSJed Brown    test:
79c4762a1bSJed Brown       suffix: baij_2_2_d
80c4762a1bSJed Brown       nsize: 2
81dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
82c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -f2 ${DATAFILESPATH}/matrices/smallbs2 -mat_type baij -mat_block_size 2
83c4762a1bSJed Brown       filter: grep -v Mat_
84c4762a1bSJed Brown 
85c4762a1bSJed Brown    test:
86c4762a1bSJed Brown       suffix: sbaij_1_1
87dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
88c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_type sbaij -mat_block_size 1
89c4762a1bSJed Brown       filter: grep -v Mat_
90c4762a1bSJed Brown 
91c4762a1bSJed Brown    test:
92c4762a1bSJed Brown       suffix: sbaij_1_2
93dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
94c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_type sbaij -mat_block_size 2
95c4762a1bSJed Brown       filter: grep -v Mat_
96c4762a1bSJed Brown 
97c4762a1bSJed Brown    test:
98c4762a1bSJed Brown       suffix: sbaij_2_1_d
99c4762a1bSJed Brown       nsize: 2
100dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
101c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -f2 ${DATAFILESPATH}/matrices/smallbs2 -mat_type sbaij -mat_block_size 1
102c4762a1bSJed Brown       filter: grep -v Mat_
103c4762a1bSJed Brown 
104c4762a1bSJed Brown    test:
105c4762a1bSJed Brown       suffix: sbaij_2_2
106c4762a1bSJed Brown       nsize: 2
107dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
108c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/small -mat_type sbaij -mat_block_size 2
109c4762a1bSJed Brown       filter: grep -v Mat_
110c4762a1bSJed Brown 
111c4762a1bSJed Brown TEST*/
112