MPI-SCATCI 2.0
An MPI version of SCATCI
Loading...
Searching...
No Matches
MemoryManager_module.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 Memory manager module
23!> \authors A Al-Refaie
24!> \date 2017
25!>
26!> This module has an extremely basic memory management system in place, this will hopefully be upgraded in due time
27!>
28!> \note 30/01/2017 - Ahmed Al-Refaie: Initial revision.
29!> \note 16/01/2019 - Jakub Benda: Unifom coding style and expanded documentation.
30!>
31module memorymanager_module
32
33 use precisn, only: longint, wp
34 use const_gbl, only: stdout
35
36 implicit none
37
38 public master_memory
39
40 private
41
42 !> \brief This is a simple class to handle memory management tracking
43 !> \authors A Al-Refaie
44 !> \date 2017
45 !>
46 type :: memorymanager
47 private
48 integer(longint) :: total_local_memory = 0 !< The total memory we have been assigned
49 integer(longint) :: available_local_memory = 0 !< Memory we have left over
50 integer(longint) :: cached_tracked_memory = 0
51 logical :: initialized = .false.
52 contains
53 private
54 procedure, public :: construct => memorymanager_ctor
55 procedure, public :: init_memory !< Initialize the memory with a new value
56 procedure, public :: track_memory !< Tracks
57 procedure, public :: free_memory
58 procedure, public :: get_available_memory
59 procedure, public :: get_scaled_available_memory
60 procedure, public :: print_memory_report
61 end type memorymanager
62
63 !> \brief This is the global memory manager.
64 !> \authors A Al-Refaie
65 !> \date 2017
66 !>
67 !> No other memory manager can exist
68 !>
69 type(memorymanager) :: master_memory
70
71contains
72
73 !> \brief Constructor, used to initialize total available memory
74 !> \authors A Al-Refaie
75 !> \date 2017
76 !>
77 !> \param[inout] this Manager object to update.
78 !> \param[in] total_memory Total available memory for the process in bytes
79 !>
80 subroutine memorymanager_ctor (this, total_memory)
81 class(memorymanager) :: this
82 integer(longint), intent(in) :: total_memory
83
84 this % total_local_memory = total_memory
85 this % available_local_memory = this % total_local_memory
86 this % initialized = .true.
87 this % available_local_memory = this % available_local_memory - this % cached_tracked_memory
88
89 end subroutine memorymanager_ctor
90
91
92 !> \brief Same as constructor, used to reinitialize te memory again (not really used)
93 !> \authors A Al-Refaie
94 !> \date 2017
95 !>
96 !> \param[inout] this Manager object to update.
97 !> \param[in] total_memory Total available memory for the process in bytes
98 !>
99 subroutine init_memory (this, total_memory)
100 class(memorymanager) :: this
101 integer(longint), intent(in) :: total_memory
102
103 this % total_local_memory = total_memory
104 this % available_local_memory = this % total_local_memory
105 this % available_local_memory = this % available_local_memory - this % cached_tracked_memory
106 this % initialized = .true.
107
108 end subroutine init_memory
109
110
111 !> \brief Tracks memory allocation, usually called after an allocation
112 !> \authors A Al-Refaie
113 !> \date 2017
114 !>
115 !> \param[inout] this Manager object to update.
116 !> \param[in] elem_size Size in byte of each element
117 !> \param[in] nelem Number of elements of size elem_size
118 !> \param[in] stat error Value after an allocate call
119 !> \param[in] array_name Assigned name of array, used to debug where a memory allocation error occured.
120 !>
121 subroutine track_memory (this, elem_size, nelem, stat, array_name)
122 class(memorymanager) :: this
123 character(len=*), intent(in) :: array_name
124 integer, intent(in) :: elem_size,nelem,stat
125 integer(longint) :: total_mem
126
127 !Check if allocation was successfull
128 if (stat /= 0) then
129 !If not print name and error stats
130 write (stdout, "('Array : ',a)") array_name
131 write (stdout, "('stat returned error ',i16,' when trying to allocate nelems ',i16,' of size ', i16)") &
132 stat, nelem, elem_size
133 !Halt the program
134 stop "[Memory manager] Memory allocation error"
135 end if
136
137 ! calculate in bytes total memory used
138 total_mem = elem_size * nelem
139
140 !Remove from available memory
141 !$OMP ATOMIC
142 this % cached_tracked_memory = this % cached_tracked_memory + total_mem
143 if (.not. this % initialized) return
144 !$OMP ATOMIC
145 this % available_local_memory = this % available_local_memory - total_mem
146
147 !If we've hit negative then throw error
148 if (this % available_local_memory < 0) then
149 call this % print_memory_report()
150 write (stdout, "('Array : ',a)") array_name
151 write (stdout, "('Overrun allowed memory space when trying to allocate nelems ',i16,' of size ', 2i16,' vs',i16)") &
152 nelem, elem_size, total_mem, this % get_available_memory()
153 stop "[Memory manager] Overrun space"
154 end if
155 end subroutine track_memory
156
157
158 !> \brief Tracks memory deallocation, usually called before a deallocation
159 !> \authors A Al-Refaie
160 !> \date 2017
161 !>
162 !> \param[inout] this Manager object to update.
163 !> \param[in] elem_size Size in byte of each element
164 !> \param[in] nelem Number of elements of size elem_size
165 !>
166 subroutine free_memory (this, elem_size, nelem)
167 class(memorymanager) :: this
168 integer, intent(in) :: elem_size, nelem
169 integer(longint) :: total_mem
170
171 ! Calculate total deallocated
172 total_mem = elem_size * nelem
173 !$OMP ATOMIC
174 this % cached_tracked_memory = this % cached_tracked_memory - total_mem
175 !Add memory back to pool
176 if (.not. this % initialized) return
177 !$OMP ATOMIC
178 this % available_local_memory = this % available_local_memory + total_mem
179
180 !If we have more memory than available then flag as there is a mismatching of track and free
181 if (this % available_local_memory > this % total_local_memory) then
182 stop "Memory mismatch"
183 end if
184
185 end subroutine free_memory
186
187
188 !> \brief Get currently available memory
189 !> \authors A Al-Refaie
190 !> \date 2017
191 !>
192 !> \result Total memory in bytes
193 !>
194 integer(longint) function get_available_memory (this)
195 class(memorymanager) :: this
196
197 get_available_memory = this % available_local_memory
198
199 end function get_available_memory
200
201
202 !> \brief Get currently available memory scaled.
203 !> \authors A Al-Refaie
204 !> \date 2017
205 !>
206 !> This is particularly useful when trying to determine how many elements we can fit into memory by giving a margin of safety
207 !>
208 !> \param[inout] this Manager object to update.
209 !> \param[in] scale_ Scale valued of available memory (0.0 - 1.0).
210 !>
211 !> \result Total memory in bytes
212 !>
213 integer(longint) function get_scaled_available_memory (this, scale_)
214 class(memorymanager) :: this
215 real(wp), intent(in) :: scale_
216 real(wp) :: s
217
218 !Clamp scale between 0.0 ad 1.0
219 s = min(abs(scale_), 1.0_wp)
220
221 !Return scaled
222 get_scaled_available_memory = int(real(this % available_local_memory) * s, longint)
223
224 end function get_scaled_available_memory
225
226
227 !> \brief Prints a very simple memory report
228 !> \authors A Al-Refaie
229 !> \date 2017
230 !>
231 subroutine print_memory_report (this)
232 class(memorymanager) :: this
233 real(wp) :: mem_total, mem_avail, percentage
234
235 !Compute total in gibibytes
236 mem_total = real(this % total_local_memory) / real(1024**3)
237
238 !Compute available in gibibytes
239 mem_avail = real(this % available_local_memory) / real(1024**3)
240
241 !Compute percentage available
242 percentage = mem_avail * 100.0 / mem_total
243 write (stdout, "('Currently we have ',f18.8,'GB / ',f18.8,'GiB available')") mem_avail, mem_total
244 write (stdout, "('Memory availability at ',f6.1,'%')") percentage
245
246 end subroutine print_memory_report
247
248end module memorymanager_module