xref: /petsc/src/mat/tests/ex107.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1 
2 static char help[] = "Test MatCreate() with MAT_STRUCTURE_ONLY .\n\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **argv)
7 {
8   Mat            mat;
9   PetscInt       m = 7,n,i,j,rstart,rend;
10   PetscErrorCode ierr;
11   PetscMPIInt    size;
12   PetscScalar    v;
13   PetscBool      struct_only=PETSC_TRUE;
14 
15   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
16   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
17   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
18 
19   ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr);
20   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
21   ierr = PetscOptionsGetBool(NULL,NULL,"-struct_only",&struct_only,NULL);CHKERRQ(ierr);
22   n    = m;
23 
24   /* ------- Assemble matrix, test MatValid() --------- */
25   ierr = MatCreate(PETSC_COMM_WORLD,&mat);CHKERRQ(ierr);
26   ierr = MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
27   ierr = MatSetFromOptions(mat);CHKERRQ(ierr);
28   if (struct_only) {
29     ierr = MatSetOption(mat,MAT_STRUCTURE_ONLY,PETSC_TRUE);CHKERRQ(ierr);
30   }
31   ierr = MatSetUp(mat);CHKERRQ(ierr);
32   ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr);
33   for (i=rstart; i<rend; i++) {
34     for (j=0; j<n; j++) {
35       v    = 10.0*i+j;
36       ierr = MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
37     }
38   }
39   ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
40   ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
41   ierr = MatView(mat,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
42 
43   /* Free data structures */
44   ierr = MatDestroy(&mat);CHKERRQ(ierr);
45   ierr = PetscFinalize();
46   return ierr;
47 }
48 
49 /*TEST
50 
51    test:
52       output_file: output/ex107.out
53 
54    test:
55       suffix: 2
56       args: -mat_type baij -mat_block_size 2 -m 10
57 
58 TEST*/
59