TABLE OF CONTENTS


ABINIT/m_profiling_abi [ Modules ]

[ Top ] [ Modules ]

NAME

 m_profiling_abi

FUNCTION

  This module is used for tracing memory allocations/deallocations
  when we compile the code with --enable-memory-profiling="yes" that,
  in turn, defines the CPP macro HAVE_MEM_PROFILE in abi_common.h
  The main entry point is abimem_init. abimem_record is interfaced via CPP macros
  defined in abi_common

COPYRIGHT

 Copyright (C) 2010-2024 ABINIT group (MG)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .

NOTES

 This module is not thread-safe. However, this should not represent
 a significant limitation since memory-tracing is only enabled in debug mode.

m_abimem/abimem_abort [ Functions ]

[ Top ] [ Functions ]

NAME

  abimem_abort

FUNCTION

  Stop the code if an error occurs.

 INPUT
  msg=Error message
  file=File name
  line=Line number

SOURCE

515 subroutine abimem_abort(msg, file, line)
516 
517 !Arguments ------------------------------------
518  integer,intent(in) :: line
519  character(len=*),intent(in) :: msg,file
520 
521 !Local variables-------------------------------
522  integer :: ierr
523 
524 !Local variables-------------------------------
525  integer :: unt_found
526  logical :: isopen
527 ! *************************************************************************
528 
529  write(std_out,*)trim(msg),", file: ", trim(file), ", line: ", line
530 
531  ! Close logfile if it's connected to flush io buffers and avoid file corruption
532  inquire(file=minfo%logfile, number=unt_found, opened=isopen)
533  if (isopen .and. (unt_found == minfo%logunt)) close(unit=minfo%logunt)
534 
535  ierr = 0
536 #ifdef HAVE_MPI
537  call MPI_ABORT(MPI_COMM_WORLD, MPI_ERR_UNKNOWN, ierr)
538 #endif
539  stop
540 
541 end subroutine abimem_abort

m_profiling_abi/abimem_basename [ Functions ]

[ Top ] [ m_profiling_abi ] [ Functions ]

NAME

 abimem_basename

FUNCTION

  Returns the final component of a pathname.

INPUTS

  string=The input string

NOTES

  * If the input string in not a valid path to a file (i.e not in the form foo/name)
    a blank strink is returned
  * We do a backward search becase we want to optimize the algorithm for Fortran strings.

SOURCE

563 pure function abimem_basename(string) result(basename)
564 
565  character(len=*),intent(in) :: string
566  character(len=LEN_TRIM(string)) :: basename
567 
568 !Local variables-------------------------------
569  integer :: ic,nch_trim,nch
570  character(len=1),parameter :: DIR_SEPARATOR = '/'
571  character(len=1),parameter :: BLANK=' '
572 !************************************************************************
573 
574  nch     =LEN     (string)
575  nch_trim=LEN_TRIM(string)
576 
577  ic = INDEX (TRIM(string), DIR_SEPARATOR, back=.TRUE.)
578  !write(*,*)'DEBUG ',TRIM(string),ic
579 
580  if (ic >= 1 .and. ic <= nch_trim-1) then ! there is stuff after the separator.
581    basename = string(ic+1:nch_trim)
582    return
583  else if (ic==0 .or. ic == nch_trim+1) then ! no separator in string or zero length string,
584    basename = TRIM(string)                   ! return trimmed string.
585    return
586  else              ! (ic == nch_trim) separator is the last char.
587    basename= BLANK  ! This is not a valid path to a file, return blank.
588    return
589  end if
590 
591 end function abimem_basename

m_profiling_abi/abimem_get_info [ Functions ]

[ Top ] [ m_profiling_abi ] [ Functions ]

NAME

 abimem_get_info

FUNCTION

  Function that returns the number of allocations and deallocations that have
  been performed in Fortran and the memory currently used

OUTPUT

  nalloc: number of allocations that have been performed in Fortran
  nfree:  number of deallocations that have been performed in Fortran
  allocmemory:  total memory used (Fortran)
  nalloc_c, nfree_c: Similar to nalloc and nfree but for the C code.

m_profiling_abi/abimem_init [ Functions ]

[ Top ] [ m_profiling_abi ] [ Functions ]

NAME

 abimem_init

FUNCTION

  Initialize memory profiling module.

 INPUT
  level = Integer selecting the operation mode:
       0 -> no file abimem.mocc is created, only memory allocation counters running
       1 -> light version. Only memory peaks are written.
       2 -> file abimem.mocc is created with full information inside.
       3 -> Write info only if allocation/deallocation is larger or smaller than limit_mb
                depending on of the sign of limit_mb
    NOTE: By default, only master node writes, use negative values to make all MPI procs write info to disk.
  [delta_time]=Interval in second for snapshots. Will write report to std_out every delta_time seconds.
  [filename] = If present, activate memory logging only inside filename (basename).
  [limit_mb]= Set memory limit in Mb if level == 3. Print allocation/deallocation only above this limit.
    Positive value to print above the threshold.
    Negative value to print below the threshold.

m_profiling_abi/abimem_record [ Functions ]

[ Top ] [ m_profiling_abi ] [ Functions ]

NAME

 abimem_record

FUNCTION

  Control the memory occupation by calculating the overall size of the allocated arrays
  At the end of the calculation a short report is printed on the screen,
  some information can be also written on disk following the needs

SOURCE

395 subroutine abimem_record(istat, vname, addr, act, isize, file, line)
396 
397 !Arguments ------------------------------------
398  integer,intent(in) :: istat,line
399  integer(i8b), intent(in) :: isize,addr
400  character(len=*), intent(in) :: vname,act,file
401 
402 !Local variables-------------------------------
403  !integer :: ierr
404  !real(dp) :: now
405  logical :: do_log, new_peak
406  character(len=500) :: msg
407 ! *************************************************************************
408 
409  ! Handle possible allocate/deallocate failures.
410  if (istat /= 0) then
411    if (isize >= 0) then
412      write(msg,*)" Problem of allocation of variable: ",trim(vname),', error code= ',istat
413      _ABORT(msg)
414    else if (isize<0) then
415      write(msg,*)" Problem of deallocation of variable ",trim(vname),', error code= ',istat
416      _ABORT(msg)
417    end if
418  end if
419 
420  ! Increase total counter
421  minfo%memory = minfo%memory + isize
422  new_peak = .False.
423  if (isize > minfo%peak) then
424    ! New peak
425    new_peak = .True.
426    minfo%peak = isize
427    minfo%peak_vname = vname
428    minfo%peak_file = file
429    minfo%peak_fileline = line
430  end if
431 
432  if (isize > 0) then
433    minfo%num_alloc = minfo%num_alloc + 1
434  else if (isize < 0) then
435    minfo%num_free = minfo%num_free + 1
436  else
437    ! This is the correct check but tests fail!
438    !if (act == "A") then
439    !  minfo%num_alloc = minfo%num_alloc + 1
440    !else if (act == "D") then
441    !  minfo%num_free = minfo%num_free + 1
442    !else
443    !  _ABORT("Wrong action: "//trim(act))
444    !end if
445  end if
446 
447  ! Selective memory tracing
448  do_log = .True.
449  if (minfo%select_file /= NONE_STRING) do_log = (minfo%select_file == file)
450  do_log = do_log .and. minfo%iwrite
451 
452  ! Snapshot (write to std_out)
453  !if (do_log .and. minfo%last_snapshot >= zero) then
454  !  now = abimem_wtime()
455  !  if (now - minfo%last_snapshot >= minfo%dt_snapshot) then
456  !    call abimem_report("", std_out)
457  !    minfo%last_snapshot = now
458  !  end if
459  !end if
460 
461  ! IMPORTANT:
462  ! Remember to change the pyton code in ~abinit/tests/pymods/memprof.py to account for changes in the format
463  if (do_log) then
464    select case (minfo%level)
465    case (0)
466      ! No action required
467 
468    case (1)
469      ! Write only if we have a new peak
470      if (new_peak) then
471        write(minfo%logunt,'(a,t60,a,1x,2(i0,1x),a,1x,2(i0,1x))') &
472          trim(vname), trim(act), addr, isize, trim(abimem_basename(file)), line, minfo%memory
473      end if
474 
475    case (2)
476      ! To be used for inspecting a variable which is not deallocated
477      write(minfo%logunt,'(a,t60,a,1x,2(i0,1x),a,1x,2(i0,1x))') &
478        trim(vname), trim(act), addr, isize, trim(abimem_basename(file)), line, minfo%memory
479 
480    case (3)
481      ! Write memory allocations larger than limit_mb
482      if ((abs(isize * b2Mb) > minfo%limit_mb .and. minfo%limit_mb > zero) .or. &
483          (abs(isize * b2Mb) < minfo%limit_mb .and. minfo%limit_mb < zero)) then
484        write(minfo%logunt,'(a,t60,a,1x,2(i0,1x),a,1x,2(i0,1x))') &
485          trim(vname), trim(act), addr, isize, trim(abimem_basename(file)), line, minfo%memory
486      end if
487 
488    case default
489      _ABORT("Invalid abimem_level")
490    end select
491  end if
492 
493 end subroutine abimem_record

m_profiling_abi/abimem_report [ Functions ]

[ Top ] [ m_profiling_abi ] [ Functions ]

NAME

 abimem_report

FUNCTION

  Print info about memory usage to unit `unt`.
  Add mallinfo values if `with_mallinfo` (default: True)

 INPUT

m_profiling_abi/abimem_set_snapshot_time [ Functions ]

[ Top ] [ m_profiling_abi ] [ Functions ]

NAME

 abimem_set_snapshot_time

FUNCTION

  Set time interface for snapshots.

 INPUT
  delta_time=Interval in second for snapshots.

 INPUT

m_profiling_abi/abimem_shutdown [ Functions ]

[ Top ] [ m_profiling_abi ] [ Functions ]

NAME

 abimem_shutdown

FUNCTION

 Perform final cleanup of the module and close files.

 INPUT

m_profiling_abi/abimem_t [ Types ]

[ Top ] [ m_profiling_abi ] [ Types ]

NAME

  minfo_t

FUNCTION

  Internal datastructure storing information on the memory allocated at run-time

SOURCE

 68  type :: abimem_t
 69 
 70    integer :: level = huge(1)
 71    ! Integer selecting the operation mode
 72    ! The initial value is set to huge so that main executables that don't call abimem_init will
 73    ! produce an Error when the first alocation is performed and abimem_record is invoked.
 74 
 75    integer(i8b) :: memory = 0
 76    ! Total memory allocated so far in bytes.
 77 
 78    integer(i8b) :: peak = 0
 79    ! Memory peak in bytes.
 80 
 81    integer :: peak_fileline = -1
 82    ! Line in peak_file
 83 
 84    integer(i8b) :: num_alloc = 0
 85    ! Total numer of allocations performed so far.
 86 
 87    integer(i8b) :: num_free = 0
 88    ! Total numer of deallocations performed so far.
 89 
 90    integer :: logunt = 99
 91    ! Unit number of logfile (hardcoded)
 92 
 93    integer :: my_rank = 0
 94    ! Rank of this processor.
 95 
 96    logical :: iwrite = .False.
 97    ! True if this MPI rank should write to file.
 98 
 99    !real(dp),private,save :: start_time
100    ! Origin of time in seconds.
101 
102    real(dp) :: last_snapshot = -one
103    ! time of the last snapshot in seconds.
104 
105    real(dp) :: dt_snapshot = -one
106    ! time between two consecutive snapshots in seconds.
107 
108    real(dp) :: limit_mb = 20_dp
109    ! Optional memory limit in Mb. used when level == 3
110 
111    character(len=slen) :: peak_vname = "_vname"
112    ! Name of the last variable for which the memory peak occurred.
113 
114    character(len=slen) :: peak_file = NONE_STRING
115    ! Name of the file in which peak occurred.
116 
117    ! Selective memory tracing
118    character(fnlen) :: select_file = NONE_STRING
119 
120    character(len=fnlen) :: logfile
121    ! File used for logging allocations/deallocations.
122 
123  end type abimem_t

m_time/abimem_wtime [ Functions ]

[ Top ] [ m_time ] [ Functions ]

NAME

  abimem_wtime

FUNCTION

  Return wall clock time in seconds since some arbitrary start.
  Call the F90 intrinsic date_and_time .

INPUTS

  (no inputs)

OUTPUT

  wall= wall clock time in seconds

SOURCE

612 function abimem_wtime() result(wall)
613 
614 !Arguments ------------------------------------
615 !scalars
616  real(dp) :: wall
617 
618 !Local variables-------------------------------
619 !scalars
620 #ifndef HAVE_MPI
621  integer,parameter :: nday(24)=(/31,28,31,30,31,30,31,31,30,31,30,31,&
622 &                                31,28,31,30,31,30,31,31,30,31,30,31/)
623  integer,save :: month_init,month_now,start=1,year_init
624  integer :: months
625  character(len=8)   :: date
626  character(len=10)  :: time
627  character(len=5)   :: zone
628 !arrays
629  integer :: values(8)
630 #endif
631 
632 ! *************************************************************************
633 
634 #ifndef HAVE_MPI
635 
636 !The following section of code is standard F90, but it is useful only if the intrinsics
637 !date_and_time is accurate at the 0.01 sec level, which is not the case for a P6 with the pghpf compiler ...
638 !Year and month initialisation
639  if(start==1)then
640    start=0
641    call date_and_time(date,time,zone,values)
642    year_init=values(1)
643    month_init=values(2)
644  end if
645 
646 !Uses intrinsic F90 subroutine Date_and_time for
647 !wall clock (not correct when a change of year happen)
648  call date_and_time(date,time,zone,values)
649 
650 !Compute first the number of seconds from the beginning of the month
651  wall=(values(3)*24.0d0+values(5))*3600.0d0+values(6)*60.0d0+values(7)+values(8)*0.001d0
652 
653 !If the month has changed, compute the number of seconds
654 !to be added. This fails if the program ran one year !!
655  month_now=values(2)
656  if(month_now/=month_init)then
657    if(year_init+1==values(1))then
658      month_now=month_now+12
659    end if
660    if(month_now<=month_init)then
661      _ABORT('Problem with month and year numbers.')
662    end if
663    do months=month_init,month_now-1
664      wall=wall+86400.0d0*nday(months)
665    end do
666  end if
667 
668 !Now take into account bissextile years (I think 2000 is bissextile, but I am not sure ...)
669  if(mod(year_init,4)==0 .and. month_init<=2 .and. month_now>2)   wall=wall+3600.0d0
670  if(mod(values(1),4)==0 .and. month_init<=14 .and. month_now>14) wall=wall+3600.0d0
671 
672 #else
673 !Use the timer provided by MPI1.
674  wall = MPI_WTIME()
675 #endif
676 
677 end function abimem_wtime