Actual source code: ex11f.F90
1: ! Test MatCreateMPIAdj() with NULL argument 'values'
3: #include <petsc/finclude/petscmat.h>
4: program main
5: use petscmat
6: implicit none
8: Mat :: mesh, dual
9: MatPartitioning :: part
10: IS :: is
11: PetscInt, parameter :: Nvertices = 6, ncells = 2
12: PetscInt :: ii(3), jj(6)
13: PetscMPIInt :: sz, rnk
14: PetscErrorCode :: ierr
16: PetscCallA(PetscInitialize(ierr))
18: PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, sz, ierr))
19: PetscCheckA(sz == 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, 'This example is for exactly two processes')
20: PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rnk, ierr))
21: ii = [0, 3, 6]
22: if (rnk == 0) then
23: jj = [0, 1, 2, 1, 2, 3]
24: else
25: jj = [1, 4, 5, 1, 3, 5]
26: end if
28: PetscCallA(MatCreateMPIAdj(PETSC_COMM_WORLD, ncells, Nvertices, ii, jj, PETSC_NULL_INTEGER_ARRAY, mesh, ierr))
29: PetscCallA(MatMeshToCellGraph(mesh, 2_PETSC_INT_KIND, dual, ierr))
30: PetscCallA(MatView(dual, PETSC_VIEWER_STDOUT_WORLD, ierr))
32: PetscCallA(MatPartitioningCreate(PETSC_COMM_WORLD, part, ierr))
33: PetscCallA(MatPartitioningSetAdjacency(part, dual, ierr))
34: PetscCallA(MatPartitioningSetFromOptions(part, ierr))
35: PetscCallA(MatPartitioningApply(part, is, ierr))
36: PetscCallA(ISView(is, PETSC_VIEWER_STDOUT_WORLD, ierr))
37: PetscCallA(ISDestroy(is, ierr))
38: PetscCallA(MatPartitioningDestroy(part, ierr))
40: PetscCallA(MatDestroy(mesh, ierr))
41: PetscCallA(MatDestroy(dual, ierr))
42: PetscCallA(PetscFinalize(ierr))
44: end program
46: !/*TEST
47: !
48: ! build:
49: ! requires: parmetis
50: !
51: ! test:
52: ! nsize: 2
53: !
54: !TEST*/