xref: /petsc/src/mat/tests/ex107.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Test MatCreate() with MAT_STRUCTURE_ONLY .\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **argv)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            mat;
9c4762a1bSJed Brown   PetscInt       m = 7,n,i,j,rstart,rend;
10c4762a1bSJed Brown   PetscMPIInt    size;
11c4762a1bSJed Brown   PetscScalar    v;
12c4762a1bSJed Brown   PetscBool      struct_only=PETSC_TRUE;
13c4762a1bSJed Brown 
14*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
155f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
162c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
17c4762a1bSJed Brown 
185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-struct_only",&struct_only,NULL));
21c4762a1bSJed Brown   n    = m;
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /* ------- Assemble matrix, test MatValid() --------- */
245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&mat));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,m,n));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(mat));
27c4762a1bSJed Brown   if (struct_only) {
285f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOption(mat,MAT_STRUCTURE_ONLY,PETSC_TRUE));
29c4762a1bSJed Brown   }
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(mat));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(mat,&rstart,&rend));
32c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
33c4762a1bSJed Brown     for (j=0; j<n; j++) {
34c4762a1bSJed Brown       v    = 10.0*i+j;
355f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES));
36c4762a1bSJed Brown     }
37c4762a1bSJed Brown   }
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(mat,PETSC_VIEWER_STDOUT_WORLD));
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   /* Free data structures */
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&mat));
44*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
45*b122ec5aSJacob Faibussowitsch   return 0;
46c4762a1bSJed Brown }
47c4762a1bSJed Brown 
48c4762a1bSJed Brown /*TEST
49c4762a1bSJed Brown 
50c4762a1bSJed Brown    test:
51c4762a1bSJed Brown       output_file: output/ex107.out
52c4762a1bSJed Brown 
53c4762a1bSJed Brown    test:
54c4762a1bSJed Brown       suffix: 2
55c4762a1bSJed Brown       args: -mat_type baij -mat_block_size 2 -m 10
56c4762a1bSJed Brown 
57c4762a1bSJed Brown TEST*/
58