Point/Line Sets

Add Cell Connectivity to Points

class PVGeo.filters.xyz.AddCellConnToPoints(**kwargs)[source]

Bases: FilterBase

This filter will add linear cell connectivity between scattered points. You have the option to add VTK_Line or VTK_PolyLine connectivity. VTK_Line connectivity makes a straight line between the points in order (either in the order by index or using a nearest neighbor calculation). The VTK_PolyLine adds a poly line connectivity between all points as one spline (either in the order by index or using a nearest neighbor calculation). Type map is specified in vtkCellType.h.

Cell Connectivity Types:

  • 4: Poly Line

  • 3: Line

RequestData(request, inInfo, outInfo)[source]

Used by pipeline to generate output data object

_connect_cells(pdi, pdo, log_time=False)[source]

Internal helper to perform the connection

set_cell_type(cell_type)[source]

Set the cell typ by the integer id as specified in vtkCellType.h

set_use_nearest_nbr(flag)[source]

Set a flag on whether to a KDTree nearest neighbor algorithm to sort the points to before adding linear connectivity.

set_use_unique_points(flag)[source]

Set a flag on whether to only use unique points

Append Cell Centers

class PVGeo.filters.xyz.AppendCellCenters(**kwargs)[source]

Bases: FilterPreserveTypeBase

RequestData(request, inInfo, outInfo)[source]

Used by pipeline to generate output

Build Surface From Points

class PVGeo.filters.xyz.BuildSurfaceFromPoints(**kwargs)[source]

Bases: FilterBase

From the sorted x, y, and z station locations in the input PolyData, create a surface to project down from the line of those points. Use the Z cells to control the size of the mesh surface

RequestData(request, inInfo, outInfo)[source]

Execute on pipeline

static create_surface(points, z_range)[source]

From the sorted x, y, and z station locations, create a surface to display a seismic recording/migration on in space. The result is defined in the X,Y,Z-z_range 3D space.

The z_range should be treated as relative coordinates to the values given on the third column of the points array. If you want the values in the z_range to be treated as the absolute coordinates, simply do not pass any Z values in the points array - if points is N by 2, then the values in z_range will be inferred as absolute.

Parameters:
  • points (np.ndarray) – array-like of the station x and y locations

  • z_range ((npts by 2-3)) – The linear space of the z

  • location. (dimension. This will be filled out for every station) –

Returns:

pyvista.UnstructuredGrid

set_z_coords(zcoords)[source]

Set the spacings for the cells in the Z direction

Parameters:

zcoords (list or np.array(floats)) – the spacings along the Z-axis

set_z_coords_str(zcoordstr)[source]

Set the spacings for the cells in the Z direction

Parameters:

zcoordstr (str) – the spacings along the Z-axis in the UBC style

Convert XYZ Units

class PVGeo.filters.xyz.ConvertUnits(conversion='meter_to_feet', **kwargs)[source]

Bases: FilterPreserveTypeBase

Convert points in an input data object from meters to feet or vice versa. This simply uses a vtkTransformFilter and scales input data object with common conversions.

RequestData(request, inInfo, outInfo)[source]

Used by pipeline to generate output

get_conversion()[source]

Get the conversion value

static lookup_conversions(get_keys=False)[source]

All Available conversions

Returns:

dictionary of conversion units

Return type:

dict

set_conversion(key)[source]

Set the conversion via a lookup table

Extract Cell Centers

class PVGeo.filters.xyz.ExtractCellCenters(**kwargs)[source]

Bases: FilterBase

RequestData(request, inInfo, outInfo)[source]

Used by pipeline to generate output

Extract Points

class PVGeo.filters.xyz.ExtractPoints[source]

Bases: FilterBase

Extracts XYZ coordinates and point/cell data from an input vtkDataSet

RequestData(request, inInfo, outInfo)[source]

Used by pipeline to generate output

Iterate Over Points

class PVGeo.filters.xyz.IterateOverPoints(dt=1.0)[source]

Bases: FilterBase

Iterate over points in a time varying manner.

RequestData(request, inInfo, outInfo)[source]

Used by pipeline to generate output

RequestInformation(request, inInfo, outInfo)[source]

Used by pipeline to set the time information

_update_time_steps()[source]

For internal use only

get_normal()[source]

Get the current normal vector

get_point()[source]

Get the current point

get_time_step_values()[source]

Use this in ParaView decorator to register timesteps

set_decimate(percent)[source]

Set the percent (1 to 100) to decimate

set_time_delta(dt)[source]

Set the time step interval in seconds

Lat Lon To UTM

class PVGeo.filters.xyz.LonLatToUTM(**kwargs)[source]

Bases: FilterPreserveTypeBase

Converts Points from Lon Lat to UTM

RequestData(request, inInfo, outInfo)[source]

Used by pipeline to generate output

__convert_2d(lon, lat, elev)

Converts 2D Lon Lat coords to 2D XY UTM points

static get_available_ellps(idx=None)[source]

Returns the available ellps

set_ellps(ellps)[source]

Set the ellipsoid type

set_zone(zone)[source]

Set the UTM zone number

Points to Tube

class PVGeo.filters.xyz.PointsToTube(num_sides=20, radius=10.0, capping=False, **kwargs)[source]

Bases: AddCellConnToPoints

Takes points from a vtkPolyData object and constructs a line of those points then builds a polygonal tube around that line with some specified radius and number of sides.

_connect_cells(pdi, pdo, log_time=False)[source]

This uses the parent’s _connect_cells() to build a tube around

set_capping(flag)[source]

Set a boolean flag on whether or not to cap the ends of the tube

set_number_of_sides(num)[source]

Set the number of sides (resolution) for the tube

set_radius(radius)[source]

Set the radius of the tube

Rotate Points

class PVGeo.filters.xyz.RotatePoints(angle=45.0, origin=None, use_corner=True)[source]

Bases: FilterBase

Rotates XYZ coordinates in vtkPolyData around an origin at a given angle on the XY plane.

RequestData(request, inInfo, outInfo)[source]

Used by pipeline to generate output.

set_origin(xo, yo)[source]

Sets the origin to perform the rotate around.

set_rotation_degrees(theta)[source]

Sets the rotational angle in degrees.

set_use_corner(flag)[source]

A flag to use a corner of the input data set as the rotational origin.

Rotation Tool

class PVGeo.filters.xyz.RotationTool(decimals=6)[source]

Bases: object

A class that holds a set of methods/tools for performing and estimating coordinate rotations.

_converge_angle(pt1, pt2)[source]

Internal use only: pts should only be a two neighboring points.

_estimate_angle_and_spacing(pts, sample=0.5)[source]

internal use only

static _get_rotation_matrix(theta)[source]

Internal helper to generate a rotation matrix given a rotation angle

static cos_between(pts)[source]

Gets the cosine between two points

static distance_between(pts)[source]

Gets the distance between two points

estimate_and_rotate(x, y, z)[source]

A method to estimate the rotation of a set of points and correct that rotation on the XY plane

static rotate(pts, theta)[source]

Rotate points around (0,0,0) given an angle on the XY plane

static rotate_around(pts, theta, origin)[source]

Rotate points around an origins given an angle on the XY plane

static rotation_matrix(vector_orig, vector_fin)[source]

Calculate the rotation matrix required to rotate from one vector to another. For the rotation of one vector to another, there are an infinite series of rotation matrices possible. Due to axial symmetry, the rotation axis can be any vector lying in the symmetry plane between the two vectors. Hence the axis-angle convention will be used to construct the matrix with the rotation axis defined as the cross product of the two vectors. The rotation angle is the arccosine of the dot product of the two unit vectors. Given a unit vector parallel to the rotation axis, w = [x, y, z] and the rotation angle a, the rotation matrix R is:

    |  1 + (1-cos(a))*(x*x-1)   -z*sin(a)+(1-cos(a))*x*y   y*sin(a)+(1-cos(a))*x*z |
R = |  z*sin(a)+(1-cos(a))*x*y   1 + (1-cos(a))*(y*y-1)   -x*sin(a)+(1-cos(a))*y*z |
    | -y*sin(a)+(1-cos(a))*x*z   x*sin(a)+(1-cos(a))*y*z   1 + (1-cos(a))*(z*z-1)  |
Parameters:
  • vector_orig (umpy array, len 3) – The unrotated vector defined in the reference frame.

  • vector_fin (numpy array, len 3) – The rotated vector defined in the reference frame.

Note

This code was adopted from printipi under the MIT license.

static sin_between(pts)[source]

Calculate the sine angle between two points