xref: /petsc/src/mat/tests/ex43.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Saves a dense matrix in a dense format (binary).\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            C;
9c4762a1bSJed Brown   PetscScalar    v;
10c4762a1bSJed Brown   PetscInt       i,j,m = 4,n = 4;
11c4762a1bSJed Brown   PetscErrorCode ierr;
12c4762a1bSJed Brown   PetscMPIInt    rank,size;
13c4762a1bSJed Brown   PetscViewer    viewer;
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
16*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
17*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
18*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   /* PART 1:  Generate matrix, then write it in binary format */
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /* Generate matrix */
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqDense(PETSC_COMM_WORLD,m,n,NULL,&C));
25c4762a1bSJed Brown   for (i=0; i<m; i++) {
26c4762a1bSJed Brown     for (j=0; j<n; j++) {
27c4762a1bSJed Brown       v    = i*m+j;
28*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
29c4762a1bSJed Brown     }
30c4762a1bSJed Brown   }
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_WRITE,&viewer));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,viewer));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
37c4762a1bSJed Brown   ierr = PetscFinalize();
38c4762a1bSJed Brown   return ierr;
39c4762a1bSJed Brown }
40c4762a1bSJed Brown 
41c4762a1bSJed Brown /*TEST
42c4762a1bSJed Brown 
43c4762a1bSJed Brown    test:
44c4762a1bSJed Brown 
45c4762a1bSJed Brown TEST*/
46