FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_ordering.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2020 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
5 
6 !----------------------------------------------------------------------
8 ! for direct solver
9 !----------------------------------------------------------------------
11  use hecmw_util
12  implicit none
13 
14  private
15  public :: hecmw_ordering_gen
16 
17  integer(kind=kint), parameter:: ORDERING_DEFAULT = 0
18  integer(kind=kint), parameter:: ORDERING_QMD = 1
19  integer(kind=kint), parameter:: ORDERING_METIS = 2
20  integer(kind=kint), parameter:: ORDERING_RCM = 3
21  integer(kind=kint), parameter:: ORDERING_NMAX = 3
22 
23  integer(kind=kint), parameter:: ORDERING_DEBUG = 0
24 
25 contains
26 
27  !======================================================================!
29  !======================================================================!
30  subroutine hecmw_ordering_gen(Neqns,Nttbr,Xadj,Adj0,Perm,Invp,opt,loglevel)
34  implicit none
35  !------
36  integer(kind=kint), intent(in):: neqns
37  integer(kind=kint), intent(in):: nttbr
38  integer(kind=kint), intent(in):: adj0(:)
39  integer(kind=kint), intent(in):: xadj(:)
40  integer(kind=kint), intent(in):: opt
41  integer(kind=kint), intent(in):: loglevel
42  integer(kind=kint), intent(out):: perm(:)
43  integer(kind=kint), intent(out):: invp(:)
44  !------
45  integer(kind=kint):: ordering
46  ordering = opt
47  if (ordering < 0 .or. ordering > ordering_nmax) then
48  stop "ERROR ordering option for direct solver out of range"
49  endif
50  if (ordering == ordering_default) then
51 #ifdef HECMW_WITH_METIS
52  ordering = ordering_metis
53 #else
54  ordering = ordering_qmd
55 #endif
56  endif
57  select case (ordering)
58  case(ordering_qmd)
59  if (loglevel > 0) write(*,*) 'Ordering method: QMD'
60  call hecmw_ordering_genqmd(neqns,nttbr,xadj,adj0,perm,invp)
61  case(ordering_metis)
62  if (loglevel > 0) write(*,*) 'Ordering method: METIS_NodeND'
63  call hecmw_ordering_metis_nodend(neqns,xadj,adj0,perm,invp)
64  case(ordering_rcm)
65  if (loglevel > 0) write(*,*) 'Ordering method: RCM'
66  call hecmw_ordering_genrcm(neqns,xadj,adj0,perm,invp)
67  end select
68  if (ordering_debug > 0) then
69  call write_nonzero_profile(neqns, xadj, adj0, perm, invp)
70  call write_perm(neqns, perm, invp)
71  endif
72  end subroutine hecmw_ordering_gen
73 
74  subroutine write_nonzero_profile(N, index, item, perm, iperm)
75  implicit none
76  integer(kind=kint), intent(in) :: N
77  integer(kind=kint), intent(in) :: index(:)
78  integer(kind=kint), intent(in) :: item(:)
79  integer(kind=kint), intent(in) :: perm(:), iperm(:)
80  integer(kind=kint), parameter :: F_ORG = 901
81  integer(kind=kint), parameter :: F_NEW = 902
82  integer(kind=kint) :: i, j, irow, jcol
83  open(f_org, file='nzprof_org.txt', status='replace')
84  do irow = 1, n
85  i = irow
86  do j = index(i), index(i+1)-1
87  jcol = item(j)
88  write(f_org,*) irow, jcol
89  end do
90  end do
91  close(f_org)
92  open(f_new, file='nzprof_new.txt', status='replace')
93  do irow = 1, n
94  i = perm(irow)
95  do j = index(i), index(i+1)-1
96  jcol = item(j)
97  write(f_new,*) irow, iperm(jcol)
98  end do
99  end do
100  close(f_new)
101  end subroutine write_nonzero_profile
102 
103  subroutine write_perm(N, perm, iperm)
104  implicit none
105  integer(kind=kint), intent(in) :: N
106  integer(kind=kint), intent(in) :: perm(:), iperm(:)
107  integer(kind=kint), parameter :: F_PERM = 903
108  integer(kind=kint) :: i
109  open(f_perm, file='perm_iperm.txt', status='replace')
110  do i = 1, n
111  write(f_perm,*) perm(i), iperm(i)
112  end do
113  close(f_perm)
114  end subroutine write_perm
115 
116 end module hecmw_ordering
hecmw_ordering_qmd
HECMW_ORDERING_QMD is a program for the minimum degree.
Definition: hecmw_ordering_qmd.f90:10
hecmw_ordering
HECMW_ORDERING is a program for fill-reducing ordering.
Definition: hecmw_ordering.F90:10
hecmw_ordering_metis::hecmw_ordering_metis_nodend
subroutine, public hecmw_ordering_metis_nodend(Neqns, Xadj, Adj0, Perm, Invp)
hecmw_ordering_metis_NodeND
Definition: hecmw_ordering_metis.F90:23
hecmw_ordering_rcm
HECMW_ORDERING_RCM is a program for fill-reducing ordering.
Definition: hecmw_ordering_rcm.f90:10
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
hecmw_ordering_rcm::hecmw_ordering_genrcm
subroutine, public hecmw_ordering_genrcm(Neqns, Xadj, Adj0, Perm, Invp)
Definition: hecmw_ordering_rcm.f90:20
hecmw_ordering_metis
HECMW_ORDERING_METIS is a program for the Nested Dissection.
Definition: hecmw_ordering_metis.F90:10
hecmw_ordering::write_nonzero_profile
subroutine write_nonzero_profile(N, index, item, perm, iperm)
Definition: hecmw_ordering.F90:75
hecmw_ordering_qmd::hecmw_ordering_genqmd
subroutine, public hecmw_ordering_genqmd(Neqns, Nttbr, Xadj, Adj0, Perm, Invp)
hecmw_ordering_GENQMD
Definition: hecmw_ordering_qmd.f90:23
hecmw_ordering::hecmw_ordering_gen
subroutine, public hecmw_ordering_gen(Neqns, Nttbr, Xadj, Adj0, Perm, Invp, opt, loglevel)
hecmw_ordering_gen
Definition: hecmw_ordering.F90:31