Merge lp:~albertog/siesta/4.0-pool into lp:siesta/4.0

Proposed by Alberto Garcia
Status: Merged
Merged at revision: 593
Proposed branch: lp:~albertog/siesta/4.0-pool
Merge into: lp:siesta/4.0
Diff against target: 109 lines (+47/-10) (has conflicts)
2 files modified
Src/matel_registry.F90 (+42/-10)
version.info (+5/-0)
Text conflict in version.info
To merge this branch: bzr merge lp:~albertog/siesta/4.0-pool
Reviewer Review Type Date Requested Status
Nick Papior Approve
Review via email: mp+364488@code.launchpad.net

Description of the change

Make radial function pool dynamically extensible.

To post a comment you must log in.
Revision history for this message
Nick Papior (nickpapior) wrote :

Nice! Great job!

LGTM!

review: Approve

Preview Diff

[H/L] Next/Prev Comment, [J/K] Next/Prev File, [N/P] Next/Prev Hunk
=== modified file 'Src/matel_registry.F90'
--- Src/matel_registry.F90 2016-01-25 16:00:16 +0000
+++ Src/matel_registry.F90 2019-03-15 09:23:30 +0000
@@ -8,7 +8,7 @@
8module m_matel_registry8module m_matel_registry
9 !9 !
10 ! This module provides a general interface to the10 ! This module provides a general interface to the
11 ! functions needed by matel11 ! functions needed by (new_)matel
12 ! 12 !
13 ! Sketch of usage (see examples in the code):13 ! Sketch of usage (see examples in the code):
14 !14 !
@@ -20,7 +20,7 @@
20 !20 !
21 ! Once registered, a function can be evaluated, its cutoff and 21 ! Once registered, a function can be evaluated, its cutoff and
22 ! angular momentum obtained, etc. So far only the functions22 ! angular momentum obtained, etc. So far only the functions
23 ! needed by matel are implemented: rcut, lcut, evaluate (former phiatm).23 ! needed by (new_)matel are implemented: rcut, lcut, evaluate (former phiatm).
24 !24 !
25 ! The registry contains some meta-data for each function (its25 ! The registry contains some meta-data for each function (its
26 ! kind as a string, and the original indexes handled by the26 ! kind as a string, and the original indexes handled by the
@@ -71,13 +71,13 @@
71 integer :: lcut = -huge(1)71 integer :: lcut = -huge(1)
72 end type registered_function_t72 end type registered_function_t
73 !73 !
74 ! There should be a mechanism to make the size74 ! This is a dynamic array, allocated and extended as needed.
75 ! of the pool dynamic.75 type(registered_function_t), dimension(:), allocatable, save :: matel_pool
76 !
77 integer, parameter :: nmax_funcs = 300
78 type(registered_function_t), dimension(nmax_funcs), save :: matel_pool
7976
80 integer, private :: nfuncs = 077 integer, parameter :: initial_size = 100 ! initial size of pool
78 integer, parameter :: delta_size = 50 ! chunk size for extensions
79 integer, private :: nmax_funcs = 0 ! current maximum size of the pool
80 integer, private :: nfuncs = 0 ! number of functions in pool
8181
82 public :: register_in_rf_pool, register_in_tf_pool82 public :: register_in_rf_pool, register_in_tf_pool
83 public :: evaluate, rcut, lcut83 public :: evaluate, rcut, lcut
@@ -93,6 +93,36 @@
93 ok = (gindex > 0 .AND. gindex <= nfuncs)93 ok = (gindex > 0 .AND. gindex <= nfuncs)
94 end function valid94 end function valid
9595
96 !> Extends the size of the pool if needed
97 subroutine check_size_of_pool(nfuncs)
98 integer, intent(in) :: nfuncs
99
100 type(registered_function_t), dimension(:), allocatable :: tmp_pool
101
102 if (nfuncs > nmax_funcs) then
103 if (.not. allocated(matel_pool)) then
104 allocate(matel_pool(initial_size))
105 nmax_funcs = initial_size
106 else
107 ! copy data to tmp array
108 allocate(tmp_pool(nmax_funcs))
109
110 ! Each element (derived type) involves just two scalars
111 ! and two pointers, so the copy overhead is small.
112 ! Note that we do not use the F2003 "move_alloc" idiom.
113
114 tmp_pool(1:nmax_funcs) = matel_pool(1:nmax_funcs)
115 deallocate(matel_pool)
116 allocate(matel_pool(nmax_funcs + delta_size))
117 matel_pool(1:nmax_funcs) = tmp_pool(1:nmax_funcs)
118 nmax_funcs = nmax_funcs + delta_size
119
120 deallocate(tmp_pool)
121
122 endif
123 endif
124 end subroutine check_size_of_pool
125
96 !126 !
97 ! This is the main entry to the registry for simple radial functions127 ! This is the main entry to the registry for simple radial functions
98 !128 !
@@ -109,7 +139,8 @@
109 n_indexes = size(indexes)139 n_indexes = size(indexes)
110140
111 nfuncs = nfuncs + 1141 nfuncs = nfuncs + 1
112 if (nfuncs > nmax_funcs) call die("Overflow in registry")142 call check_size_of_pool(nfuncs)
143
113 gindex = nfuncs144 gindex = nfuncs
114145
115 allocate(matel_pool(gindex)%rf)146 allocate(matel_pool(gindex)%rf)
@@ -139,7 +170,8 @@
139 integer, intent(out) :: gindex ! global index170 integer, intent(out) :: gindex ! global index
140171
141 nfuncs = nfuncs + 1172 nfuncs = nfuncs + 1
142 if (nfuncs > nmax_funcs) call die("Overflow in registry")173 call check_size_of_pool(nfuncs)
174
143 gindex = nfuncs175 gindex = nfuncs
144176
145 allocate(matel_pool(gindex)%tf)177 allocate(matel_pool(gindex)%tf)
146178
=== modified file 'version.info'
--- version.info 2019-02-20 19:40:46 +0000
+++ version.info 2019-03-15 09:23:30 +0000
@@ -1,1 +1,6 @@
1<<<<<<< TREE
1siesta-4.0--5922siesta-4.0--592
3=======
4siesta-4.0--591--pool-593
5
6>>>>>>> MERGE-SOURCE

Subscribers

People subscribed via source and target branches