257 lines
8.6 KiB
Plaintext
257 lines
8.6 KiB
Plaintext
|
{
|
||
|
"cells": [
|
||
|
{
|
||
|
"cell_type": "markdown",
|
||
|
"metadata": {},
|
||
|
"source": [
|
||
|
"# Surface Wave Velocity (Single Station Approach)"
|
||
|
]
|
||
|
},
|
||
|
{
|
||
|
"cell_type": "code",
|
||
|
"execution_count": null,
|
||
|
"metadata": {
|
||
|
"jupyter": {
|
||
|
"outputs_hidden": false
|
||
|
}
|
||
|
},
|
||
|
"outputs": [],
|
||
|
"source": [
|
||
|
"# Preparation: load packages, set some basic options \n",
|
||
|
"%matplotlib inline\n",
|
||
|
"from obspy import UTCDateTime\n",
|
||
|
"from obspy.clients.fdsn import Client\n",
|
||
|
"from obspy.clients.fdsn import URL_MAPPINGS\n",
|
||
|
"import numpy as np\n",
|
||
|
"import matplotlib.pylab as plt\n",
|
||
|
"from movingWindow import movingWindowAnalysis, hann # import the moving window analysis routines\n",
|
||
|
"from multipleFilter import multipleFilterAnalysis, normalizeMFT # import the multiple filter analysis routines\n",
|
||
|
"plt.style.use('ggplot')\n",
|
||
|
"plt.rcParams['figure.figsize'] = 15, 4\n",
|
||
|
"plt.rcParams['lines.linewidth'] = 0.5"
|
||
|
]
|
||
|
},
|
||
|
{
|
||
|
"cell_type": "markdown",
|
||
|
"metadata": {},
|
||
|
"source": [
|
||
|
"## Data from data archive using obspy webfdsn client\n",
|
||
|
"FDSN stands for the International Federal of Digital Seismic Networks (www.fdsn.org). Citation from their web site: \"The International Federation of Digital Seismograph Networks (FDSN) is a global organization. Its membership is comprised of groups responsible for the installation and maintenance of seismographs either within their geographic borders or globally. Membership in the FDSN is open to all organizations that operate more than one broadband station. Members agree to coordinate station siting and provide free and open access to their data. This cooperation helps scientists all over the world to further the advancement of earth science and particularly the study of global seismic activity.\"\n",
|
||
|
"BGR (Bundesanstalt für Geowissenschaften und Rohstoffe) operates broadband stations in Germany and operates a data archive. It is member of the FDSN. That is why we can freely download data from them.\n",
|
||
|
"There are other data archives around the world which are also member of the FDSN and allow opn access to their data (see OBSPY documentation "
|
||
|
]
|
||
|
},
|
||
|
{
|
||
|
"cell_type": "code",
|
||
|
"execution_count": null,
|
||
|
"metadata": {
|
||
|
"jupyter": {
|
||
|
"outputs_hidden": false
|
||
|
}
|
||
|
},
|
||
|
"outputs": [],
|
||
|
"source": [
|
||
|
"for key in sorted(URL_MAPPINGS.keys()): # eine Liste der Archive\n",
|
||
|
" print(\"{0:<7} {1}\".format(key, URL_MAPPINGS[key]))"
|
||
|
]
|
||
|
},
|
||
|
{
|
||
|
"cell_type": "markdown",
|
||
|
"metadata": {},
|
||
|
"source": [
|
||
|
"Get waveforms for the Tohoku earthquake for station TNS of the German Regional Seismic Network."
|
||
|
]
|
||
|
},
|
||
|
{
|
||
|
"cell_type": "code",
|
||
|
"execution_count": null,
|
||
|
"metadata": {
|
||
|
"jupyter": {
|
||
|
"outputs_hidden": false
|
||
|
}
|
||
|
},
|
||
|
"outputs": [],
|
||
|
"source": [
|
||
|
"client = Client(\"BGR\") # data archive from which data are fetched\n",
|
||
|
"t1 = UTCDateTime(\"2018-01-14T09:19:00.000\") # desired start time of data\n",
|
||
|
"stc = client.get_waveforms(\"GR\", \"TNS\", \"\", \"LHZ\", t1, # use fdsn web client to download data\n",
|
||
|
" t1 + 7200., attach_response=True)"
|
||
|
]
|
||
|
},
|
||
|
{
|
||
|
"cell_type": "code",
|
||
|
"execution_count": null,
|
||
|
"metadata": {
|
||
|
"jupyter": {
|
||
|
"outputs_hidden": false
|
||
|
}
|
||
|
},
|
||
|
"outputs": [],
|
||
|
"source": [
|
||
|
"stc.detrend('linear') # take out linear trend\n",
|
||
|
"stc.plot()\n",
|
||
|
"print(stc)"
|
||
|
]
|
||
|
},
|
||
|
{
|
||
|
"cell_type": "code",
|
||
|
"execution_count": null,
|
||
|
"metadata": {
|
||
|
"jupyter": {
|
||
|
"outputs_hidden": false
|
||
|
}
|
||
|
},
|
||
|
"outputs": [],
|
||
|
"source": [
|
||
|
"st = stc.copy() # proceed with copy to avoid repeated downloading\n",
|
||
|
"ta = UTCDateTime(\"2018-01-14T10:05:00.000000Z\") # cut to surface wave train\n",
|
||
|
"te = UTCDateTime(\"2018-01-14T10:22:00.000000Z\")\n",
|
||
|
"st.trim(ta,te)\n",
|
||
|
"st.taper(0.05,\"hann\") # apply a taper to bring signal to zero at ends\n",
|
||
|
"st.plot() # dispersion is now well visible\n",
|
||
|
"print(st)\n",
|
||
|
"tr = st[0]"
|
||
|
]
|
||
|
},
|
||
|
{
|
||
|
"cell_type": "markdown",
|
||
|
"metadata": {},
|
||
|
"source": [
|
||
|
"## Moving window approach"
|
||
|
]
|
||
|
},
|
||
|
{
|
||
|
"cell_type": "code",
|
||
|
"execution_count": null,
|
||
|
"metadata": {
|
||
|
"jupyter": {
|
||
|
"outputs_hidden": false
|
||
|
}
|
||
|
},
|
||
|
"outputs": [],
|
||
|
"source": [
|
||
|
"wlen = 256.0 # window length in sec\n",
|
||
|
"npts = tr.stats.npts # number of samples in the trace\n",
|
||
|
"dt = tr.stats.delta # sample interval\n",
|
||
|
"nwin = int(wlen/dt) # nuber of samples in window\n",
|
||
|
"nwin = int(pow(2, np.ceil(np.log(nwin)/np.log(2)))) # find next higher power of 2 of nwin\n",
|
||
|
"shift = nwin//8\n",
|
||
|
"tshift = shift*dt\n",
|
||
|
"nseg,mwa = movingWindowAnalysis(tr.data,hann,nwin,shift,2.0) # compute spectrogram, exponent=2.0\n",
|
||
|
"freq = np.fft.rfftfreq(nwin,dt) # Fourier frequencies\n",
|
||
|
"print(\"Frequency range [Hz]: \",freq[0],freq[-1])\n",
|
||
|
"print(\"Number of frequency samples: \",len(freq))\n",
|
||
|
"print(\"Window length in samples: \",nwin)\n",
|
||
|
"print(\"Window shift in samples: \",shift)\n",
|
||
|
"print(\"Number of segments: \",nseg)"
|
||
|
]
|
||
|
},
|
||
|
{
|
||
|
"cell_type": "code",
|
||
|
"execution_count": null,
|
||
|
"metadata": {
|
||
|
"jupyter": {
|
||
|
"outputs_hidden": false
|
||
|
}
|
||
|
},
|
||
|
"outputs": [],
|
||
|
"source": [
|
||
|
"# plot the result\n",
|
||
|
"st.plot() # plot seismogram again\n",
|
||
|
"fmax = 0.08 # maximum frequency to be plotted\n",
|
||
|
"fmin = 0.02\n",
|
||
|
"nf = len(freq)-1 # number of frequencies\n",
|
||
|
"jfmax = int(fmax/freq[-1]*nf) # max index for plotting \n",
|
||
|
"jfmin = int(fmin/freq[-1]*nf) # min index for plotting \n",
|
||
|
"extent = (0.5*wlen,0.5*wlen+nseg*tshift,fmin,fmax) # extent of matrix in true time and frequency\n",
|
||
|
"asp = 0.5*nseg*tshift/(fmax-fmin)\n",
|
||
|
"fg = plt.figure(figsize = [12.5,12])\n",
|
||
|
"plt.imshow(mwa[jfmin:jfmax,:],extent = extent, aspect = asp, \n",
|
||
|
" origin = \"lower\", interpolation = \"bicubic\") # plotting with interpolation\n",
|
||
|
"plt.xlabel('Window center time [s]')\n",
|
||
|
"plt.ylabel('Frequency [Hz]')\n",
|
||
|
"plt.show()"
|
||
|
]
|
||
|
},
|
||
|
{
|
||
|
"cell_type": "markdown",
|
||
|
"metadata": {},
|
||
|
"source": [
|
||
|
"Observe the change of frequency content with time from 0.03 Hz initially to 0.06 Hz towards the end of the record. Given a epicentral distance $d$, we could transform the time axis into a group velocity axis via $U=d/t$."
|
||
|
]
|
||
|
},
|
||
|
{
|
||
|
"cell_type": "markdown",
|
||
|
"metadata": {},
|
||
|
"source": [
|
||
|
"## Multiple filtering approach"
|
||
|
]
|
||
|
},
|
||
|
{
|
||
|
"cell_type": "code",
|
||
|
"execution_count": null,
|
||
|
"metadata": {
|
||
|
"jupyter": {
|
||
|
"outputs_hidden": false
|
||
|
}
|
||
|
},
|
||
|
"outputs": [],
|
||
|
"source": [
|
||
|
"alfa = 50.0 # filter width parameter\n",
|
||
|
"cfreq = np.linspace(fmin,fmax,100) # array of filter center frequencies\n",
|
||
|
"mfa = multipleFilterAnalysis(tr.data,alfa,cfreq,dt,1)\n",
|
||
|
"mfa = normalizeMFT(mfa,'freq',2.0)"
|
||
|
]
|
||
|
},
|
||
|
{
|
||
|
"cell_type": "code",
|
||
|
"execution_count": null,
|
||
|
"metadata": {
|
||
|
"jupyter": {
|
||
|
"outputs_hidden": false
|
||
|
}
|
||
|
},
|
||
|
"outputs": [],
|
||
|
"source": [
|
||
|
"st.plot()\n",
|
||
|
"fg = plt.figure(figsize = [12.5,12])\n",
|
||
|
"extent = (0.,npts*dt,fmin,fmax) # extent of matrix in true time and frequency\n",
|
||
|
"asp = 0.5*npts*dt/(fmax-fmin)\n",
|
||
|
"plt.imshow(mfa,extent = extent, aspect = asp, \n",
|
||
|
" origin = \"lower\", interpolation = \"bicubic\") # plotting\n",
|
||
|
"plt.xlabel('Time [s]')\n",
|
||
|
"plt.ylabel('Frequency [Hz]')\n",
|
||
|
"plt.show()"
|
||
|
]
|
||
|
},
|
||
|
{
|
||
|
"cell_type": "code",
|
||
|
"execution_count": null,
|
||
|
"metadata": {},
|
||
|
"outputs": [],
|
||
|
"source": []
|
||
|
}
|
||
|
],
|
||
|
"metadata": {
|
||
|
"kernelspec": {
|
||
|
"display_name": "Python 3",
|
||
|
"language": "python",
|
||
|
"name": "python3"
|
||
|
},
|
||
|
"language_info": {
|
||
|
"codemirror_mode": {
|
||
|
"name": "ipython",
|
||
|
"version": 3
|
||
|
},
|
||
|
"file_extension": ".py",
|
||
|
"mimetype": "text/x-python",
|
||
|
"name": "python",
|
||
|
"nbconvert_exporter": "python",
|
||
|
"pygments_lexer": "ipython3",
|
||
|
"version": "3.6.7"
|
||
|
}
|
||
|
},
|
||
|
"nbformat": 4,
|
||
|
"nbformat_minor": 4
|
||
|
}
|