Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

114

115

116

117

118

119

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

145

146

147

148

149

150

151

152

153

154

155

156

157

158

159

160

161

162

163

164

165

166

167

168

169

170

171

172

173

174

175

176

177

178

179

180

181

182

183

184

185

186

187

188

189

190

191

192

193

194

195

196

197

198

199

200

201

202

203

204

205

206

207

208

209

210

211

212

213

214

215

216

217

218

219

220

221

# -*- coding: utf-8 -*- 

#------------------------------------------------------------------- 

# Filename: spectrogram.py 

#  Purpose: Plotting spectrogram of Seismograms. 

#   Author: Christian Sippl, Moritz Beyreuther 

#    Email: sippl@geophysik.uni-muenchen.de 

# 

# Copyright (C) 2008-2012 Christian Sippl 

#--------------------------------------------------------------------- 

""" 

Plotting spectrogram of seismograms. 

 

:copyright: 

    The ObsPy Development Team (devs@obspy.org) 

:license: 

    GNU General Public License (GPL) 

    (http://www.gnu.org/licenses/gpl.txt) 

""" 

 

from matplotlib import mlab 

from matplotlib.colors import Normalize 

from obspy.core.util import getMatplotlibVersion 

import math as M 

import matplotlib.pyplot as plt 

import numpy as np 

 

 

MATPLOTLIB_VERSION = getMatplotlibVersion() 

 

 

def _nearestPow2(x): 

    """ 

    Find power of two nearest to x 

 

    >>> _nearestPow2(3) 

    2.0 

    >>> _nearestPow2(15) 

    16.0 

 

    :type x: Float 

    :param x: Number 

    :rtype: Int 

    :return: Nearest power of 2 to x 

    """ 

    a = M.pow(2, M.ceil(np.log2(x))) 

    b = M.pow(2, M.floor(np.log2(x))) 

    if abs(a - x) < abs(b - x): 

        return a 

    else: 

        return b 

 

 

def spectrogram(data, samp_rate, per_lap=0.9, wlen=None, log=False, 

                outfile=None, fmt=None, axes=None, dbscale=False, 

                mult=8.0, cmap=None, zorder=None, title=None, show=True, 

                sphinx=False, clip=[0.0, 1.0]): 

    """ 

    Computes and plots spectrogram of the input data. 

 

    :param data: Input data 

    :type samp_rate: float 

    :param samp_rate: Samplerate in Hz 

    :type per_lap: float 

    :param per_lap: Percentage of overlap of sliding window, ranging from 0 

        to 1. High overlaps take a long time to compute. 

    :type wlen: int or float 

    :param wlen: Window length for fft in seconds. If this parameter is too 

        small, the calculation will take forever. 

    :type log: bool 

    :param log: Logarithmic frequency axis if True, linear frequency axis 

        otherwise. 

    :type outfile: String 

    :param outfile: String for the filename of output file, if None 

        interactive plotting is activated. 

    :type fmt: String 

    :param fmt: Format of image to save 

    :type axes: :class:`matplotlib.axes.Axes` 

    :param axes: Plot into given axes, this deactivates the fmt and 

        outfile option. 

    :type dbscale: bool 

    :param dbscale: If True 10 * log10 of color values is taken, if False the 

        sqrt is taken. 

    :type mult: float 

    :param mult: Pad zeros to lengh mult * wlen. This will make the spectrogram 

        smoother. Available for matplotlib > 0.99.0. 

    :type cmap: :class:`matplotlib.colors.Colormap` 

    :param cmap: Specify a custom colormap instance 

    :type zorder: float 

    :param zorder: Specify the zorder of the plot. Only of importance if other 

        plots in the same axes are executed. 

    :type title: String 

    :param title: Set the plot title 

    :type show: bool 

    :param show: Do not call `plt.show()` at end of routine. That way, further 

        modifications can be done to the figure before showing it. 

    :type sphinx: bool 

    :param sphinx: Internal flag used for API doc generation, default False 

    :type clip: [float, float] 

    :param clip: adjust colormap to clip at lower and/or upper end. The given 

        percentages of the amplitude range (linear or logarithmic depending 

        on option `dbscale`) are clipped. 

    """ 

    # enforce float for samp_rate 

    samp_rate = float(samp_rate) 

 

    # set wlen from samp_rate if not specified otherwise 

    if not wlen: 

        wlen = samp_rate / 100. 

 

    npts = len(data) 

    # nfft needs to be an integer, otherwise a deprecation will be raised 

    #XXX add condition for too many windows => calculation takes for ever 

    nfft = int(_nearestPow2(wlen * samp_rate)) 

    if nfft > npts: 

        nfft = int(_nearestPow2(npts / 8.0)) 

 

    if mult != None: 

        mult = int(_nearestPow2(mult)) 

        mult = mult * nfft 

    nlap = int(nfft * float(per_lap)) 

 

    data = data - data.mean() 

    end = npts / samp_rate 

 

    # Here we call not plt.specgram as this already produces a plot 

    # matplotlib.mlab.specgram should be faster as it computes only the 

    # arrays 

    # XXX mlab.specgram uses fft, would be better and faster use rfft 

    if MATPLOTLIB_VERSION >= [0, 99, 0]: 

        specgram, freq, time = mlab.specgram(data, Fs=samp_rate, NFFT=nfft, 

                                              pad_to=mult, noverlap=nlap) 

    else: 

        specgram, freq, time = mlab.specgram(data, Fs=samp_rate, 

                                                NFFT=nfft, noverlap=nlap) 

    # db scale and remove zero/offset for amplitude 

    if dbscale: 

        specgram = 10 * np.log10(specgram[1:, :]) 

    else: 

        specgram = np.sqrt(specgram[1:, :]) 

    freq = freq[1:] 

 

    vmin, vmax = clip 

    if vmin < 0 or vmax > 1 or vmin >= vmax: 

        msg = "Invalid parameters for clip option." 

        raise ValueError(msg) 

    _range = float(specgram.max() - specgram.min()) 

    vmin = specgram.min() + vmin * _range 

    vmax = specgram.min() + vmax * _range 

    norm = Normalize(vmin, vmax, clip=True) 

 

    if not axes: 

        fig = plt.figure() 

        ax = fig.add_subplot(111) 

    else: 

        ax = axes 

 

    # calculate half bin width 

    halfbin_time = (time[1] - time[0]) / 2.0 

    halfbin_freq = (freq[1] - freq[0]) / 2.0 

 

    if log: 

        # pcolor expects one bin more at the right end 

        freq = np.concatenate((freq, [freq[-1] + 2 * halfbin_freq])) 

        time = np.concatenate((time, [time[-1] + 2 * halfbin_time])) 

        # center bin 

        time -= halfbin_time 

        freq -= halfbin_freq 

        # pcolormesh issue was fixed in matplotlib r5716 (2008-07-07) 

        # inbetween tags 0.98.2 and 0.98.3 

        # see: 

        #  - http://matplotlib.svn.sourceforge.net/viewvc/... 

        #    matplotlib?revision=5716&view=revision 

        #  - http://matplotlib.sourceforge.net/_static/CHANGELOG 

        if MATPLOTLIB_VERSION >= [0, 98, 3]: 

            # Log scaling for frequency values (y-axis) 

            ax.set_yscale('log') 

            # Plot times 

            ax.pcolormesh(time, freq, specgram, cmap=cmap, zorder=zorder, 

                          norm=norm) 

        else: 

            X, Y = np.meshgrid(time, freq) 

            ax.pcolor(X, Y, specgram, cmap=cmap, zorder=zorder, norm=norm) 

            ax.semilogy() 

    else: 

        # this method is much much faster! 

        specgram = np.flipud(specgram) 

        # center bin 

        extent = (time[0] - halfbin_time, time[-1] + halfbin_time, \ 

                  freq[0] - halfbin_freq, freq[-1] + halfbin_freq) 

        ax.imshow(specgram, interpolation="nearest", extent=extent, 

                  cmap=cmap, zorder=zorder) 

 

    # set correct way of axis, whitespace before and after with window 

    # length 

    ax.axis('tight') 

    ax.set_xlim(0, end) 

    ax.grid(False) 

 

    if axes: 

        return ax 

 

    ax.set_xlabel('Time [s]') 

    ax.set_ylabel('Frequency [Hz]') 

    if title: 

        ax.set_title(title) 

 

    if not sphinx: 

        # ignoring all NumPy warnings during plot 

        temp = np.geterr() 

        np.seterr(all='ignore') 

        plt.draw() 

        np.seterr(**temp) 

    if outfile: 

        if fmt: 

            fig.savefig(outfile, format=fmt) 

        else: 

            fig.savefig(outfile) 

    elif show: 

        plt.show() 

    else: 

        return fig