obspy.core.event.base.ResourceIdentifier¶
- class ResourceIdentifier(id=None, prefix='smi:local', referred_object=None)[source]¶
Bases: builtins.object
Unique identifier of any resource so it can be referred to.
In QuakeML many elements and types can have a unique id that other elements use to refer to it. This is called a ResourceIdentifier and it is used for the same purpose in the obspy.core.event classes.
In QuakeML it has to be of the following regex form:
(smi|quakeml):[\w\d][\w\d\-\.\*\(\)_~']{2,}/[\w\d\-\.\*\(\)_~'] [\w\d\-\.\*\(\)\+\?_~'=,;#/&]*
e.g.
- smi:sub.website.org/event/12345678
- quakeml:google.org/pick/unique_pick_id
smi stands for “seismological meta-information”.
In this class it can be any hashable object, e.g. most immutable objects like numbers and strings.
Parameters: - id (str, optional) A unique identifier of the element it refers to. It is not verified, that it actually is unique. The user has to take care of that. If no resource_id is given, uuid.uuid4() will be used to create one which assures uniqueness within one Python run. If no fixed id is provided, the ID will be built from prefix and a random uuid hash. The random hash part can be regenerated by the referred object automatically if it gets changed.
- prefix (str, optional) An optional identifier that will be put in front of any automatically created resource id. The prefix will only have an effect if id is not specified (for a fixed ID string). Makes automatically generated resource ids more reasonable. By default “smi:local” is used which ensures a QuakeML conform resource identifier.
- referred_object (Python object, optional) The object this instance refers to. All instances created with the same resource_id will be able to access the object as long as at least one instance actual has a reference to it.
General Usage
>>> ResourceIdentifier('2012-04-11--385392') ResourceIdentifier(id="2012-04-11--385392") >>> # If 'id' is not specified it will be generated automatically. >>> ResourceIdentifier() ResourceIdentifier(id="smi:local/...") >>> # Supplying a prefix will simply prefix the automatically generated ID >>> ResourceIdentifier(prefix='event') ResourceIdentifier(id="event/...")
ResourceIdentifiers can, and oftentimes should, carry a reference to the object they refer to. This is a weak reference which means that if the object get deleted or runs out of scope, e.g. gets garbage collected, the reference will cease to exist.
>>> from obspy.core.event import Event >>> event = Event() >>> import sys >>> ref_count = sys.getrefcount(event) >>> res_id = ResourceIdentifier(referred_object=event) >>> # The reference does not changed the reference count of the object. >>> print(ref_count == sys.getrefcount(event)) True >>> # It actually is the same object. >>> print(event is res_id.get_referred_object()) True >>> # Deleting it, or letting the garbage collector handle the object will >>> # invalidate the reference. >>> del event >>> print(res_id.get_referred_object()) None
The most powerful ability (and reason why one would want to use a resource identifier class in the first place) is that once a ResourceIdentifier with an attached referred object has been created, any other ResourceIdentifier instances with the same ID can retrieve that object. This works across all ResourceIdentifiers that have been instantiated within one Python run. This enables, e.g. the resource references between the different QuakeML elements to work in a rather natural way.
>>> event_object = Event() >>> obj_id = id(event_object) >>> res_id = "obspy.org/event/test" >>> ref_a = ResourceIdentifier(res_id) >>> # The object is refers to cannot be found yet. Because no instance that >>> # an attached object has been created so far. >>> print(ref_a.get_referred_object()) None >>> # This instance has an attached object. >>> ref_b = ResourceIdentifier(res_id, referred_object=event_object) >>> ref_c = ResourceIdentifier(res_id) >>> # All ResourceIdentifiers will refer to the same object. >>> assert(id(ref_a.get_referred_object()) == obj_id) >>> assert(id(ref_b.get_referred_object()) == obj_id) >>> assert(id(ref_c.get_referred_object()) == obj_id)
The id can be converted to a valid QuakeML ResourceIdentifier by calling the convert_id_to_quakeml_uri() method. The resulting id will be of the form:
smi:authority_id/prefix/id
>>> res_id = ResourceIdentifier(prefix='origin') >>> res_id.convert_id_to_quakeml_uri(authority_id="obspy.org") >>> res_id ResourceIdentifier(id="smi:obspy.org/origin/...") >>> res_id = ResourceIdentifier(id='foo') >>> res_id.convert_id_to_quakeml_uri() >>> res_id ResourceIdentifier(id="smi:local/foo") >>> # A good way to create a QuakeML compatibly ResourceIdentifier from >>> # scratch is >>> res_id = ResourceIdentifier(prefix='pick') >>> res_id.convert_id_to_quakeml_uri(authority_id='obspy.org') >>> res_id ResourceIdentifier(id="smi:obspy.org/pick/...") >>> # If the given ID is already a valid QuakeML >>> # ResourceIdentifier, nothing will happen. >>> res_id = ResourceIdentifier('smi:test.org/subdir/id') >>> res_id ResourceIdentifier(id="smi:test.org/subdir/id") >>> res_id.convert_id_to_quakeml_uri() >>> res_id ResourceIdentifier(id="smi:test.org/subdir/id")
ResourceIdentifiers are considered identical if the IDs are the same.
>>> # Create two different resource identifiers. >>> res_id_1 = ResourceIdentifier() >>> res_id_2 = ResourceIdentifier() >>> assert(res_id_1 != res_id_2) >>> # Equalize the IDs. NEVER do this. This just an example. >>> res_id_2.id = res_id_1.id = "smi:local/abcde" >>> assert(res_id_1 == res_id_2)
ResourceIdentifier instances can be used as dictionary keys.
>>> dictionary = {} >>> res_id = ResourceIdentifier(id="foo") >>> dictionary[res_id] = "bar1" >>> # The same ID can still be used as a key. >>> dictionary["foo"] = "bar2" >>> items = sorted(dictionary.items(), key=lambda kv: kv[1]) >>> for k, v in items: ... print(repr(k), v) ResourceIdentifier(id="foo") bar1 ...'foo' bar2
Attributes
__dict__ __doc__ str(object=’‘) -> str __module__ str(object=’‘) -> str __weakref__ list of weak references to the object (if defined) id Unique identifier of the current instance. prefix resource_id uuid Public Methods
convertIDToQuakeMLURI convert_id_to_quakeml_uri Converts the current ID to a valid QuakeML URI. copy Returns a copy of the ResourceIdentifier. getQuakeMLURI getReferredObject get_quakeml_uri Returns the ID as a valid QuakeML URI if possible. Does not get_referred_object Returns the object associated with the resource identifier. regenerate_uuid Regenerates the uuid part of the ID. Does nothing for resource setReferredObject set_referred_object Sets the object the ResourceIdentifier refers to. Private Methods
Warning
Private methods are mainly for internal/developer use and their API might change without notice.
_repr_pretty_ Special Methods
__del__ __dir__ default dir() implementation __eq__ __format__ default object formatter __hash__ Uses the same hash as the resource id. __init__ __ne__ __new__ Create and return a new object. __reduce__ helper for pickle __reduce_ex__ helper for pickle __repr__ __sizeof__ size of object in memory, in bytes __str__ __subclasshook__ Abstract classes can override this to customize issubclass().