amfe.assembly module

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.

class 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:
  • V (ndarray, shape (N_unconstr, n)) – Reduction basis
  • S (ndarray, shape (N_unconstr, m)) – Snapshots gathered as column vectors.
  • verbose (print dots for the) –
Returns:

  • G (ndarray, shape (n*m, no_of_elements)) – Contribution matrix of internal forces. The columns form the internal force contributions on the basis V for the m snapshots gathered in S.
  • b (ndarray, shape (n*m, )) – summed force contribution

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:
  • u (ndarray) – nodal displacement of the nodes in Voigt-notation
  • t (float) – time
Returns:

  • K (sparse.csr_matrix) – unconstrained assembled stiffness matrix in sparse matrix csr format.
  • f (ndarray) – unconstrained assembled force vector

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:
  • V (ndarray, shape: (N_unconstr, n_red)) – unconstrained reduction basis
  • idxs (ndarray, shape (n_ele_hyper, )) – indices of the elements in the active element set.
  • xi (ndarray, shape (n_ele_hyper, )) – array of weights for the hyper reduced system
  • u (ndarray, shape: (n_red,)) – unconstrained displacement field
  • t (ndarray) – current time
Returns:

  • K (ndarray, shape: (n_red, n_red)) – reduced stiffness matrix
  • f_int (ndarray, shape (n_red,)) – reduced internal force vector

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:
  • V (ndarray, shape: (N_unconstr, n_red)) – unconstrained reduction basis
  • idxs (ndarray, shape (n_ele_hyper, )) – indices of the elements in the active element set.
  • xi (ndarray, shape (n_ele_hyper, )) – array of weights for the hyper reduced system
  • u (ndarray, shape: (n_red,)) – unconstrained displacement field
  • t (ndarray) – current time
Returns:

  • K (ndarray, shape: (n_red, n_red)) – reduced stiffness matrix
  • f_int (ndarray, shape (n_red,)) – reduced internal force vector

assemble_k_and_f_neumann(u=None, t=0)[source]

Assembles the stiffness matrix and the force of the Neumann skin elements.

Parameters:
  • u (ndarray) – nodal displacement of the nodes in Voigt-notation
  • t (float) – time
Returns:

  • K (sparse.csr_matrix) – unconstrained assembled stiffness matrix in sparse matrix csr format.
  • f (ndarray) – unconstrained assembled force vector

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:
  • V (ndarray, shape: (N_unconstr, n_red)) – unconstrained reduction basis
  • u (ndarray, shape: (n_red,)) – unconstrained displacement field
  • t (ndarray) – current time
Returns:

  • K (ndarray, shape: (n_red, n_red)) – reduced stiffness matrix
  • f_int (ndarray, shape (n_red,)) – reduced internal force vector

assemble_k_f_S_E(u, t)[source]

Assemble the stiffness matrix with stress recovery of the given mesh and element.

Parameters:
  • u (ndarray) – nodal displacement of the nodes in Voigt-notation
  • t (float) – time
Returns:

  • K (sparse.csr_matrix) – unconstrained assembled stiffness matrix in sparse matrix csr format.
  • f (ndarray) – unconstrained assembled force vector
  • S (ndarray) – unconstrained assembled stress tensor
  • E (ndarray) – unconstrained assembled strain tensor

assemble_m(u=None, t=0)[source]

Assembles the mass matrix of the given mesh and element.

Parameters:
  • u (ndarray or None, optional) – nodal displacement of the nodes in Voigt-notation. Default: None.
  • t (float, optional) – time. Default: 0.
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:
  • u (ndarray) – Initial solution or solution from previous iteration
  • t (float) – Time
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...

preallocate_hyper_csr(idxs)[source]

Preallocate a sparse matrix for a reduced mesh.

Parameters:idxs (ndarray, shape (n_ele_hyper, )) – indices of the elements in the active element set.
Returns:
Return type:None