sfepy.discrete.iga.extmods.igac module¶
- class sfepy.discrete.iga.extmods.igac.CNURBSContext¶
- R¶
- bf¶
- bfg¶
- bufBN¶
- cprint()¶
- dR_dx¶
- dR_dxi¶
- e_coors_max¶
- evaluate(coors, diff=False, **kwargs)¶
- iel¶
- sfepy.discrete.iga.extmods.igac.eval_bernstein_basis(funs, ders, x, degree)¶
- sfepy.discrete.iga.extmods.igac.eval_in_tp_coors(variable, indices, ref_coors, control_points, weights, degrees, cs, conn)¶
Evaluate a field variable (if given) or the NURBS geometry in the given tensor-product reference coordinates. The field variable is defined by its DOFs - the coefficients of the NURBS basis.
- Parameters:
- variablearray
The DOF values of the variable with n_c components, shape (:, n_c).
- indiceslist of arrays
The indices of knot spans for each axis, defining the Bezier element numbers.
- ref_coorslist of arrays
The reference coordinates in [0, 1] for each knot span for each axis, defining the reference coordinates in the Bezier elements given by indices.
- control_pointsarray
The NURBS control points.
- weightsarray
The NURBS weights.
- degreessequence of ints or int
The basis degrees in each parametric dimension.
- cslist of lists of 2D arrays
The element extraction operators in each parametric dimension.
- connarray
The connectivity of the global NURBS basis.
- Returns:
- outarray
The field variable values or NURBS geometry coordinates for the given reference coordinates.
- sfepy.discrete.iga.extmods.igac.eval_mapping_data_in_qp(qps, control_points, weights, degrees, cs, conn, cells=None)¶
Evaluate data required for the isogeometric domain reference mapping in the given quadrature points. The quadrature points are the same for all Bezier elements and should correspond to the Bernstein basis degree.
- Parameters:
- qpsarray
The quadrature points coordinates with components in [0, 1] reference element domain.
- control_pointsarray
The NURBS control points.
- weightsarray
The NURBS weights.
- degreessequence of ints or int
The basis degrees in each parametric dimension.
- cslist of lists of 2D arrays
The element extraction operators in each parametric dimension.
- connarray
The connectivity of the global NURBS basis.
- cellsarray, optional
If given, use only the given Bezier elements.
- Returns:
- bfsarray
The NURBS shape functions in the physical quadrature points of all elements.
- bfgsarray
The NURBS shape functions derivatives w.r.t. the physical coordinates in the physical quadrature points of all elements.
- detsarray
The Jacobians of the mapping to the unit reference element in the physical quadrature points of all elements.
- sfepy.discrete.iga.extmods.igac.eval_variable_in_qp(variable, qps, control_points, weights, degrees, cs, conn, cells=None)¶
Evaluate a field variable in the given quadrature points. The quadrature points are the same for all Bezier elements and should correspond to the Bernstein basis degree. The field variable is defined by its DOFs - the coefficients of the NURBS basis.
- Parameters:
- variablearray
The DOF values of the variable with n_c components, shape (:, n_c).
- qpsarray
The quadrature points coordinates with components in [0, 1] reference element domain.
- control_pointsarray
The NURBS control points.
- weightsarray
The NURBS weights.
- degreessequence of ints or int
The basis degrees in each parametric dimension.
- cslist of lists of 2D arrays
The element extraction operators in each parametric dimension.
- connarray
The connectivity of the global NURBS basis.
- cellsarray, optional
If given, use only the given Bezier elements.
- Returns:
- coorsarray
The physical coordinates of the quadrature points of all elements.
- valsarray
The field variable values in the physical quadrature points.
- detsarray
The Jacobians of the mapping to the unit reference element in the physical quadrature points.
- sfepy.discrete.iga.extmods.igac.is_nurbs(weights)¶
Return True if some weights are not one.