Line data Source code
1 : ! ***********************************************************************
2 : !
3 : ! Copyright (C) 2010 Ed Brown & The MESA Team
4 : !
5 : ! This program is free software: you can redistribute it and/or modify
6 : ! it under the terms of the GNU Lesser General Public License
7 : ! as published by the Free Software Foundation,
8 : ! either version 3 of the License, or (at your option) any later version.
9 : !
10 : ! This program is distributed in the hope that it will be useful,
11 : ! but WITHOUT ANY WARRANTY; without even the implied warranty of
12 : ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
13 : ! See the GNU Lesser General Public License for more details.
14 : !
15 : ! You should have received a copy of the GNU Lesser General Public License
16 : ! along with this program. If not, see <https://www.gnu.org/licenses/>.
17 : !
18 : ! ***********************************************************************
19 :
20 : module nuclide_set_mod
21 :
22 : use chem_def
23 : use const_def, only: dp
24 :
25 : implicit none
26 :
27 : contains
28 :
29 1381026 : function rank_in_set(iso, set)
30 : character(len=iso_name_length), intent(in) :: iso
31 : type(nuclide_set), dimension(:), intent(in) :: set
32 : integer :: rank_in_set
33 : integer :: low, high, mid, i
34 : integer, parameter :: max_cycles = 20
35 1381026 : low = 1
36 1381026 : high = size(set)
37 1381026 : if (adjustl(iso) < adjustl(set(low)% nuclide) .or. adjustl(iso) > adjustl(set(high)% nuclide)) then
38 20 : rank_in_set = nuclide_not_found
39 20 : return
40 : end if
41 19050759 : do i = 1, max_cycles
42 19050759 : if (high-low <=1) then
43 1381006 : if (adjustl(iso) == adjustl(set(high)% nuclide)) then
44 1375821 : rank_in_set = set(high)% rank
45 5185 : else if (adjustl(iso) == adjustl(set(low)% nuclide)) then
46 20 : rank_in_set = set(low)% rank
47 : else
48 : rank_in_set = nuclide_not_found
49 : end if
50 1381006 : return
51 : end if
52 17669753 : mid = (high+low)/2
53 17669753 : if (adjustl(iso) <= adjustl(set(mid)% nuclide)) then
54 : high = mid
55 8463299 : else if (adjustl(iso) > adjustl(set(mid)% nuclide)) then
56 8463299 : low = mid
57 : end if
58 : end do
59 0 : rank_in_set = nuclide_not_found
60 0 : end function rank_in_set
61 :
62 :
63 5 : subroutine sort_nuclide_set(set)
64 : type(nuclide_set), dimension(:), intent(inout) :: set
65 : integer :: n, i, ir, j, l
66 : type(nuclide_set) :: ts
67 :
68 5 : n = size(set)
69 5 : if (size(set) < 2) return
70 5 : l = n/2+1
71 5 : ir = n
72 58910 : do
73 58915 : if (l > 1) then
74 19640 : l = l-1
75 19640 : ts = set(l)
76 : else
77 39275 : ts = set(ir)
78 39275 : set(ir) = set(1)
79 39275 : ir = ir-1
80 39275 : if (ir == 1) then
81 5 : set(1) = ts
82 5 : return
83 : end if
84 : end if
85 58910 : i = l
86 58910 : j = l+l
87 : do
88 500880 : if (j > ir) exit
89 441970 : if (j < ir) then
90 441895 : if (compare_lt(set(j), set(j+1))) j = j+1
91 : end if
92 500880 : if (compare_lt(ts, set(j))) then
93 422520 : set(i) = set(j)
94 422520 : i = j
95 422520 : j = j+j
96 : else
97 19450 : j = ir + 1
98 : end if
99 : end do
100 58910 : set(i) = ts
101 : end do
102 :
103 : contains
104 :
105 883865 : logical function compare_lt(a, b)
106 : type(nuclide_set), intent(in) :: a, b
107 883865 : compare_lt = (adjustl(a% nuclide) < adjustl(b% nuclide))
108 883865 : end function compare_lt
109 :
110 : end subroutine sort_nuclide_set
111 :
112 : end module nuclide_set_mod
|