CONGEN 5.0
Configuration generation for SCATCI
Loading...
Searching...
No Matches
congen_bstree.f90
Go to the documentation of this file.
1! Copyright 2019
2!
3! For a comprehensive list of the developers that contributed to these codes
4! see the UK-AMOR website.
5!
6! This file is part of UKRmol-in (UKRmol+ suite).
7!
8! UKRmol-in is free software: you can redistribute it and/or modify
9! it under the terms of the GNU General Public License as published by
10! the Free Software Foundation, either version 3 of the License, or
11! (at your option) any later version.
12!
13! UKRmol-in is distributed in the hope that it will be useful,
14! but WITHOUT ANY WARRANTY; without even the implied warranty of
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16! GNU General Public License for more details.
17!
18! You should have received a copy of the GNU General Public License
19! along with UKRmol-in (in source/COPYING). Alternatively, you can also visit
20! <https://www.gnu.org/licenses/>.
21
22!> \brief Determinant binary search tree
23!> \authors J Benda
24!> \date 2018 - 2019
25!>
26!> This module contains a special data structure, \ref det_tree, used (e.g.) by CONGEN to carry out quick look-ups in the
27!> storage of determinants. It is based on the type "bstree" from libAMOR.
28!>
30
31 use containers, only: bstree
32
33 implicit none
34
35 !> \brief Determinant binary search tree
36 !> \authors J Benda
37 !> \date 2019
38 !>
39 !> This is an extension of the plain integer binary search tree to one that operates on the determinant
40 !> storage. The relation used in the ordering is the lexicographical order of the determinants as provided
41 !> by the overriden subroutine "compare".
42 !>
43 type, extends(bstree) :: det_tree
44 integer, pointer :: det(:) => null()
45 integer, pointer :: ndo(:) => null()
46 integer :: nelt = 0
47 contains
48 procedure, public :: init => init_det_tree
49 procedure, public :: compare => compare_determinats
50 procedure, public :: locate_det => locate_in_det_tree
51 procedure, public :: dtwrite => write_determinant
53 end type det_tree
54
55contains
56
57 !> \brief Initialize determinant tree
58 !> \authors J Benda
59 !> \date 2019
60 !>
61 !> Stores pointer to the determinant storage, so that it can be used when comparing determinant indices
62 !> in the binary search tree subroutines.
63 !>
64 !> \param this Tree object to initialize
65 !> \param ndo Determinant storage
66 !> \param nelt Size of a determinant
67 !>
68 subroutine init_det_tree (this, ndo, nelt)
69
70 class(det_tree), intent(inout) :: this
71 integer, pointer, intent(in) :: ndo(:)
72 integer, intent(in) :: nelt
73
74 this % ndo => ndo
75 this % nelt = nelt
76
77 allocate (this % det(nelt))
78
79 end subroutine init_det_tree
80
81
82 !> \brief Finalize determinant tree
83 !> \authors J Benda
84 !> \date 2019
85 !>
86 !> Releases all allocated memory.
87 !>
88 !> \param this Tree object to finalize
89 !>
90 subroutine final_det_tree (this)
91
92 type(det_tree), intent(inout) :: this
93
94 if (associated(this % det)) then
95 deallocate (this % det)
96 end if
97
98 end subroutine final_det_tree
99
100
101 !> \brief Find a specific determinant in the storage
102 !> \authors J Benda
103 !> \date 2019
104 !>
105 !> Return index in the determinant storage, where the given determinant is located. Return -1 when there is no such
106 !> determinant.
107 !>
108 !> \param this Tree object to search
109 !> \param det Determinant to find (length nelt)
110 !>
111 integer function locate_in_det_tree (this, det) result (res)
112
113 class(det_tree), intent(inout) :: this
114 integer, dimension(this % nelt), intent(in) :: det
115
116 this % det = det
117
118 res = this % locate(-1) ! negative id signalizes unary comparison, see "compare_determinats"
119
120 end function locate_in_det_tree
121
122
123 !> \brief Write determinat
124 !> \authors J Benda
125 !> \date 2019
126 !>
127 !> Used in debuggind output of the whole tree.
128 !>
129 !> \param this Binary search tree
130 !> \param lu Unit for writing
131 !> \param id Tree node id
132 !>
133 subroutine write_determinant (this, lu, id)
134
135 class(det_tree), intent(in) :: this
136 integer, intent(in) :: lu, id
137
138 integer :: i
139
140 write (lu, '("[")', advance = 'no')
141
142 do i = 1, this % nelt
143 write (lu, '(1x,I0)', advance = 'no') this % ndo((id - 1) * this % nelt + i)
144 end do
145
146 write (lu, '(1x,"]")')
147
148 end subroutine write_determinant
149
150
151 !> \brief Lexicographically compare two determinats
152 !> \authors J Benda
153 !> \date 2018 - 2019
154 !>
155 !> Binary tree needs a notion of order of its elements to be able to work with the data
156 !> so efficiently. In case of arrays of numbers, the typical order is the "lexicographical"
157 !> order. The array A is less than B if for the first position I for which A(I) /= B(I)
158 !> holds that A(I) < B(I).
159 !>
160 !> If any of the given indices is negative, the reference determinant stored within the
161 !> predicate type is used instead.
162 !>
163 !> \param this Determinant tree object to use
164 !> \param i Index of the first determinant
165 !> \param j Index of the second determinant
166 !> \param data Additional payload required by base class interface. Not used in congen.
167 !>
168 !> \return -1 if the first determinant is lexicographically less,
169 !> 0 if they are equal,
170 !> +1 if the first is greater.
171 !>
172 integer function compare_determinats (this, i, j, data) result (verdict)
173
174 use iso_c_binding, only: c_ptr
175
176 class(det_tree), intent(in) :: this
177 integer, intent(in) :: i, j
178 type(c_ptr), optional, intent(in) :: data
179
180 integer, pointer :: d1(:), d2(:)
181 integer :: k
182
183 ! get first determinant view
184 if (i >= 0) then
185 d1 => this % ndo((i - 1) * this % nelt + 1 : i * this % nelt)
186 else
187 d1 => this % det(1:this % nelt)
188 end if
189
190 ! get second determinant view
191 if (j >= 0) then
192 d2 => this % ndo((j - 1) * this % nelt + 1 : j * this % nelt)
193 else
194 d2 => this % det(1:this % nelt)
195 end if
196
197 ! compare determinants element by element
198 do k = 1, this % nelt
199
200 ! is first less?
201 if (d1(k) < d2(k)) then
202 verdict = -1
203 return
204 end if
205
206 ! is first greater?
207 if (d1(k) > d2(k)) then
208 verdict = 1
209 return
210 end if
211
212 end do
213
214 ! they are equal
215 verdict = 0
216
217 end function compare_determinats
218
219end module congen_bstree
Determinant binary search tree.
subroutine init_det_tree(this, ndo, nelt)
Initialize determinant tree.
subroutine write_determinant(this, lu, id)
Write determinat.
integer function compare_determinats(this, i, j, data)
Lexicographically compare two determinats.
integer function locate_in_det_tree(this, det)
Find a specific determinant in the storage.
subroutine final_det_tree(this)
Finalize determinant tree.
Determinant binary search tree.