Basic assembly module for the finite element code. Assumes to have all elements in the inertial frame.
Provides an Assembly class which knows the mesh. It can assemble the vector of nonlinear forces, the mass matrix and the tangential stiffness matrix. Some parts of the code – mostly the indexing of the sparse matrices – are substituted by fortran routines, as they allow for a huge speedup.
amfe.assembly.
Assembly
(mesh)[source]¶Bases: object
Class for the more fancy assembly of meshes with non-heterogeneous elements.
C_csr
¶scipy.sparse.csr.csr_matrix – Matrix containing the sparsity pattern of the problem
element_indices
¶list – Ragged list containing the global indices for the local variables of an element. The entry [i,j] gives the index in the global vector of element i with dof j
mesh
¶amfe.Mesh – Mesh-Class the Assembly is associated with
neumann_indices
¶list – Ragged list equivalently to element_indices for the neumann boundary skin elements.
nodes_voigt
¶np.ndarray – vector of all nodal coordinates in voigt-notation. Dimension is (ndofs_total, )
elements_on_node
¶np.ndarray – array containing the number of adjacent elements per node. Is necessary for stress recovery.
assemble_g_and_b
(V, S, verbose=False)[source]¶Assembles the element contributin matrix G for the given basis V and the snapshots S.
This function is needed for cubature bases Hyper reduction methods like the ECSW.
Parameters: |
|
---|---|
Returns: |
|
Note
This assembly works on the unconstrained variables, so both, V and S have to contain the unassembled nodal displacements.
assemble_k_and_f
(u, t)[source]¶Assembles the stiffness matrix of the given mesh and element.
Parameters: |
|
---|---|
Returns: |
|
assemble_k_and_f_DEIM
(E_tilde, proj_list, V, u_red=None, t=0, symmetric=False)[source]¶Assemble the DEIM matrix and vector. This function is suited for large systems, as only hte proj list is used.
assemble_k_and_f_hyper
(V, idxs, xi, u, t)[source]¶Assembly routine for a hyper reduced ECSW system.
Parameters: |
|
---|---|
Returns: |
|
See also
assemble_k_and_f_hyper_no_inplace
(V, idxs, xi, u, t)[source]¶Assembly routine for a hyper reduced ECSW system.
The difference to the function assemble_k_and_f_hyper is, that here the sparse matrix is assembled first and then the reduction is performed.
Parameters: |
|
---|---|
Returns: |
|
See also
assemble_k_and_f_neumann
(u=None, t=0)[source]¶Assembles the stiffness matrix and the force of the Neumann skin elements.
Parameters: |
|
---|---|
Returns: |
|
assemble_k_and_f_red
(V, u, t)[source]¶Assembly routine for reduces systems. Note, that V has to be unconstrained and u is the reduced displacement field.
Parameters: |
|
---|---|
Returns: |
|
assemble_k_f_S_E
(u, t)[source]¶Assemble the stiffness matrix with stress recovery of the given mesh and element.
Parameters: |
|
---|---|
Returns: |
|
assemble_m
(u=None, t=0)[source]¶Assembles the mass matrix of the given mesh and element.
Parameters: |
|
---|---|
Returns: | M – unconstrained assembled mass matrix in sparse matrix csr-format. |
Return type: | sparse.csr_matrix |
Examples
TODO
compute_c_deim
()[source]¶Assembles the matrix C_deim which assembles the unassembled system C_matrix: f_a = C_deim @ f_u
None
compute_element_indices
()[source]¶Compute the element indices which are the global dofs of every element.
The element_indices are a list, where every element of the list denotes the dofs of the element in the correct order.
Parameters: | None – |
---|---|
Returns: | |
Return type: | None |
f_nl_unassembled
(u=None, t=0)[source]¶Computes the unassemebled nonlinear force for DEIM.
Parameters: |
|
---|---|
Returns: | f – Unassembled nonlinear internal force |
Return type: | ndarray |
preallocate_csr
()[source]¶Compute the sparsity pattern of the assembled matrices and store an empty matrix in self.C_csr.
The matrix self.C_csr serves as a ‘blueprint’ matrix which is filled in the assembly process.
Parameters: | None – |
---|---|
Returns: | |
Return type: | None |
Notes
This preallocation routine can take some while for large matrices. Furthermore it is not implemented memory-efficient, so for large systems and low RAM this might be an issue...