amfe.hyper_red subpackage

training_set_generation

Training set generation

amfe.hyper_red.training_set_generation.compute_nskts(mechanical_system, F_ext_max=None, no_of_moments=4, no_of_static_cases=8, load_factor=2, no_of_force_increments=20, no_of_procs=None, norm='impedance', verbose=True, force_basis='krylov')[source]

Compute the Nonlinear Stochastic Krylov Training Sets (NSKTS).

NSKTS can be used as training sets for Hyper Reduction of nonlinear systems.

Parameters:
  • mechanical_system (instance of amfe.MechanicalSystem) – Mechanical System to which the NSKTS should be computed
  • F_ext_max (ndarray, optional) – Maximum external force. If None is given, 10 random samples of the external force in the time range t = [0,1] are computed and the maximum value is taken. Default value is None.
  • no_of_moments (int, optional) – Number of moments consiedered in the Krylov force subspace. Default value is 4.
  • no_of_static_cases (int, optional) – Number of stochastic static cases which are solved. Default value is 8.
  • load_factor (int, optional) – Load amplification factor with which the maximum external force is multiplied. Default value is 2.
  • no_of_force_increments (int, optional) – Number of force increments for nonlinear solver. Default value is 20.
  • no_of_procs ({int, None}, optional) – Number of processes which are started parallel. For None no_of_static_cases processes are run.
  • norm (str {'impedance', 'eucledian', 'kinetic'}, optional) – Norm which will be used to scale the higher order moments for the Krylov force subspace. Default value is ‘impedance’.
  • verbose (bool, optional) – Flag for setting verbose output. Default value is True.
  • force_basis (str {'krylov', 'modal'}, optional) – Type of force basis used. Either krylov meaning the classical NSKTS or modal meaning the forces producing vibration modes.
Returns:

  • nskts_arr (ndarray) – Nonlinear Stochastic Krylov Training Sets. Every column in nskts_arr represents one NSKTS displacement field.
  • Reference
  • ———
  • Todo

amfe.hyper_red.training_set_generation.krylov_force_subspace(M, K, b, omega=0, no_of_moments=3, orth='euclidean')[source]

Compute a krylov force subspace for the computation of snapshots needed in hyper reduction.

The Krylov force basis is given as

..code:

[b, M @ inv(K - omega**2) @ b, ...,
 (M @ inv(K - omega**2))**(no_of_moments-1) @ b]
Parameters:
  • M (ndarray) – Mass matrix of the system.
  • K (ndarray) – Stiffness matrix of the system.
  • b (ndarray) – input vector of external forcing.
  • omega (float, optional) – frequency for the frequency shift of the stiffness. Default value 0.
  • no_of_moments (int, optional) – number of moments matched. Default value 3.
  • orth (str, {'euclidean', 'impedance', 'kinetic'} optional) –

    flag for setting orthogonality of returnd Krylov basis vectors.

    • ‘euclidean’ : V.T @ V = eye
    • ‘impedance’ : V.T @ inv(K) @ V = eye
    • ‘kinetic’ : V.T @ inv(K).T @ M @ inv(K) @ V = eye
Returns:

V – Krylov force basis where vectors V[:,i] give the basis vectors.

Return type:

ndarray

amfe.hyper_red.training_set_generation.modal_force_subspace(M, K, no_of_modes=3, orth='euclidean')[source]

Force subspace spanned by the forces producing the vibration modes.

DEIM

Discrete Empirical Interpolation Method (DEIM) hyper-reduction

class amfe.hyper_red.deim.DEIMSystem(*args, **kwargs)[source]

Bases: amfe.mechanical_system.ReducedSystem

Hyper-reduction technique inherited from ReducedSystem class. Currently only DEIM based on POD is implemented.

This class can handle UDEIM, SUDEIM and their symmetric variants and also the other variants of the above based on choosing one dof or nodal dofs or the entire elements dofs (refered to as collocaiton variants).

DEIM_type

str, – String indicating the DEIM type.

E_tilde

ndarray, ndim (no_of_active_elements) – Indices of all active elements

oblique_proj

ndarray, ndim (no_of_dofs, no_of_force_modes) – oblique projection matrix formed by V.T @ U @ inv(P.T @ U) This matrix is compact and carries the oblique projection of the nonlinear force after the collocation.

K0_deim

ndarray, ndim (nred, nred) – Reduced linear stiffness matrix which corrects the DEIM nonlinear force and tangential stiffness matrix. It is computed as the difference of the linear stiffness matrix minus the linear effects coming from the DEIM procedure. Hence, this matrix has to be added to the nakedly assembled DEIM tangential stiffness matrix to obtain the tangential stiffness matrix of DEIM sytem.

proj_list

ndarray, ndim (no_of_active_elements, no_of_dofs, – no_of_dofs_per_element) This tensor is used as a list. It contains all oblique projectors of the active elements onto the kinematic subspace. The i-th element has the projector proj_list[i] as oblique projector to give the contribution of the nonlinear force.

Note

The key approximation in DEIM is the approximation of the full nonlinear force vector f_int through the unassembled force vector f_u:

\[\begin{split}f \approx U (P^T U)^{-1} P^T f_u \\ f_{red} \approx \underbrace{V^T U (P^T U)^{-1}P^T}_{X} f_u\end{split}\]

It is important to note, that the nonlinear force given here is a split-off part of the nonlinear restoring force. Here the nolinear force is the force without the linear part.

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

Compute tangential stiffness matrix and nonlinear force vector with the help of the computed interpolation and selective evaluation based on the selected DEIM elements

Parameters:
  • u (ndarray) – Initial solution or solution from previous iteration
  • t (float) – Time
Returns:

  • K (ndarray) – Hyper-reduced stiffness matrix of reduced system
  • f_int (ndarray) – Hyper-reduced internal force of reduced system

collocate_DEIM(U_f)[source]

Build the collocation matrix P for the DEIM algorithm

collocate_UDEIM(U_f_u)[source]

Run collocation with UDEIM algorithms

export_paraview(filename, field_list=None)[source]

Export the produced results to ParaView via XDMF format.

force_basis(U_snapshots_unconstr, no_of_force_modes, unassembled=True)[source]

Compute the unassembled or assembled force modes based on an SVD of the force snapshots generated by the displacement vectors.

Parameters:
  • U_snapshots (ndarray) – Training snapshots of the system (Unconstrained full solution)
  • no_of_force_modes (int, optional) – First ‘m’ number of force modes to be chosen
  • unassembled (bool, optional) – Flag setting, if unassembled (True) or assembled (False) force basis is computed
Returns:

U_f – Force mode either unassembled or assembled depending on the option

Return type:

ndarray

preprocess_DEIM(DEIM_type='unassem-deim-dof')[source]

Setting the flags necessary for the force basis and collocation procedure.

Parameters:DEIM_type (str, optional) – string indicating the technique used for DEIM reduction
reduce_mesh(U_snapshots, no_of_force_modes=10, DEIM_type='unassem-deim-dof')[source]

Reduce the mesh using an interpolation gathered from the snapshots.

Parameters:
  • U_snapshots (ndarray) – Constrained training snapshots of the system
  • no_of_force_modes (int, optional) – First ‘m’ number of force modes to be chosen
  • DEIM_type (str, optional) –

    String indicating the technique used for DEIM reduction. The string can be composed of different keywords indicating the DEIM technique:

    • collocation_type: {‘dof’, ‘ele’, ‘node’, ‘x’, ‘y’, ‘z’}
      The element is collocated only in the selected direction (‘dof’), all dofs of the whole element are collocated (‘ele’), all dofs of the node are chosen (‘node’) or a direction is additionally chosen (‘x’, ‘y’, ‘z’)
    • symmetric: {‘symm’}
      If this keyword is in the DEIM_type, the symmetric DEIM technique is applied.
    • surrogate: {‘surr’}
      If this keyword is there, a surrogate technique for avoiding a large SVD is used. This is the so-called surrogate DEIM.

    The string can then be composed in arbitrary order, ie. ‘dof-symm-surr’ to make a surrogate symemtric DEIM with dof collocation.

Returns:

Return type:

None

References

DEIM: Chaturantabut, S. and Sorensen, D.C., 2010. Nonlinear model reduction via discrete empirical interpolation. SIAM Journal on Scientific Computing, 32(5), pp.2737-2764.

UDEIM, SUDEIM: Tiso, P. and Rixen, D.J., 2013. Discrete empirical interpolation method for finite element structural dynamics. In Topics in Nonlinear Dynamics, Volume 1 (pp. 203-212). Springer New York.

symmetric DEIM: Chaturantabut, S., Beattie, C. and Gugercin, S., 2016. Structure-Preserving Model Reduction for Nonlinear Port-Hamiltonian Systems. arXiv preprint arXiv:1601.00527.

symmetric UDEIM: Ravichandran, T.K., Investigation of accuracy, speed and stability of hyper-reduction techniques for nonlinear FE. Masters Thesis. TU Delft, Delft University of Technology, 2016.

amfe.hyper_red.deim.reduce_mechanical_system_deim(mechanical_system, V, overwrite=False, assembly='indirect')[source]

Reduce the given mechanical system with the linear basis V and the given weights.

Parameters:
  • mechanical_system (instance of MechanicalSystem) – Mechanical system which will be transformed to a ReducedSystem.
  • V (ndarray, shape (N_constrained, n_red)) – Reduction Basis for the reduced system
  • overwrite (bool, optional) – switch, if mechanical system should be overwritten (is less memory intensive for large systems) or not.
Returns:

reduced_system – Reduced system with same properties of the mechanical system and reduction basis V

Return type:

instance of ReducedSystem

Example

ECSW

TODO: Write introduction to ECSW

class amfe.hyper_red.ecsw.ECSWSystem(**kwargs)[source]

Bases: amfe.mechanical_system.ReducedSystem

Hyper Reduced system using ECSW for the redcution.

K(u=None, t=0)[source]
K_and_f(u=None, t=0)[source]
export_paraview(filename, field_list=None)[source]
f_int(u, t=0)[source]
reduce_mesh(W_red, tau=0.001, verbose=True)[source]

Compute a reduced mesh using a sparse NNLS solver for gaining the weights.

Parameters:
  • W_red (ndarray, shape(n_red, no_of_snapshots)) – Snapshot training matrix for which the energy equality is ensured.
  • tau (float, optional) – tolerance for fitting the best solution
Returns:

  • weight_indices (ndarray) – indices of members in the reduced mesh
  • weights (ndarray) – weights for reduced mesh
  • stats (ndarray) – statistics about the convergence of the sparse NNLS solver. The first column shows the size of the active set, the secont column the residual. If verbose=False, an empty array is returned.

Note

The indices and weights of the reduced mesh are also internally saved.

amfe.hyper_red.ecsw.reduce_mechanical_system_ecsw(mechanical_system, V, overwrite=False, assembly='indirect')[source]

Reduce the given mechanical system with the linear basis V and the given weights.

Parameters:
  • mechanical_system (instance of MechanicalSystem) – Mechanical system which will be transformed to a ReducedSystem.
  • V (ndarray, shape (N_constrained, n_red)) – Reduction Basis for the reduced system
  • overwrite (bool, optional) – switch, if mechanical system should be overwritten (is less memory intensive for large systems) or not.
  • assembly (str {'direct', 'indirect'}) – flag setting, if direct or indirect assembly is done. For larger reduction bases, the indirect method is much faster.
Returns:

reduced_system – Reduced system with same properties of the mechanical system and reduction basis V

Return type:

instance of ReducedSystem

Example

amfe.hyper_red.ecsw.sparse_nnls(G, b, tau, conv_stats=False, verbose=True)[source]

Run the sparse NNLS-solver in order to find a sparse vector xi satisfying

\[|| G \xi - b ||_2 \leq \tau ||b||_2 \quad\text{with}\quad \min||\xi||_0\]
Parameters:
  • G (ndarray, shape: (n*m, no_of_elements)) – force contribution matrix
  • b (ndarray, shape: (n*m)) – force contribution vector
  • tau (float) – tolerance
  • conv_info (bool) – Flag for setting, that more detailed output is produced with convergence information.
Returns:

  • indices (ndarray, shape: (k,)) – The indices of the non-zero elements.
  • xi_red (ndarray, shape: (k,)) – The values of the non-zero elements.
  • stats (ndarray) – Infos about the convergence of the system. The first column shows the size of the active set, the second column the residual. If conv_info is set to False, an empty array is returned.

References

[R1]C. L. Lawson and R. J. Hanson. Solving least squares problems, volume 15. SIAM, 1995.
[R2]T. Chapman, P. Avery, P. Collins, and C. Farhat. Accelerated mesh sampling for the hyper reduction of nonlinear computational models. International Journal for Numerical Methods in Engineering, 2016.