obspy.signal.rotate.rotate2zne¶
- rotate2zne(data_1, azimuth_1, dip_1, data_2, azimuth_2, dip_2, data_3, azimuth_3, dip_3, inverse=False)[source]¶
Rotates an arbitrarily oriented three-component vector to ZNE.
Each components orientation is described with a azimuth and a dip. The azimuth is defined as the degrees from North, clockwise and the dip is the defined as the number of degrees, down from horizontal. Both definitions are according to the SEED standard.
The three components need not be orthogonal to each other but the components have to be linearly independent. The function performs a full base change to orthogonal Vertical, North, and East orientations.
Parameters: - data_1 Data component 1.
- azimuth_1 The azimuth of component 1.
- dip_1 The dip of component 1.
- data_2 Data component 2.
- azimuth_2 The azimuth of component 2.
- dip_2 The dip of component 2.
- data_3 Data component 3.
- azimuth_3 The azimuth of component 3.
- dip_3 The dip of component 3.
- inverse (bool) If True, the data arrays will be converted from ZNE to whatever coordinate system the azimuths and dips specify. In that case data_1, data_2, data_3 have to be data arrays for Z, N, and E and the dips and azimuths specify where to transform to.
Return type: Tuple of three NumPy arrays.
Returns: The three rotated components, oriented in Z, N, and E if inverse is False. Otherwise they will be oriented as specified by the dips and azimuths.
An input of ZNE yields an output of ZNE
>>> rotate2zne(np.arange(3), 0, -90, np.arange(3) * 2, 0, 0, np.arange(3) * 3, 90, 0) (array([ 0., 1., 2.]), array([ 0., 2., 4.]), array([ 0., 3., 6.]))
An input of ZSE yields an output of ZNE
>>> rotate2zne(np.arange(3), 0, -90, np.arange(3) * 2, 180, 0, np.arange(3) * 3, 90, 0) (array([ 0., 1., 2.]), array([ 0., -2., -4.]), array([ 0., 3., 6.]))
Mixed up components should get rotated to ZNE.
>>> rotate2zne(np.arange(3), 0, 0, np.arange(3) * 2, 90, 0, np.arange(3) * 3, 0, -90) (array([ 0., 3., 6.]), array([ 0., 1., 2.]), array([ 0., 2., 4.]))