1 !***********************************************************************
2 !* GNU Lesser General Public License
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 !> @defgroup mpp_memutils_mod mpp_memutils_mod
21 !> @brief Routines to initialize and report on memory usage during the model run.
23 !> @addtogroup mpp_memutils_mod
25 module mpp_memutils_mod
27 use mpp_mod, only: mpp_min, mpp_max, mpp_sum, mpp_pe, mpp_root_pe
28 use mpp_mod, only: mpp_error, FATAL, stderr, mpp_npes
34 public :: mpp_print_memuse_stats, mpp_mem_dump
35 public :: mpp_memuse_begin, mpp_memuse_end
38 logical :: memuse_started = .false.
42 !#######################################################################
43 !> Initialize the memory module, and record the initial memory use.
44 subroutine mpp_memuse_begin
45 if(memuse_started) then
46 call mpp_error(FATAL, "mpp_memutils_mod: mpp_memuse_begin was already called")
48 memuse_started = .true.
50 call mpp_mem_dump(begin_memuse)
51 end subroutine mpp_memuse_begin
53 !#######################################################################
54 !> End the memory collection, and report on total memory used during the
55 !! execution of the model run.
56 subroutine mpp_memuse_end(text, unit)
57 character(len=*), intent(in) :: text !< Text to include in memory use statement
58 integer, intent(in), optional :: unit !< Fortran unit number where memory report should go.
60 real :: m, mmin, mmax, mavg, mstd, end_memuse
63 if(.NOT.memuse_started) then
64 call mpp_error(FATAL, "mpp_memutils_mod: mpp_memuse_begin must be called before calling mpp_memuse_being")
66 memuse_started = .false.
68 call mpp_mem_dump(end_memuse)
70 mu = stderr(); if( PRESENT(unit) )mu = unit
71 m = end_memuse - begin_memuse
72 mmin = m; call mpp_min(mmin)
73 mmax = m; call mpp_max(mmax)
74 mavg = m; call mpp_sum(mavg); mavg = mavg/mpp_npes()
75 mstd = (m-mavg)**2; call mpp_sum(mstd); mstd = sqrt( mstd/mpp_npes() )
76 if( mpp_pe().EQ.mpp_root_pe() )write( mu,'(a64,4es11.3)' ) &
77 'Memory(MB) used in '//trim(text)//'=', mmin, mmax, mstd, mavg
80 end subroutine mpp_memuse_end
82 !#######################################################################
83 !> Print the current memory high water mark to stderr, or the unit
85 subroutine mpp_print_memuse_stats(text, unit)
86 character(len=*), intent(in) :: text !< Text to include in memory print statement
87 integer, intent(in), optional :: unit !< Fortran unit number where print statement should go.
89 real :: m, mmin, mmax, mavg, mstd
92 mu = stderr(); if( PRESENT(unit) )mu = unit
95 mmin = m; call mpp_min(mmin)
96 mmax = m; call mpp_max(mmax)
97 mavg = m; call mpp_sum(mavg); mavg = mavg/mpp_npes()
98 mstd = (m-mavg)**2; call mpp_sum(mstd); mstd = sqrt( mstd/mpp_npes() )
99 if( mpp_pe().EQ.mpp_root_pe() )write( mu,'(a64,4es11.3)' ) &
100 'Memuse(MB) at '//trim(text)//'=', mmin, mmax, mstd, mavg
103 end subroutine mpp_print_memuse_stats
105 !#######################################################################
107 !> \brief Return the memory high water mark in MiB
109 !! Query the system for the memory high water mark, return the result in MiB.
110 subroutine mpp_mem_dump(memuse)
111 real, intent(out) :: memuse !< Memory, high water mark, in MiB
114 integer(KIND=c_size_t) function getpeakrss() bind(c, name="getpeakrss")
115 use, intrinsic :: iso_c_binding
116 end function getpeakrss
119 ! Value of Bytes to Mebibytes
120 real, parameter :: B_to_MiB = 1048576.0
122 ! Get the max memory use, convert to MiB
123 memuse = real(getpeakrss())/B_to_MiB
126 end subroutine mpp_mem_dump
127 end module mpp_memutils_mod
129 ! close documentation grouping