sfepy.mesh.splinebox module

class sfepy.mesh.splinebox.SplineBox(bbox, coors, nsg=None, field=None)[source]

B-spline geometry parametrization. The geometry can be modified by moving spline control points.

static create_spb(bbox, coors, degree=3, nsg=None)[source]
evaluate(cp_values=None, outside=True)[source]

Evaluate the new position of the mesh coordinates.

Parameters:
cp_valuesarray

The actual control point values. If None, use self.control_values.

outsidebool

If True, return also the coordinates outside the spline box.

Returns:
new_coorsarray

The new position of the mesh coordinates.

evaluate_derivative(cpoint, dirvec)[source]

Evaluate derivative of the spline in a given control point and direction.

Parameters:
cpointint, list

The position (index or grid indicies) of the spline control point.

dirvecarray

The directional vector.

Returns:
diffarray

The derivative field.

static gen_cp_idxs(ncp)[source]
get_box_matrix()[source]
Returns:
mtx2D array

The matrix containing the coefficients of b-spline basis functions.

get_control_points(init=False)[source]

Get the spline control points coordinates.

Returns:
cpt_coorsarray

The coordinates of the spline control points.

initbool

If True, return the initial state.

get_coors_shape()[source]

Get the shape of the coordinates.

move_control_point(cpoint, val)[source]

Change shape of spline parametrization.

Parameters:
cpointint, list

The position (index or grid indicies) of the spline control point.

valarray

Displacement.

set_control_points(cpt_coors, add=False)[source]

Set the spline control points position.

Parameters:
cpt_coorsarray

The coordinates of the spline control points.

addbool

If True, coors += cpt_coors

write_control_net(filename, deform_by_values=True)[source]

Write the SplineBox shape to the VTK file.

Parameters:
filenamestr

The VTK file name.

class sfepy.mesh.splinebox.SplineRegion2D(spl_bnd, coors, rho=1000.0)[source]

B-spline geometry parametrization. The boundary of the SplineRegion2D is defined by BSpline curves.

static create_spb(spl_bnd, coors, rho=10)[source]

Initialize SplineBox knots, control points, base functions, …

static define_control_points(cp_bnd_coors, ncp)[source]

Find positions of “inner” control points depending on boundary splines.

find_ts(coors)[source]

Function finds parameters (t, s) corresponding to given points (coors).

static points_in_poly(points, poly, tol=1e-06)[source]

Find which points are located inside the polygon.