seismic.ASDFdatabase package¶
Submodules¶
seismic.ASDFdatabase.ClientUtils module¶
seismic.ASDFdatabase.FederatedASDFDataSet module¶
- Description:
Wrapper Class for providing fast access to data contained within a set of ASDF files
References:
CreationDate: 12/12/18 Developer: rakib.hassan@ga.gov.au
- Revision History:
LastUpdate: 12/12/18 RH LastUpdate: 2020-04-10 Fei Zhang clean up + added example run for the script
- class seismic.ASDFdatabase.FederatedASDFDataSet.FederatedASDFDataSet(asdf_source, logger=None, single_item_read_limit_in_mb=1024)[source]¶
Bases:
object
Initializer for FederatedASDFDataSet.
- Parameters
asdf_source – Path to a text file containing a list of ASDF files. Entries can be commented out with ‘#’
logger – logger instance
- get_closest_stations(lon, lat, nn=1)[source]¶
- Parameters
lon – longitude (degree)
lat – latitude (degrees)
nn – number of closest stations to fetch
- Returns
A tuple containing a list of closest ‘network.station’ names and a list of distances (in ascending order) in kms
- get_global_time_range(network, station, location=None, channel=None)[source]¶
- Parameters
network – network code
station – station code
location – location code (optional)
channel – channel code (optional)
- Returns
tuple containing min and max times as UTCDateTime objects. If no matching records are found min is set to 2100-01-01T00:00:00.000000Z and max is set to 1900-01-01T00:00:00.000000Z
- get_stations(starttime, endtime, network=None, station=None, location=None, channel=None)[source]¶
- Parameters
starttime – start time string in UTCDateTime format; can also be an instance of obspy.UTCDateTime
endtime – end time string in UTCDateTime format; can also be an instance of obspy.UTCDateTime
network – network code (optional)
station – station code (optional)
location – location code (optional)
channel – channel code (optional)
- Returns
a list containing [net, sta, loc, cha, lon, lat] in each row
- get_waveform_count(network, station, location, channel, starttime, endtime)[source]¶
Count the number of traces within the given parameters of network, station, etc.. and date range. This is a fast method of determing whether any trace data exists in a given time period, if you don’t actually need the waveform data itself.
- Parameters
network – network code
station – station code
location – location code
channel – channel code
starttime – start time string in UTCDateTime format; can also be an instance of obspy.UTCDateTime
endtime – end time string in UTCDateTime format; can also be an instance of obspy.UTCDateTime
- Returns
The number of streams containing waveform data over the time-range provided
- get_waveforms(network, station, location, channel, starttime, endtime, trace_count_threshold=200)[source]¶
- Parameters
network – network code
station – station code
location – location code
channel – channel code
starttime – start time string in UTCDateTime format; can also be an instance of obspy.UTCDateTime
endtime – end time string in UTCDateTime format; can also be an instance of obspy.UTCDateTime
trace_count_threshold – returns an empty Stream if the number of traces within the time-range provided exceeds the threshold (default 200). This is particularly useful for filtering out data from bad stations, e.g. those from the AU.Schools network
- Returns
an obspy.Stream containing waveform data over the time-rage provided
- local_net_sta_list()[source]¶
This function provides an iterator over the entire data volume contained in all the ASDF files listed in the text file during instantiation. When FederatedASDFDataSet is instantiated in an MPI-parallel environment, meta-data for the entire data volume are equally partitioned over all processors – in such instances, this function provides an iterator over the data allocated to a given processor. This functionality underpins parallel operations, e.g. picking arrivals.
- Returns
tuples containing [net, sta, start_time, end_time]; start- and end-times are instances of obspy.UTCDateTime
- property unique_coordinates¶
- Returns
dictionary containing [lon, lat] coordinates indexed by ‘net.sta’
seismic.ASDFdatabase.asdf2mseed module¶
- Description:
Small utility for exporting mseed files from an asdf file in parallel.
References:
CreationDate: 06/08/18 Developer: rakib.hassan@ga.gov.au
- Revision History:
LastUpdate: 04/12/17 RH LastUpdate: dd/mm/yyyy Who Optional description
- seismic.ASDFdatabase.asdf2mseed.dump_traces(ds, sn_list, start_date, end_date, length, output_folder)[source]¶
Dump mseed traces from an ASDF file in parallel
- Parameters
ds – ASDF Dataset
sn_list – station list to process
start_date – start date
end_date – end date
length – length of each mseed file
output_folder – output folder
seismic.ASDFdatabase.asdf_preprocess module¶
- Description:
Reads waveforms from an ASDF file, optionally applies instrument response correction, resamples and outputs them to another ASDF file. This preprocessing is crucial for large-scale studies involving > 10000 Green’s Functions, e.g. in ambient noise tomography. This approach significantly reduces IO bottlenecks and computational costs associated with having to apply instrument response corrections on data from a given station in alternative workflows.
References:
CreationDate: 18/07/19
Developer: rakib.hassan@ga.gov.au
- Revision History:
LastUpdate: 18/07/19 RH LastUpdate: dd/mm/yyyy Who Optional description
- seismic.ASDFdatabase.asdf_preprocess.create_station_asdf(input_asdf, output_folder, resample_rate, instrument_response_inventory, instrument_response_output, water_level)[source]¶
seismic.ASDFdatabase.create_small_chunks module¶
seismic.ASDFdatabase.plot_data_quality module¶
- Description:
Reads waveform data from a FederatedASDFDataSet and generates data-quality plots into a multi-page PDF file
For each station, the program read-in all the waveform data over the period specified (usually over a year), and then perform statistics analysis over the waveform data for subsequent plotting.
The first page of the PDF shows the stations on a basemap with marker colored according to the number of usable days (“good” data available in the day, not gap, not all zeros) And the color is determined by c=mapper.to_rgba(days): Relatively speaking, black/blue/cold-ish colors indicate higher percentage of good data in the station. And red/yellow/warm-ish colors indicate higher percentage of problematic data.
The following pages of the PDF will be plotting the daily-averaged seismic wave data with gaps and all-zeros days are shaded.
CreationDate: 19/09/19 Developer: rakib.hassan@ga.gov.au
- Revision History:
LastUpdate: 19/09/2019 RH LastUpdate: 27/05/2020 FZ Refactoring and docs run examples etc.
- Todo:
The script currently have a low pylint score 4.3/10. Need to refactor to comply with Python standard pep-8: add required docstrings. use smaller function to modularize better. Consider to generate an additional page to executive summarise info for user.
- seismic.ASDFdatabase.plot_data_quality.process_data(rank, fds, stations, start_time, end_time)[source]¶
- Parameters
rank – processor rank in parallel runs
fds – FederatedASDFDataSer instance
stations – list containing tuples (net, sta, loc, cha, lon, lat)
start_time – start-time in UTCDateTime format
end_time – end-time in UTCDateTime format
- Returns
a list containing UTCDateTimes marking the start of each day and their corresponding means in each row
seismic.ASDFdatabase.query_input_yes_no module¶
From http://code.activestate.com/recipes/577058/
- seismic.ASDFdatabase.query_input_yes_no.query_yes_no(question, default='yes')[source]¶
Ask a yes/no question via input() and return their answer.
“question” is a string that is presented to the user. “default” is the presumed answer if the user just hits <Enter>. It must be “yes” (the default), “no” or None (meaning an answer is required of the user).
The “answer” return value is one of “yes” or “no”.
seismic.ASDFdatabase.sc3toasdf module¶
- Description:
Reads waveforms (within a given time-range) from a Seiscomp3 server and dumps out ASDF files, along with a json file containing associated metadata
References:
CreationDate: 13/09/18 Developer: rakib.hassan@ga.gov.au
- Revision History:
LastUpdate: 22/08/18 RH LastUpdate: dd/mm/yyyy Who Optional description
seismic.ASDFdatabase.seisds module¶
- class seismic.ASDFdatabase.seisds.SeisDB(json_file=False, generate_numpy_index=True)[source]¶
Bases:
object
- get_unique_information()[source]¶
Method to retreive the unique channels and tags within an ASDF file from the JSON Database :return: (unique channels, unique tags)
Method to test if a channel code is related to a given net/stn/tag/loc :param chan: Channel Code, String :param net: Network Code, String :param sta: Station Code, String :param sta: Location Code, String :return: True/False, Bool