xref: /petsc/src/mat/tests/ex107.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1c4762a1bSJed Brown static char help[] = "Test MatCreate() with MAT_STRUCTURE_ONLY .\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   Mat         mat;
8c4762a1bSJed Brown   PetscInt    m = 7, n, i, j, rstart, rend;
9c4762a1bSJed Brown   PetscMPIInt size;
10c4762a1bSJed Brown   PetscScalar v;
11c4762a1bSJed Brown   PetscBool   struct_only = PETSC_TRUE;
12c4762a1bSJed Brown 
13327415f7SBarry Smith   PetscFunctionBeginUser;
14*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
16be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
17c4762a1bSJed Brown 
189566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_COMMON));
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-struct_only", &struct_only, NULL));
21c4762a1bSJed Brown   n = m;
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /* ------- Assemble matrix, test MatValid() --------- */
249566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &mat));
259566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat, PETSC_DECIDE, PETSC_DECIDE, m, n));
269566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(mat));
271baa6e33SBarry Smith   if (struct_only) PetscCall(MatSetOption(mat, MAT_STRUCTURE_ONLY, PETSC_TRUE));
289566063dSJacob Faibussowitsch   PetscCall(MatSetUp(mat));
299566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(mat, &rstart, &rend));
30c4762a1bSJed Brown   for (i = rstart; i < rend; i++) {
31c4762a1bSJed Brown     for (j = 0; j < n; j++) {
32c4762a1bSJed Brown       v = 10.0 * i + j;
339566063dSJacob Faibussowitsch       PetscCall(MatSetValues(mat, 1, &i, 1, &j, &v, INSERT_VALUES));
34c4762a1bSJed Brown     }
35c4762a1bSJed Brown   }
369566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
379566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
389566063dSJacob Faibussowitsch   PetscCall(MatView(mat, PETSC_VIEWER_STDOUT_WORLD));
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   /* Free data structures */
419566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
429566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
43b122ec5aSJacob Faibussowitsch   return 0;
44c4762a1bSJed Brown }
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /*TEST
47c4762a1bSJed Brown 
48c4762a1bSJed Brown    test:
49c4762a1bSJed Brown       output_file: output/ex107.out
50c4762a1bSJed Brown 
51c4762a1bSJed Brown    test:
52c4762a1bSJed Brown       suffix: 2
53c4762a1bSJed Brown       args: -mat_type baij -mat_block_size 2 -m 10
54c4762a1bSJed Brown 
55c4762a1bSJed Brown TEST*/
56