def cylindrical_to_cartesian(grid, vec=None):
    """
    Take a grid defined in cylindrical coordinates :math:`(r, \theta, z)` and
    transform it to cartesian coordinates.
    """
    grid = np.atleast_2d(grid)

    if vec is None:
        return np.hstack([
            mkvc(grid[:, 0] * np.cos(grid[:, 1]), 2),
            mkvc(grid[:, 0] * np.sin(grid[:, 1]), 2),
            mkvc(grid[:, 2], 2),
        ])

    if len(vec.shape) == 1 or vec.shape[1] == 1:
        vec = vec.reshape(grid.shape, order="F")

    x = vec[:, 0] * np.cos(grid[:, 1]) - vec[:, 1] * np.sin(grid[:, 1])
    y = vec[:, 0] * np.sin(grid[:, 1]) + vec[:, 1] * np.cos(grid[:, 1])

    newvec = [x, y]
    if grid.shape[1] == 3:
        z = vec[:, 2]
        newvec += [z]

    return np.vstack(newvec).T
def cartesian_to_cylindrical(grid, vec=None):
    """
    Take a grid defined in cartesian coordinates and transform it to cyl
    coordinates
    """
    if vec is None:
        vec = grid

    vec = np.atleast_2d(vec)
    grid = np.atleast_2d(grid)

    theta = np.arctan2(grid[:, 1], grid[:, 0])

    return np.hstack([
        mkvc(np.cos(theta) * vec[:, 0] + np.sin(theta) * vec[:, 1], 2),
        mkvc(-np.sin(theta) * vec[:, 0] + np.cos(theta) * vec[:, 1], 2),
        mkvc(vec[:, 2], 2),
    ])
def interpolation_matrix(locs, x, y=None, z=None):
    """Local interpolation computed for each receiver point in turn

    :param numpy.ndarray loc: Location of points to interpolate to
    :param numpy.ndarray x: Tensor of 1st dimension of grid.
    :param numpy.ndarray y: Tensor of 2nd dimension of grid. None by default.
    :param numpy.ndarray z: Tensor of 3rd dimension of grid. None by default.
    :rtype: scipy.sparse.csr_matrix
    :return: Interpolation matrix

    .. plot::

        import discretize
        import numpy as np
        import matplotlib.pyplot as plt
        locs = np.random.rand(50)*0.8+0.1
        x = np.linspace(0, 1, 7)
        dense = np.linspace(0, 1, 200)
        fun = lambda x: np.cos(2*np.pi*x)
        Q = discretize.utils.interpolation_matrix(locs, x)
        plt.plot(x, fun(x), 'bs-')
        plt.plot(dense, fun(dense), 'y:')
        plt.plot(locs, Q*fun(x), 'mo')
        plt.plot(locs, fun(locs), 'rx')
        plt.show()

    """

    npts = locs.shape[0]
    locs = locs.astype(float)
    x = x.astype(float)
    if y is None and z is None:
        shape = [x.size]
        inds, vals = _interpmat1D(mkvc(locs), x)
    elif z is None:
        y = y.astype(float)
        shape = [x.size, y.size]
        inds, vals = _interpmat2D(locs, x, y)
    else:
        y = y.astype(float)
        z = z.astype(float)
        shape = [x.size, y.size, z.size]
        inds, vals = _interpmat3D(locs, x, y, z)

    I = np.repeat(range(npts), 2 ** len(shape))
    J = sub2ind(shape, inds)
    Q = sp.csr_matrix((vals, (I, J)), shape=(npts, np.prod(shape)))
    return Q
예제 #4
0
def rotate_points_from_normals(xyz, v0, v1, x0=np.r_[0.0, 0.0, 0.0]):
    """Rotate a set of xyz locations about a specified point.

    Rotate a grid of Cartesian points about a location x0 according to the
    rotation defined from vector v0 to v1.

    Let :math:`\\mathbf{x}` represent an input xyz location, let :math:`\\mathbf{x_0}` be
    the origin of rotation, and let :math:`\\mathbf{R}` denote the rotation matrix from
    vector v0 to v1. Where :math:`\\mathbf{x'}` is the new xyz location, this function
    outputs the following operation for all input locations:

    .. math::
        \\mathbf{x'} = \\mathbf{R (x - x_0)} + \\mathbf{x_0}

    Parameters
    ----------
    xyz : (n, 3) numpy.ndarray
        locations to rotate
    v0 : (3) numpy.ndarray
        Starting orientation direction
    v1 : (3) numpy.ndarray
        Finishing orientation direction
    x0 : (3) numpy.ndarray, optional
        The origin of rotation.

    Returns
    -------
    (n, 3) numpy.ndarray
        The rotated xyz locations.
    """

    # Compute rotation matrix between v0 and v1
    R = rotation_matrix_from_normals(v0, v1)

    if xyz.shape[1] != 3:
        raise ValueError("Grid of xyz points should be n x 3")
    if len(x0) != 3:
        raise ValueError("x0 should have length 3")

    # Define origin
    X0 = np.ones([xyz.shape[0], 1]) * mkvc(x0)

    return (xyz - X0).dot(R.T) + X0  # equivalent to (R*(xyz - X0)).T + X0
def rotate_points_from_normals(XYZ, n0, n1, x0=np.r_[0.0, 0.0, 0.0]):
    """
    rotates a grid so that the vector n0 is aligned with the vector n1

    :param numpy.array n0: vector of length 3, should have norm 1
    :param numpy.array n1: vector of length 3, should have norm 1
    :param numpy.array x0: vector of length 3, point about which we perform the rotation
    :rtype: numpy.array, 3x3
    :return: rotation matrix which rotates the frame so that n0 is aligned with n1
    """

    R = rotation_matrix_from_normals(n0, n1)

    if XYZ.shape[1] != 3:
        raise ValueError("Grid XYZ should be 3 wide")
    if len(x0) != 3:
        raise ValueError("x0 should have length 3")

    X0 = np.ones([XYZ.shape[0], 1]) * mkvc(x0)

    return (XYZ - X0).dot(R.T) + X0  # equivalent to (R*(XYZ - X0)).T + X0
예제 #6
0
def face_info(xyz, A, B, C, D, average=True, normalize_normals=True, **kwargs):
    """Returns normal surface vector and area for a given set of faces.

    Let *xyz* be an (n, 3) array denoting a set of vertex locations.
    Now let vertex locations *a, b, c* and *d* define a quadrilateral
    (regular or irregular) in 2D or 3D space. For this quadrilateral,
    we organize the vertices as follows:

    CELL VERTICES::

            a -------Vab------- b
           /                   /
          /                   /
        Vda       (X)       Vbc
        /                   /
       /                   /
      d -------Vcd------- c

    where the normal vector *(X)* is pointing into the page. For a set
    of quadrilaterals whose vertices are indexed in arrays *A, B, C* and *D* ,
    this function returns the normal surface vector(s) and the area
    for each quadrilateral.

    At each vertex, there are 4 cross-products that can be used to compute the
    vector normal the surface defined by the quadrilateral. In 3D space however,
    the vertices indexed may not define a quadrilateral exactly and thus the normal vectors
    computed at each vertex might not be identical. In this case, you may choose output
    the normal vector at *a, b, c* and *d* or compute
    the average normal surface vector as follows:

    .. math::
        \\bf{n} = \\frac{1}{4} \\big (
        \\bf{v_{ab} \\times v_{da}} +
        \\bf{v_{bc} \\times v_{ab}} +
        \\bf{v_{cd} \\times v_{bc}} +
        \\bf{v_{da} \\times v_{cd}} \\big )


    For computing the surface area, we assume the vertices define a quadrilateral.

    Parameters
    ----------
    xyz : (n, 3) numpy.ndarray
        The x,y, and z locations for all verticies
    A : (n_face) numpy.ndarray
        Vector containing the indicies for the **a** vertex locations
    B : (n_face) numpy.ndarray
        Vector containing the indicies for the **b** vertex locations
    C : (n_face) numpy.ndarray
        Vector containing the indicies for the **c** vertex locations
    D : (n_face) numpy.ndarray
        Vector containing the indicies for the **d** vertex locations
    average : bool, optional
        If *True*, the function returns the average surface
        normal vector for each surface. If *False* , the function will
        return the normal vectors computed at the *A, B, C* and *D*
        vertices in a cell array {nA,nB,nC,nD}.
    normalize_normal : bool, optional
        If *True*, the function will normalize the surface normal
        vectors. This is applied regardless of whether the *average* parameter
        is set to *True* or *False*. If *False*, the vectors are not normalized.

    Returns
    -------
    N : (n_face) numpy.ndarray or (4) list of (n_face) numpy.ndarray
        Normal vector(s) for each surface. If *average* = *True*, the function
        returns an ndarray with the average surface normal vectos. If *average* = *False* ,
        the function returns a cell array {nA,nB,nC,nD} containing the normal vectors
        computed using each vertex of the surface.
    area : (n_face) numpy.ndarray
        The surface areas.

    Examples
    --------
    Here we define a set of vertices for a tensor mesh. We then
    index 4 vertices for an irregular quadrilateral. The
    function *face_info* is used to compute the normal vector
    and the surface area.

    >>> from discretize.utils import face_info
    >>> from discretize import TensorMesh
    >>> import numpy as np
    >>> import matplotlib.pyplot as plt
    >>> import matplotlib as mpl
    >>> mpl.rcParams.update({"font.size": 14})

    Define Corners of a uniform cube.

    >>> h = [1, 1]
    >>> mesh = TensorMesh([h, h, h])
    >>> xyz = mesh.nodes

    Choose the face indices,

    >>> A = np.array([0])
    >>> B = np.array([4])
    >>> C = np.array([26])
    >>> D = np.array([18])

    Compute average surface normal vector (normalized),

    >>> nvec, area = face_info(xyz, A, B, C, D)
    >>> nvec, area
    (array([[-0.70710678,  0.70710678,  0.        ]]), array([4.24264069]))

    Plot surface for example 1 on mesh

    .. collapse:: Expand to show scripting for plot

        >>> fig = plt.figure(figsize=(7, 7))
        >>> ax = fig.gca(projection='3d')
        >>> mesh.plot_grid(ax=ax)
        >>> k = [0, 4, 26, 18, 0]
        >>> xyz_quad = xyz[k, :]
        >>> ax.plot(xyz_quad[:, 0], xyz_quad[:, 1], xyz_quad[:, 2], 'r')
        >>> ax.text(-0.25, 0., 3., 'Area of the surface: {:g} $m^2$'.format(area[0]))
        >>> ax.text(-0.25, 0., 2.8, 'Normal vector: ({:.2f}, {:.2f}, {:.2f})'.format(
        ...     nvec[0, 0], nvec[0, 1], nvec[0, 2])
        ... )
        >>> plt.show()

    In our second example, the vertices are unable to define a flat
    surface in 3D space. However, we will demonstrate the *face_info*
    returns the average normal vector and an approximate surface area.

    Define the face indicies
    >>> A = np.array([0])
    >>> B = np.array([5])
    >>> C = np.array([26])
    >>> D = np.array([18])

    Compute average surface normal vector

    >>> nvec, area = face_info(xyz, A, B, C, D)
    >>> nvec, area
    (array([[-0.4472136 ,  0.89442719,  0.        ]]), array([2.23606798]))

    Plot surface for example 2 on mesh

    .. collapse:: Expand to show scripting for plot

        >>> fig = plt.figure(figsize=(7, 7))
        >>> ax = fig.gca(projection='3d')
        >>> mesh.plot_grid(ax=ax)
        >>> k = [0, 5, 26, 18, 0]
        >>> xyz_quad = xyz[k, :]
        >>> ax.plot(xyz_quad[:, 0], xyz_quad[:, 1], xyz_quad[:, 2], 'g')
        >>> ax.text(-0.25, 0., 3., 'Area of the surface: {:g} $m^2$'.format(area[0]))
        >>> ax.text(-0.25, 0., 2.8, 'Average normal vector: ({:.2f}, {:.2f}, {:.2f})'.format(
        ...     nvec[0, 0], nvec[0, 1], nvec[0, 2])
        ... )
        >>> plt.show()
    """
    if "normalizeNormals" in kwargs:
        warnings.warn(
            "The normalizeNormals keyword argument has been deprecated, please use normalize_normals. "
            "This will be removed in discretize 1.0.0",
            FutureWarning,
        )
        normalize_normals = kwargs["normalizeNormals"]
    if not isinstance(average, bool):
        raise TypeError("average must be a boolean")
    if not isinstance(normalize_normals, bool):
        raise TypeError("normalize_normals must be a boolean")

    AB = xyz[B, :] - xyz[A, :]
    BC = xyz[C, :] - xyz[B, :]
    CD = xyz[D, :] - xyz[C, :]
    DA = xyz[A, :] - xyz[D, :]

    def cross(X, Y):
        return np.c_[X[:, 1] * Y[:, 2] - X[:, 2] * Y[:, 1],
                     X[:, 2] * Y[:, 0] - X[:, 0] * Y[:, 2],
                     X[:, 0] * Y[:, 1] - X[:, 1] * Y[:, 0], ]

    nA = cross(AB, DA)
    nB = cross(BC, AB)
    nC = cross(CD, BC)
    nD = cross(DA, CD)

    length = lambda x: np.sqrt(x[:, 0]**2 + x[:, 1]**2 + x[:, 2]**2)
    normalize = lambda x: x / np.kron(np.ones(
        (1, x.shape[1])), mkvc(length(x), 2))
    if average:
        # average the normals at each vertex.
        N = (nA + nB + nC + nD) / 4  # this is intrinsically weighted by area
        # normalize
        N = normalize(N)
    else:
        if normalize_normals:
            N = [normalize(nA), normalize(nB), normalize(nC), normalize(nD)]
        else:
            N = [nA, nB, nC, nD]

    # Area calculation
    #
    # Approximate by 4 different triangles, and divide by 2.
    # Each triangle is one half of the length of the cross product
    #
    # So also could be viewed as the average parallelogram.
    #
    # TODO: This does not compute correctly for concave quadrilaterals
    area = (length(nA) + length(nB) + length(nC) + length(nD)) / 4

    return N, area
예제 #7
0
def cartesian_to_cylindrical(grid, vec=None):
    """
    Transform a grid or a vector from Cartesian coordinates :math:`(x, y, z)` to
    cylindrical coordinates :math:`(r, \\theta, z)`.

    Parameters
    ----------
    grid : (n, 3) array_like
        Location points defined in Cartesian coordinates :math:`(x, y z)`.
    vec : (n, 3) array_like, optional
        Vector defined in Cartesian coordinates. This also accepts a flattened array
        with the same total elements in column major order.

    Returns
    -------
    (n, 3) numpy.ndarray
        If `vec` is ``None``, this returns the transformed `grid` array, otherwise
        this is the transformed `vec` array.

    Examples
    --------
    Here, we convert a series of vectors in 3D space from Cartesian coordinates
    to cylindrical coordinates.

    >>> from discretize.utils import cartesian_to_cylindrical
    >>> import numpy as np

    Create set of vectors in Cartesian coordinates

    >>> r = np.ones(9)
    >>> phi = np.linspace(0, 2*np.pi, 9)
    >>> z = np.linspace(-4., 4., 9)
    >>> x = r*np.cos(phi)
    >>> y = r*np.sin(phi)
    >>> u = np.c_[x, y, z]
    >>> u
    array([[ 1.00000000e+00,  0.00000000e+00, -4.00000000e+00],
           [ 7.07106781e-01,  7.07106781e-01, -3.00000000e+00],
           [ 6.12323400e-17,  1.00000000e+00, -2.00000000e+00],
           [-7.07106781e-01,  7.07106781e-01, -1.00000000e+00],
           [-1.00000000e+00,  1.22464680e-16,  0.00000000e+00],
           [-7.07106781e-01, -7.07106781e-01,  1.00000000e+00],
           [-1.83697020e-16, -1.00000000e+00,  2.00000000e+00],
           [ 7.07106781e-01, -7.07106781e-01,  3.00000000e+00],
           [ 1.00000000e+00, -2.44929360e-16,  4.00000000e+00]])

    Compute equivalent set of vectors in cylindrical coordinates

    >>> v = cartesian_to_cylindrical(u)
    >>> v
    array([[ 1.00000000e+00,  0.00000000e+00, -4.00000000e+00],
           [ 1.00000000e+00,  7.85398163e-01, -3.00000000e+00],
           [ 1.00000000e+00,  1.57079633e+00, -2.00000000e+00],
           [ 1.00000000e+00,  2.35619449e+00, -1.00000000e+00],
           [ 1.00000000e+00,  3.14159265e+00,  0.00000000e+00],
           [ 1.00000000e+00, -2.35619449e+00,  1.00000000e+00],
           [ 1.00000000e+00, -1.57079633e+00,  2.00000000e+00],
           [ 1.00000000e+00, -7.85398163e-01,  3.00000000e+00],
           [ 1.00000000e+00, -2.44929360e-16,  4.00000000e+00]])
    """
    grid = as_array_n_by_dim(grid, 3)
    theta = np.arctan2(grid[:, 1], grid[:, 0])
    if vec is None:
        return np.c_[np.linalg.norm(grid[:, :2], axis=-1), theta, grid[:, 2]]
    vec = as_array_n_by_dim(vec, 3)

    return np.hstack([
        mkvc(np.cos(theta) * vec[:, 0] + np.sin(theta) * vec[:, 1], 2),
        mkvc(-np.sin(theta) * vec[:, 0] + np.cos(theta) * vec[:, 1], 2),
        mkvc(vec[:, 2], 2),
    ])
예제 #8
0
def cylindrical_to_cartesian(grid, vec=None):
    """
    Transform a grid or a vector from cylindrical coordinates :math:`(r, \\theta, z)` to
    Cartesian coordinates :math:`(x, y, z)`. :math:`\\theta` is given in radians.

    Parameters
    ----------
    grid : (n, 3) array_like
        Location points defined in cylindrical coordinates :math:`(r, \\theta, z)`.
    vec : (n, 3) array_like, optional
        Vector defined in cylindrical coordinates :math:`(r, \\theta, z)` at the
        locations grid. Will also except a flattend array in column major order with the
        same number of elements.

    Returns
    -------
    (n, 3) numpy.ndarray
        If `vec` is ``None``, this returns the transformed `grid` array, otherwise
        this is the transformed `vec` array.

    Examples
    --------
    Here, we convert a series of vectors in 3D space from cylindrical coordinates
    to Cartesian coordinates.

    >>> from discretize.utils import cylindrical_to_cartesian
    >>> import numpy as np

    Construct original set of vectors in cylindrical coordinates

    >>> r = np.ones(9)
    >>> phi = np.linspace(0, 2*np.pi, 9)
    >>> z = np.linspace(-4., 4., 9)
    >>> u = np.c_[r, phi, z]
    >>> u
    array([[ 1.        ,  0.        , -4.        ],
           [ 1.        ,  0.78539816, -3.        ],
           [ 1.        ,  1.57079633, -2.        ],
           [ 1.        ,  2.35619449, -1.        ],
           [ 1.        ,  3.14159265,  0.        ],
           [ 1.        ,  3.92699082,  1.        ],
           [ 1.        ,  4.71238898,  2.        ],
           [ 1.        ,  5.49778714,  3.        ],
           [ 1.        ,  6.28318531,  4.        ]])

    Create equivalent set of vectors in Cartesian coordinates

    >>> v = cylindrical_to_cartesian(u)
    >>> v
    array([[ 1.00000000e+00,  0.00000000e+00, -4.00000000e+00],
           [ 7.07106781e-01,  7.07106781e-01, -3.00000000e+00],
           [ 6.12323400e-17,  1.00000000e+00, -2.00000000e+00],
           [-7.07106781e-01,  7.07106781e-01, -1.00000000e+00],
           [-1.00000000e+00,  1.22464680e-16,  0.00000000e+00],
           [-7.07106781e-01, -7.07106781e-01,  1.00000000e+00],
           [-1.83697020e-16, -1.00000000e+00,  2.00000000e+00],
           [ 7.07106781e-01, -7.07106781e-01,  3.00000000e+00],
           [ 1.00000000e+00, -2.44929360e-16,  4.00000000e+00]])
    """
    grid = np.atleast_2d(grid)

    if vec is None:
        return np.hstack([
            mkvc(grid[:, 0] * np.cos(grid[:, 1]), 2),
            mkvc(grid[:, 0] * np.sin(grid[:, 1]), 2),
            mkvc(grid[:, 2], 2),
        ])
    vec = np.asanyarray(vec)
    if len(vec.shape) == 1 or vec.shape[1] == 1:
        vec = vec.reshape(grid.shape, order="F")

    x = vec[:, 0] * np.cos(grid[:, 1]) - vec[:, 1] * np.sin(grid[:, 1])
    y = vec[:, 0] * np.sin(grid[:, 1]) + vec[:, 1] * np.cos(grid[:, 1])

    newvec = [x, y]
    if grid.shape[1] == 3:
        z = vec[:, 2]
        newvec += [z]

    return np.vstack(newvec).T
예제 #9
0
def interpolation_matrix(locs, x, y=None, z=None):
    """
    Generate interpolation matrix which maps a tensor quantity to a set of locations.

    This function generates a sparse matrix for interpolating tensor quantities to a set
    of specified locations. It uses nD linear interpolation. The user may generate the
    interpolation matrix for tensor quantities that live on 1D, 2D or 3D tensors. This
    functionality is frequently used to interpolate quantites from cell centers or nodes
    to specified locations.

    In higher dimensions the ordering of the output has the 1st dimension changing the
    quickest.

    Parameters
    ----------
    locs : (n, dim) numpy.ndarray
        The locations for the interpolated values. Here *n* is
        the number of locations and *dim* is the dimension (1, 2 or 3)
    x : (nx) numpy.ndarray
        Vector defining the locations of the tensor along the x-axis
    y : (ny) numpy.ndarray, optional
        Vector defining the locations of the tensor along the y-axis. Required if
        ``dim`` is 2.
    z : (nz) numpy.ndarray, optional
        Vector defining the locations of the tensor along the z-axis. Required if
        ``dim`` is 3.

    Returns
    -------
    (n, nx * ny * nz) scipy.sparse.csr_matrix
        A sparse matrix which interpolates the tensor quantity on cell centers or nodes
        to the set of specified locations.

    Examples
    --------
    Here is a 1D example where a function evaluated on a regularly spaced grid
    is interpolated to a set of random locations. To compare the accuracy, the
    function is evaluated at the set of random locations.

    >>> from discretize.utils import interpolation_matrix
    >>> from discretize import TensorMesh
    >>> import numpy as np
    >>> import matplotlib.pyplot as plt
    >>> np.random.seed(14)

    Create an interpolation matrix

    >>> locs = np.random.rand(50)*0.8+0.1
    >>> x = np.linspace(0, 1, 7)
    >>> dense = np.linspace(0, 1, 200)
    >>> fun = lambda x: np.cos(2*np.pi*x)
    >>> Q = interpolation_matrix(locs, x)

    Plot original function and interpolation

    .. collapse:: Expand to show scripting for plot

        >>> fig1 = plt.figure(figsize=(5, 3))
        >>> ax = fig1.add_axes([0.1, 0.1, 0.8, 0.8])
        >>> ax.plot(dense, fun(dense), 'k:', lw=3)
        >>> ax.plot(x, fun(x), 'ks', markersize=8)
        >>> ax.plot(locs, Q*fun(x), 'go', markersize=4)
        >>> ax.plot(locs, fun(locs), 'rs', markersize=4)
        >>> ax.legend(
        ...     [
        ...         'True Function',
        ...         'True (discrete loc.)',
        ...         'Interpolated (computed)',
        ...         'True (interp. loc.)'
        ...     ],
        ...     loc='upper center'
        ... )
        >>> plt.show()

    Here, demonstrate a similar example on a 2D mesh using a 2D Gaussian distribution.
    We interpolate the Gaussian from the nodes to cell centers and examine the relative
    error.

    >>> hx = np.ones(10)
    >>> hy = np.ones(10)
    >>> mesh = TensorMesh([hx, hy], x0='CC')
    >>> def fun(x, y):
    ...     return np.exp(-(x**2 + y**2)/2**2)

    Define the the value at the mesh nodes,

    >>> nodes = mesh.nodes
    >>> val_nodes = fun(nodes[:, 0], nodes[:, 1])

    >>> centers = mesh.cell_centers
    >>> A = interpolation_matrix(
    ...     centers, mesh.nodes_x, mesh.nodes_y
    ... )
    >>> val_interp = A.dot(val_nodes)

    Plot the interpolated values, along with the true values at cell centers,

    .. collapse:: Expand to show scripting for plot

        >>> val_centers = fun(centers[:, 0], centers[:, 1])
        >>> fig = plt.figure(figsize=(11,3.3))
        >>> clim = (0., 1.)
        >>> ax1 = fig.add_subplot(131)
        >>> ax2 = fig.add_subplot(132)
        >>> ax3 = fig.add_subplot(133)
        >>> mesh.plot_image(val_centers, ax=ax1, clim=clim)
        >>> mesh.plot_image(val_interp, ax=ax2, clim=clim)
        >>> mesh.plot_image(val_centers-val_interp, ax=ax3, clim=clim)
        >>> ax1.set_title('Analytic at Centers')
        >>> ax2.set_title('Interpolated from Nodes')
        >>> ax3.set_title('Relative Error')
        >>> plt.show()
    """

    npts = locs.shape[0]
    locs = locs.astype(float)
    x = x.astype(float)
    if y is None and z is None:
        shape = [x.size]
        inds, vals = _interpmat1D(mkvc(locs), x)
    elif z is None:
        y = y.astype(float)
        shape = [x.size, y.size]
        inds, vals = _interpmat2D(locs, x, y)
    else:
        y = y.astype(float)
        z = z.astype(float)
        shape = [x.size, y.size, z.size]
        inds, vals = _interpmat3D(locs, x, y, z)

    I = np.repeat(range(npts), 2 ** len(shape))
    J = sub2ind(shape, inds)
    Q = sp.csr_matrix((vals, (I, J)), shape=(npts, np.prod(shape)))
    return Q
예제 #10
0
def face_info(xyz, A, B, C, D, average=True, normalize_normals=True, **kwargs):
    """
    function [N] = face_info(y,A,B,C,D)

       Returns the averaged normal, area, and edge lengths for a given set of faces.

       If average option is FALSE then N is a cell array {nA,nB,nC,nD}


    Input:
       xyz          - X,Y,Z vertex vector
       A,B,C,D      - vert index of the face (counter clockwize)

    Options:
       average      - [true]/false, toggles returning all normals or the average

    Output:
       N            - average face normal or {nA,nB,nC,nD} if average = false
       area         - average face area
       edgeLengths  - exact edge Lengths, 4 column vector [AB, BC, CD, DA]

    see also testFaceNormal testFaceArea

    @author Rowan Cockett

    Last modified on: 2013/07/26

    """
    if "normalizeNormals" in kwargs:
        warnings.warn(
            "The normalizeNormals keyword argument has been deprecated, please use normalize_normals. "
            "This will be removed in discretize 1.0.0",
            DeprecationWarning,
        )
        normalize_normals = kwargs["normalizeNormals"]
    if not isinstance(average, bool):
        raise TypeError("average must be a boolean")
    if not isinstance(normalize_normals, bool):
        raise TypeError("normalize_normals must be a boolean")
    # compute normal that is pointing away from you.
    #
    #    A -------A-B------- B
    #    |                   |
    #    |                   |
    #   D-A       (X)       B-C
    #    |                   |
    #    |                   |
    #    D -------C-D------- C

    AB = xyz[B, :] - xyz[A, :]
    BC = xyz[C, :] - xyz[B, :]
    CD = xyz[D, :] - xyz[C, :]
    DA = xyz[A, :] - xyz[D, :]

    def cross(X, Y):
        return np.c_[X[:, 1] * Y[:, 2] - X[:, 2] * Y[:, 1],
                     X[:, 2] * Y[:, 0] - X[:, 0] * Y[:, 2],
                     X[:, 0] * Y[:, 1] - X[:, 1] * Y[:, 0], ]

    nA = cross(AB, DA)
    nB = cross(BC, AB)
    nC = cross(CD, BC)
    nD = cross(DA, CD)

    length = lambda x: np.sqrt(x[:, 0]**2 + x[:, 1]**2 + x[:, 2]**2)
    normalize = lambda x: x / np.kron(np.ones(
        (1, x.shape[1])), mkvc(length(x), 2))
    if average:
        # average the normals at each vertex.
        N = (nA + nB + nC + nD) / 4  # this is intrinsically weighted by area
        # normalize
        N = normalize(N)
    else:
        if normalize_normals:
            N = [normalize(nA), normalize(nB), normalize(nC), normalize(nD)]
        else:
            N = [nA, nB, nC, nD]

    # Area calculation
    #
    # Approximate by 4 different triangles, and divide by 2.
    # Each triangle is one half of the length of the cross product
    #
    # So also could be viewed as the average parallelogram.
    #
    # TODO: This does not compute correctly for concave quadrilaterals
    area = (length(nA) + length(nB) + length(nC) + length(nD)) / 4

    return N, area