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