10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-06-02 11:25:26 +02:00

Forgot memory.irp.f

This commit is contained in:
Anthony Scemama 2018-12-14 09:37:28 +01:00
parent 69c6361cd7
commit fe33840e18

98
src/Utils/memory.irp.f Normal file
View File

@ -0,0 +1,98 @@
subroutine resident_memory(value)
implicit none
BEGIN_DOC
! Returns the current used memory in gigabytes used by the current process.
END_DOC
integer :: iunit
integer, external :: getUnitAndOpen
character*(32) :: key
double precision, intent(out) :: value
value = 0.d0
iunit = getUnitAndOpen('/proc/self/status','r')
do
read(iunit,*,err=10,end=20) key, value
if (trim(key) == 'VmRSS:') then
exit
endif
10 continue
end do
20 continue
close(iunit)
value = value / (1024.d0*1024.d0)
end function
subroutine total_memory(value)
implicit none
BEGIN_DOC
! Returns the current used memory in gigabytes used by the current process.
END_DOC
integer :: iunit
integer, external :: getUnitAndOpen
character*(32) :: key
double precision, intent(out) :: value
iunit = getUnitAndOpen('/proc/self/status','r')
do
read(iunit,*,err=10,end=20) key, value
if (trim(key) == 'VmSize:') then
exit
endif
10 continue
end do
20 continue
close(iunit)
value = value / (1024.d0*1024.d0)
end function
double precision function memory_of_double(n)
implicit none
BEGIN_DOC
! Computes the memory required for n double precision elements in gigabytes.
END_DOC
integer, intent(in) :: n
double precision, parameter :: f = 8.d0 / (1024.d0*1024.d0*1024.d0)
memory_of_double = dble(n) * f
end function
double precision function memory_of_int(n)
implicit none
BEGIN_DOC
! Computes the memory required for n double precision elements in gigabytes.
END_DOC
integer, intent(in) :: n
double precision, parameter :: f = 4.d0 / (1024.d0*1024.d0*1024.d0)
memory_of_int = dble(n) * f
end function
subroutine check_mem(rss_in,routine)
implicit none
BEGIN_DOC
! Checks if n gigabytes can be allocated. If not, the program exits
END_DOC
double precision, intent(in) :: rss_in
character*(*) :: routine
double precision :: rss
!$OMP CRITICAL
call resident_memory(rss)
rss += rss_in
if (int(rss)+1 > qp_max_mem) then
print *, 'Not enough memory: aborting in ', routine
stop -1
endif
!$OMP END CRITICAL
end
subroutine print_memory_usage()
implicit none
BEGIN_DOC
! Prints the memory usage in the output
END_DOC
double precision :: rss, mem
call resident_memory(rss)
call total_memory(mem)
write(*,'(A,F14.6,A,F14.6,A)') &
'.. >>>>> [ RES MEM : ', rss , &
' GB ] [ VIRT MEM : ', mem, ' GB ] <<<<< ..'
end