xref: /petsc/src/mat/tutorials/ex15f.F90 (revision 9b88ac225e01f016352a5f4cd90e158abe5f5675)
1c4762a1bSJed Brown#include <petsc/finclude/petscmat.h>
2*c5e229c2SMartin Diehlprogram main
3c4762a1bSJed Brown  use petscmat
4c4762a1bSJed Brown
5c4762a1bSJed Brown  implicit none
6c4762a1bSJed Brown
7c4762a1bSJed Brown  Mat             :: A
8c4762a1bSJed Brown  MatPartitioning :: part
9c4762a1bSJed Brown  IS              :: is
10c4762a1bSJed Brown  PetscInt        :: r, myStart, myEnd
11c4762a1bSJed Brown  PetscInt        :: N = 10
12c4762a1bSJed Brown  PetscErrorCode  :: ierr
13c4762a1bSJed Brown  PetscScalar, pointer, dimension(:) :: vals
14c4762a1bSJed Brown  PetscInt, pointer, dimension(:) :: cols
15c4762a1bSJed Brown  PetscBool :: flg
16c4762a1bSJed Brown  PetscInt, parameter :: one = 1, two = 2, three = 3
17c4762a1bSJed Brown
18d8606c27SBarry Smith  PetscCallA(PetscInitialize(ierr))
19c4762a1bSJed Brown
20dcb3e689SBarry Smith  PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-N', N, flg, ierr))
21d8606c27SBarry Smith  PetscCallA(MatCreate(PETSC_COMM_WORLD, A, ierr))
22d8606c27SBarry Smith  PetscCallA(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, N, N, ierr))
23d8606c27SBarry Smith  PetscCallA(MatSetFromOptions(A, ierr))
245d83a8b1SBarry Smith  PetscCallA(MatSeqAIJSetPreallocation(A, three, PETSC_NULL_INTEGER_ARRAY, ierr))
255d83a8b1SBarry Smith  PetscCallA(MatMPIAIJSetPreallocation(A, three, PETSC_NULL_INTEGER_ARRAY, two, PETSC_NULL_INTEGER_ARRAY, ierr))
26c4762a1bSJed Brown
27c4762a1bSJed Brown  !/* Create a linear mesh */
28d8606c27SBarry Smith  PetscCallA(MatGetOwnershipRange(A, myStart, myEnd, ierr))
29c4762a1bSJed Brown
30c4762a1bSJed Brown  do r = myStart, myEnd - 1
31c4762a1bSJed Brown    if (r == 0) then
32c4762a1bSJed Brown      allocate (vals(2))
33c4762a1bSJed Brown      vals = 1.0
34c4762a1bSJed Brown      allocate (cols(2), source=[r, r + 1])
355d83a8b1SBarry Smith      PetscCallA(MatSetValues(A, one, [r], two, cols, vals, INSERT_VALUES, ierr))
36c4762a1bSJed Brown      deallocate (cols)
37c4762a1bSJed Brown      deallocate (vals)
38c4762a1bSJed Brown    else if (r == N - 1) then
39c4762a1bSJed Brown      allocate (vals(2))
40c4762a1bSJed Brown      vals = 1.0
41c4762a1bSJed Brown      allocate (cols(2), source=[r - 1, r])
425d83a8b1SBarry Smith      PetscCallA(MatSetValues(A, one, [r], two, cols, vals, INSERT_VALUES, ierr))
43c4762a1bSJed Brown      deallocate (cols)
44c4762a1bSJed Brown      deallocate (vals)
45c4762a1bSJed Brown    else
46c4762a1bSJed Brown      allocate (vals(3))
47c4762a1bSJed Brown      vals = 1.0
48c4762a1bSJed Brown      allocate (cols(3), source=[r - 1, r, r + 1])
495d83a8b1SBarry Smith      PetscCallA(MatSetValues(A, one, [r], three, cols, vals, INSERT_VALUES, ierr))
50c4762a1bSJed Brown      deallocate (cols)
51c4762a1bSJed Brown      deallocate (vals)
52c4762a1bSJed Brown    end if
53c4762a1bSJed Brown  end do
54d8606c27SBarry Smith  PetscCallA(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY, ierr))
55d8606c27SBarry Smith  PetscCallA(MatAssemblyend(A, MAT_FINAL_ASSEMBLY, ierr))
56d8606c27SBarry Smith  PetscCallA(MatPartitioningCreate(PETSC_COMM_WORLD, part, ierr))
57d8606c27SBarry Smith  PetscCallA(MatPartitioningSetAdjacency(part, A, ierr))
58d8606c27SBarry Smith  PetscCallA(MatPartitioningSetFromOptions(part, ierr))
59d8606c27SBarry Smith  PetscCallA(MatPartitioningApply(part, is, ierr))
60d8606c27SBarry Smith  PetscCallA(ISView(is, PETSC_VIEWER_STDOUT_WORLD, ierr))
61d8606c27SBarry Smith  PetscCallA(ISDestroy(is, ierr))
62d8606c27SBarry Smith  PetscCallA(MatPartitioningDestroy(part, ierr))
63d8606c27SBarry Smith  PetscCallA(MatDestroy(A, ierr))
64d8606c27SBarry Smith  PetscCallA(PetscFinalize(ierr))
65c4762a1bSJed Brown
66c4762a1bSJed Brownend program
67c4762a1bSJed Brown
68c4762a1bSJed Brown!/*TEST
69c4762a1bSJed Brown!
70c4762a1bSJed Brown!   test:
71c4762a1bSJed Brown!      nsize: 3
72c4762a1bSJed Brown!      requires: parmetis
73c4762a1bSJed Brown!      args: -mat_partitioning_type parmetis
74c4762a1bSJed Brown!      output_file: output/ex15_1.out
75c4762a1bSJed Brown!
76c4762a1bSJed Brown!   test:
77c4762a1bSJed Brown!      suffix: 2
78c4762a1bSJed Brown!      nsize: 3
79c4762a1bSJed Brown!      requires: ptscotch
80c4762a1bSJed Brown!      args: -mat_partitioning_type ptscotch
81c4762a1bSJed Brown!      output_file: output/ex15_2.out
82c4762a1bSJed Brown!
83c4762a1bSJed Brown!   test:
84c4762a1bSJed Brown!      suffix: 3
85c4762a1bSJed Brown!      nsize: 4
86c4762a1bSJed Brown!      requires: party
87c4762a1bSJed Brown!      args: -mat_partitioning_type party
88c4762a1bSJed Brown!      output_file: output/ex15_3.out
89c4762a1bSJed Brown!
90c4762a1bSJed Brown!   test:
91c4762a1bSJed Brown!      suffix: 4
92c4762a1bSJed Brown!      nsize: 3
93c4762a1bSJed Brown!      requires: chaco
94c4762a1bSJed Brown!      args: -mat_partitioning_type chaco
95c4762a1bSJed Brown!      output_file: output/ex15_4.out
96c4762a1bSJed Brown!
97c4762a1bSJed Brown!   test:
98c4762a1bSJed Brown!      suffix: 5
99c4762a1bSJed Brown!      nsize: 3
100c4762a1bSJed Brown!      requires: parmetis
101c4762a1bSJed Brown!      args: -mat_partitioning_type hierarch -mat_partitioning_hierarchical_nfineparts 3 -mat_partitioning_nparts 10 -N 100
102c4762a1bSJed Brown!      output_file: output/ex15_5.out
103c4762a1bSJed Brown!
104c4762a1bSJed Brown!TEST*/
105