def _compute_jacobian_numerically(form, force): r"""Compute the Jacobian matrix. The actual computation of the Jacobían matrix d_X/dX where X contains the form diagram coordinates in *Fortran* order (first all x-coordinates, then all y-coordinates) and _X contains the force diagram coordinates in *Fortran* order (first all _x-coordinates, then all _y-coordinates) The Jacobian is computed numerically using forward differences. Parameters ---------- form : compas_ags.diagrams.formdiagram.FormDiagram The form diagram. force : compas_bi_ags.diagrams.forcediagram.ForceDiagram The force diagram. Returns ------- jacobian Jacobian matrix ( 2 * _vcount, 2 * vcount) """ # Update force diagram form_update_q_from_qind(form) force_update_from_form(force, form) # Get current coordinates xy = array(form.xy(), dtype=float64).reshape((-1, 2)) X = np.vstack( (np.asmatrix(xy[:, 0]).transpose(), np.asmatrix(xy[:, 1]).transpose())) _X = _comp_perturbed_force_coordinates_from_form(form, force, X) vcount = form.number_of_vertices() _vcount = force.number_of_vertices() # -------------------------------------------------------------------------- # Jacobian from forward differencing # -------------------------------------------------------------------------- jacobian = np.zeros((_vcount * 2, vcount * 2)) nv = len(X) for i in range(nv): du = np.zeros(X.shape) du[i] = 1 _Xi_pert = _comp_perturbed_force_coordinates_from_form( form, force, X + du * 1e-10) dX = np.divide(_Xi_pert - _X, 1e-10) jacobian[:, i] = dX.A1 return jacobian
def _comp_perturbed_force_coordinates_from_form(form, force, X): r"""Compute force diagram coordinates based on X. Used to compute perturbed force diagram coordinates by sending in perturbed form diagram coordinates. Useful for computing the Jacobian matrix numerically. Parameters ---------- form : compas_ags.formdiagram.FormDiagram The form diagram. force : compas_ags.forcediagram.ForceDiagram The force diagram. X Perturbed form diagram coordinates in *Fortran* order (first all x-coordinates, then all y-coordinates). Returns ------- _X Perturbed force diagram coordinates in *Fortran* order (first all _x-coordinates, then all _y-coordinates). """ # Save current coordinates xyo = array(form.xy(), dtype=float64).reshape((-1, 2)) _xyo = array(force.xy(), dtype=float64) # Perturb form and update force nx = xyo.shape[0] ny = xyo.shape[0] xy = array(form.xy(), dtype=float64).reshape((-1, 2)) xy[:, 0] = X[:nx].T xy[:, 1] = X[nx:(nx + ny)].T _update_coordinates(form, xy) form_update_q_from_qind(form) force_update_from_form(force, form) _xy = array(force.xy(), dtype=float64) _X = np.vstack( (np.asmatrix(_xy[:, 0]).transpose(), np.asmatrix(_xy[:, 1]).transpose())) # Revert to current coordinates _update_coordinates(form, xyo) _update_coordinates(force, _xyo) return _X
(5, 18), ] form = FormDiagram.from_vertices_and_edges(vertices, edges) force = ForceDiagram.from_formdiagram(form) index_uv = form.index_uv() ind = [3, 6, 10, 13, 16] for index in ind: u, v = index_uv[index] form.edge[u][v]['is_ind'] = True form.edge[u][v]['q'] = 1. graphstatics.form_update_q_from_qind(form) graphstatics.force_update_from_form(force, form) force.vertex[7]['x'] = 0 force.vertex[7]['y'] = 0 force.vertex[8]['x'] = 0 force.vertex[8]['y'] = 2.5 force.vertex[13]['x'] = 0 force.vertex[13]['y'] = -2.5 force.vertex[6]['x'] = -2 force.vertex[6]['y'] = 2.5 force.vertex[1]['x'] = -2 force.vertex[1]['y'] = -2.5
def compute_jacobian(form, force): r"""Compute the Jacobian matrix. The actual computation of the Jacobian matrix :math:`\partial \mathbf{X}^* / \partial \mathbf{X}` where :math:`\mathbf{X}` contains the form diagram coordinates in *Fortran* order (first all :math:`\mathbf{x}`-coordinates, then all :math:`\mathbf{y}`-coordinates) and :math:`\mathbf{X}^*` contains the force diagram coordinates in *Fortran* order (first all :math:`\mathbf{x}^*`-coordinates, then all :math:`\mathbf{y}^*`-coordinates). Parameters ---------- form : compas_ags.diagrams.formdiagram.FormDiagram The form diagram. force : compas_bi_ags.diagrams.forcediagram.ForceDiagram The force diagram. Returns ------- jacobian Jacobian matrix (2 * _vcount, 2 * vcount) """ # Update force diagram based on form form_update_q_from_qind(form) force_update_from_form(force, form) # -------------------------------------------------------------------------- # form diagram # -------------------------------------------------------------------------- vcount = form.number_of_vertices() k_i = form.key_index() leaves = [k_i[key] for key in form.leaves()] free = list(set(range(form.number_of_vertices())) - set(leaves)) vicount = len(free) edges = [(k_i[u], k_i[v]) for u, v in form.edges()] xy = array(form.xy(), dtype=float64).reshape((-1, 2)) ecount = len(edges) C = connectivity_matrix(edges, 'array') E = equilibrium_matrix(C, xy, free, 'array') uv = C.dot(xy) u = np.asmatrix(uv[:, 0]).transpose() v = np.asmatrix(uv[:, 1]).transpose() Ct = C.transpose() Cti = Ct[free, :] q = array(form.q(), dtype=float64).reshape((-1, 1)) Q = np.diag(np.asmatrix(q).getA1()) Q = np.asmatrix(Q) independent_edges = list(form.edges_where({'is_ind': True})) independent_edges_idx = [edges.index(i) for i in independent_edges] dependent_edges_idx = list(set(range(ecount)) - set(independent_edges_idx)) Ed = E[:, dependent_edges_idx] Eid = E[:, independent_edges_idx] qid = q[independent_edges_idx] # -------------------------------------------------------------------------- # force diagram # -------------------------------------------------------------------------- _vcount = force.number_of_vertices() _edges = force.ordered_edges(form) _L = laplacian_matrix(_edges, normalize=False, rtype='array') _C = connectivity_matrix(_edges, 'array') _Ct = _C.transpose() _Ct = np.asmatrix(_Ct) _k_i = force.key_index() _known = [_k_i[force.anchor()]] # -------------------------------------------------------------------------- # Jacobian # -------------------------------------------------------------------------- jacobian = np.zeros((_vcount * 2, vcount * 2)) for j in range(2): # Loop for x and y idx = list(range(j * vicount, (j + 1) * vicount)) for i in range(vcount): dXdxi = np.diag(Ct[i, :]) dxdxi = np.transpose(np.asmatrix(Ct[i, :])) dEdXi = np.zeros((vicount * 2, ecount)) dEdXi[idx, :] = np.asmatrix(Cti) * np.asmatrix( dXdxi) # Always half the matrix 0 depending on j (x/y) dEdXi_d = dEdXi[:, dependent_edges_idx] dEdXi_id = dEdXi[:, independent_edges_idx] EdInv = np.linalg.inv(np.asmatrix(Ed)) dEdXiInv = -EdInv * (np.asmatrix(dEdXi_d) * EdInv) dqdXi_d = -dEdXiInv * (Eid * np.asmatrix(qid)) - EdInv * ( dEdXi_id * np.asmatrix(qid)) dqdXi = np.zeros((ecount, 1)) dqdXi[dependent_edges_idx] = dqdXi_d dqdXi[independent_edges_idx] = 0 dQdXi = np.asmatrix(np.diag(dqdXi[:, 0])) d_XdXiTop = np.zeros((_L.shape[0])) d_XdXiBot = np.zeros((_L.shape[0])) if j == 0: d_XdXiTop = solve_with_known( _L, np.array(_Ct * (dQdXi * u + Q * dxdxi)).flatten(), d_XdXiTop, _known) d_XdXiBot = solve_with_known( _L, np.array(_Ct * (dQdXi * v)).flatten(), d_XdXiBot, _known) elif j == 1: d_XdXiTop = solve_with_known( _L, np.array(_Ct * (dQdXi * u)).flatten(), d_XdXiTop, _known) d_XdXiBot = solve_with_known( _L, np.array(_Ct * (dQdXi * v + Q * dxdxi)).flatten(), d_XdXiBot, _known) d_XdXi = np.hstack((d_XdXiTop, d_XdXiBot)) jacobian[:, i + j * vcount] = d_XdXi return jacobian