# http://pyrocko.org - GPLv3 # # The Pyrocko Developers, 21st Century # ---|P------/S----------~Lg----------
'displacement': 'Displacement [m]', 'velocity': 'Velocity [m/s]', 'acceleration': 'Acceleration [m/s²]' }
store_id, source, radius, azi_begin=0., azi_end=360., dazi=1., interpolation='multilinear', components='RTZ', quantity='displacement'):
quantity='displacement', interpolation=interpolation, store_id=store_id)
Target( lat=source.lat, lon=source.lon, north_shift=coords[0, iazi] + source.north_shift, east_shift=coords[1, iazi] + source.east_shift, azimuth={ 'R': azi, 'T': azi+90., 'Z': 0. }[channel], dip=dips[channel], codes=('', 'S%01d' % iazi, '', channel), **target_kwargs) for iazi, azi in enumerate(azimuths) for channel in components]
response, fmin=None, fmax=None, component='R', envelope=False):
continue
elif fmin: tr.highpass(4, fmin) elif fmax: tr.lowpass(4, fmax)
azimuth = 360.0 - azimuth azi = azimuth * r2d alt = angle_altitude * r2d
x, y = num.gradient(array) slope = num.pi/2. - num.arctan(num.sqrt(x*x + y*y)) aspect = num.arctan2(-x, y)
shaded = num.sin(alt)*num.sin(slope) \ + num.cos(alt)*num.cos(slope)*num.cos((azi - num.pi/2.) - aspect)
return (shaded + 1.)/2.
seismogram_array, rgba_map, shad_lim=(.4, .98), contrast=1., blend_mode='multiply'): # Light source from somewhere above - psychologically the best choice # from upper left
# convolution of two 2D arrays
# if there are strong artifical edges in the data, shades get # dominated by them. Cutting off the largest and smallest 2% of # # shades helps
# # normalize shading
# # reduce range to balance gray color
rgba_map[:, :3] = 1. - ((1. - rgba_map[:, :3])*(shad[:, num.newaxis]))
engine, source, store_id, distance=300*km, azi_begin=0., azi_end=360., dazi=1., phases={'P': 'first{stored:any_P}-10%', 'S': 'last{stored:any_S}+50'}, quantity='displacement', envelope=False, component='R', fmin=0.01, fmax=0.1, hillshade=True, cmap=None, plot_mt='full', show_phases=True, show_description=True, reverse_time=False, show_nucleations=True, axes=None, nthreads=0): ''' Plot the directivity and radiation characteristics of source models.
Synthetic seismic traces (R, T or Z) are forward-modelled at a defined radius, covering the full or partial azimuthal range and projected on a polar plot. Difference in the amplitude are enhanced by hillshading the data.
:param engine: Forward modelling engine :type engine: :py:class:`~pyrocko.gf.seismosizer.Engine` :param source: Parametrized source model :type source: :py:class:`~pyrocko.gf.seismosizer.Source` :param store_id: Store ID used for forward modelling :type store_id: str :param distance: Distance in [m] :type distance: float :param azi_begin: Begin azimuth in [deg] :type azi_begin: float :param azi_end: End azimuth in [deg] :type azi_end: float :param dazi: Delta azimuth, bin size [deg] :type dazi: float :param phase_begin: Start time of the window :type phase_begin: :py:class:`~pyrocko.gf.meta.Timing` :param phase_end: End time of the window :type phase_end: :py:class:`~pyrocko.gf.meta.Timing` :param quantity: Seismogram quantity, default ``displacement`` :type quantity: str :param envelope: Plot envelop instead of seismic trace :type envelope: bool :param component: Forward modelled component, default ``R``. Choose from `RTZ` :type component: str :param fmin: Bandpass lower frequency [Hz], default ``0.01`` :type fmin: float :param fmax: Bandpass upper frequency [Hz], default ``0.1`` :type fmax: float :param hillshade: Enable hillshading, default ``True`` :type hillshade: bool :param cmap: Matplotlit colormap to use, default ``seismic``. When ``envelope`` is ``True`` the default colormap will be ``Reds``. :type cmap: str :param plot_mt: Plot a centered moment tensor, default ``full``. Choose from ``full, deviatoric, dc or False`` :type plot_mt: str, bool :param show_phases: Show annotations, default ``True`` :type show_phases: bool :param show_description: Show desciption, default ``True`` :type show_description: bool :param reverse_time: Reverse time axis. First phases arrive at the center, default ``False`` :type reverse_time: bool :param show_nucleations: Show nucleation piercing points on the moment tensor, default ``True`` :type show_nucleations: bool :param axes: Give axes to plot into :type axes: :py:class:`matplotlib.axes.Axes` :param nthreads: Number of threads used for forward modelling, default ``0`` - all available cores :type nthreads: int '''
else: fig = axes.figure ax = axes
elif cmap is None: cmap = 'seismic'
store_id, source, distance, azi_begin, azi_end, dazi, components='R', quantity=quantity)
resp, fmin, fmax, component=component, envelope=envelope)
nx = float(source.nucleation_x[0]) ny = float(source.nucleation_y[0])
except AttributeError: logger.warn('Old matplotlib version: cannot set label positions')
if v < tbegin or v > tend: return '' return '%g s' % v
ax.set_rlim(times[0] - .3*duration, times[-1]) else:
mt, ax, beachball_type=plot_mt, size=mt_size, size_units='axes', color_t=(0.7, 0.4, 0.4), position=(.5, .5), linewidth=1.)
and show_nucleations: nucleation_x = source.nucleation_x nucleation_y = source.nucleation_y
transform=ax.transAxes, zorder=10)
azimuths * d2r, times, data, cmap=cmw.cmap, shading='gouraud', zorder=0)
data, color, shad_lim=(.85, 1.), blend_mode='multiply')
label_theta*d2r, time, label, ha='left', color='k', fontsize='small')
'Component {component:s}\n' 'Distance {distance:g} km').format( component=component, distance=distance / km)
fmin=fmin, fmax=fmax) elif fmin: description += '\nHighpass {fmin:g} Hz'.format(fmin=fmin) elif fmax: description += '\nLowpass {fmax:g} Hz'.format(fmax=fmax) -.05, -.05, description, fontsize='small', ha='left', va='bottom', transform=ax.transAxes)
cmw, ax=ax, orientation='vertical', shrink=.8, pad=0.11)
if __name__ == '__main__': engine = LocalEngine(store_superdirs=['.'], use_config=True)
rect_source = RectangularSource( depth=2.6*km, strike=240., dip=76.6, rake=-.4, anchor='top',
nucleation_x=-.57, nucleation_y=-.59, velocity=2070.,
length=27*km, width=9.4*km, slip=1.4)
resp = plot_directivity( engine, rect_source, 'crust2_ib', dazi=5, component='R', quantity='displacement', envelope=True) |