xref: /petsc/src/mat/tests/ex199.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 
2 static char help[] = "Tests the different MatColoring implementatons.\n\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **args)
7 {
8   Mat            C;
9   PetscViewer    viewer;
10   char           file[128];
11   PetscBool      flg;
12   MatColoring    ctx;
13   ISColoring     coloring;
14   PetscMPIInt    size;
15 
16   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
17   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
18 
19   PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg));
20   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must use -f filename to load sparse matrix");
21   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer));
22   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
23   PetscCall(MatLoad(C,viewer));
24   PetscCall(PetscViewerDestroy(&viewer));
25 
26   PetscCall(MatColoringCreate(C,&ctx));
27   PetscCall(MatColoringSetFromOptions(ctx));
28   PetscCall(MatColoringApply(ctx,&coloring));
29   PetscCall(MatColoringTest(ctx,coloring));
30   if (size == 1) {
31     /* jp, power and greedy have bug -- need to be fixed */
32     PetscCall(MatISColoringTest(C,coloring));
33   }
34 
35   /* Free data structures */
36   PetscCall(ISColoringDestroy(&coloring));
37   PetscCall(MatColoringDestroy(&ctx));
38   PetscCall(MatDestroy(&C));
39   PetscCall(PetscFinalize());
40   return 0;
41 }
42 
43 /*TEST
44 
45    test:
46       nsize: {{3}}
47       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
48       args: -f ${DATAFILESPATH}/matrices/arco1 -mat_coloring_type {{ jp power natural greedy}} -mat_coloring_distance {{ 1 2}}
49 
50    test:
51       suffix: 2
52       nsize: {{1 2}}
53       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
54       args: -f ${DATAFILESPATH}/matrices/arco1 -mat_coloring_type {{  sl lf id }} -mat_coloring_distance 2
55       output_file: output/ex199_1.out
56 
57 TEST*/
58