1 | !======================================================================!
|
2 | ! !
|
3 | ! Software Name : HEC-MW Library for PC-cluster !
|
4 | ! Version : 2.7 !
|
5 | ! !
|
6 | ! Last Update : 2014/01/25 !
|
7 | ! Category : Linear Solver !
|
8 | ! !
|
9 | ! Written by Kazuya Goto (PExProCS LLC) !
|
10 | ! !
|
11 | ! Contact address : IIS,The University of Tokyo RSS21 project !
|
12 | ! !
|
13 | ! "Structural Analysis System for General-purpose Coupling !
|
14 | ! Simulations Using High End Computing Middleware (HEC-MW)" !
|
15 | ! !
|
16 | !======================================================================!
|
17 |
|
18 | module m_hecmw_matrix_ordering_MC
|
19 | use hecmw_util
|
20 | implicit none
|
21 |
|
22 | private
|
23 | public :: hecmw_matrix_ordering_MC
|
24 |
|
25 | contains
|
26 |
|
27 | subroutine hecmw_matrix_ordering_MC(N, indexL, itemL, indexU, itemU, &
|
28 | perm_cur, ncolor_in, ncolor_out, COLORindex, perm, iperm)
|
29 | implicit none
|
30 | integer(kind=kint), intent(in) :: N
|
31 | integer(kind=kint), intent(in) :: indexL(0:), indexU(0:)
|
32 | integer(kind=kint), intent(in) :: itemL(:), itemU(:)
|
33 | integer(kind=kint), intent(in) :: perm_cur(:)
|
34 | integer(kind=kint), intent(in) :: ncolor_in
|
35 | integer(kind=kint), intent(out) :: ncolor_out
|
36 | integer(kind=kint), intent(out) :: COLORindex(0:)
|
37 | integer(kind=kint), intent(out) :: perm(:), iperm(:)
|
38 | integer(kind=kint), allocatable :: iwk(:)
|
39 | integer(kind=kint) :: nn_color, cntall, cnt, color
|
40 | integer(kind=kint) :: i, inode, j, jnode
|
41 | allocate(iwk(N))
|
42 | iwk = 0
|
43 | nn_color = N / ncolor_in
|
44 | cntall = 0
|
45 | COLORindex(0) = 0
|
46 | do color=1,N
|
47 | cnt = 0
|
48 | do i=1,N
|
49 | inode = perm_cur(i)
|
50 | if (iwk(inode) > 0 .or. iwk(inode) == -1) cycle
|
51 | ! if (iwk(inode) == 0)
|
52 | iwk(inode) = color
|
53 | cntall = cntall + 1
|
54 | perm(cntall) = inode
|
55 | cnt = cnt + 1
|
56 | if (cnt == nn_color) exit
|
57 | if (cntall == N) exit
|
58 | ! mark all connected and uncolored nodes
|
59 | do j = indexL(inode-1)+1, indexL(inode)
|
60 | jnode = itemL(j)
|
61 | if (iwk(jnode) == 0) iwk(jnode) = -1
|
62 | end do
|
63 | do j = indexU(inode-1)+1, indexU(inode)
|
64 | jnode = itemU(j)
|
65 | if (jnode > N) cycle
|
66 | if (iwk(jnode) == 0) iwk(jnode) = -1
|
67 | end do
|
68 | end do
|
69 | COLORindex(color) = cntall
|
70 | if (cntall == N) then
|
71 | ncolor_out = color
|
72 | exit
|
73 | end if
|
74 | ! unmark all marked nodes
|
75 | do i=1,N
|
76 | if (iwk(i) == -1) iwk(i) = 0
|
77 | end do
|
78 | end do
|
79 | deallocate(iwk)
|
80 | ! make iperm
|
81 | do i=1,N
|
82 | iperm(perm(i)) = i
|
83 | end do
|
84 | end subroutine hecmw_matrix_ordering_MC
|
85 |
|
86 | end module m_hecmw_matrix_ordering_MC
|