fe_eig

Contents

fe_eig#

Purpose

Computation of normal modes associated to a second order undamped model.

Syntax

def = fe_eig(model,EigOpt)
def = fe_eig({m,k,mdof},EigOpt)
def = fe_eig({m,k,T,mdof},EigOpt)
[phi, wj] = fe_eig(m,k)
[phi, wj, kd] = fe_eig(m,k,EigOpt,imode)

Description

The normal modeshapes phi=φ and frequencies wj=sqrt(diag(Ω2)) are solution of the undamped eigenvalue problem (see section 5.2)

 
    (10.13)

and verify the two orthogonality conditions

 
    (10.14)

The outputs are the data structure def (which is more appropriate for use with high level functions feplot, nor2ss, ... since it keeps track of the signification of its content, frequencies in def.data are then in Hz) or the modeshapes (columns of phi) and frequencies wj in rad/s. Note how you provide {m,k,mdof} in a cell array to obtain a def structure without having a model.

The optional output kd corresponds to the factored stiffness matrix. It should be used with methods that do not renumber DOFs.

fe_eig implements various algorithms to solve this problem for modes and frequencies. Many options are available and it is important that you read the notes below to understand how to properly use them. The option vector EigOpt can be supplied explicitely or set using model=stack_set(model, 'info','EigOpt',EigOpt). Its format is

[method nm Shift Print Thres] (default values are [2 0 0 0 1e-5])

  • method
    • 2 full matrix solution. Cannot be used for large models, used by default when the number of searched modes exceed 25% of the matrix size.
    • 5 default Lanczos solver is an iterative solver with problem size scalability and higher robustness with convergence checks. To turn off convergence check add 2000 to the option (2105, 2005, ...), otherwise a maximum of 5 convergence iterations is performed. You can tune this value by setting the 9th value of the opt vector to the desired number. opt=[5 100 1e3 1 0 0 0 0 1];.
    • 6 IRA/Sorensen solver. Faster than 5 but less robust, issues are known for multiple modes, very close frequencies, or when computing a large number of modes.
    • 50 Callback to let the user specify an external solver method using setpref('SDT','ExternalEig').
    • The other methods are left for reference but should not be used,
      • 105, 106, 104 same as 5, 6, 4 methods but no initial DOF renumbering. This is useless with the default ofact('methodspfmex') which renumbers at factorization time.
      • 0 SVD based full matrix solution
      • 1 subspace iteration which allows to compute the lowest modes of a large problem where sparse mass and stiffness matrices are used.
      • 3 Same as 5 but using ofact('methodlu').
      • 4 Same as 5 but using ofact('methodchol').
  • nm number of modes to be returned. A non-integer or negative nm is used as the desired fmax in Hz for iterative solvers (this is limited to 12 modes with method 5). The easiest way to handle fmax at the moment is to call fe_def SubDef after fe_eigṪhe sample syntax is then def=fe_eig(model,[5 50 1e3]); def=fe_def('subdef',def,find(def.data(:,1)<=fmax));. One thus has to estimate the relevant number of modes necessary beforehand.
  • shift value of mass shift (should be non-zero for systems with rigid body modes, see notes below). The subspace iteration method supports iterations without mass shift for structures with rigid body modes. This method is used by setting the shift value to Inf.
  • print level of printout (0 none, 11 maximum)
  • thres threshold for convergence of modes (default 1e-5 for the subspace iteration and Lanczos methods)

Finally, a set of vectors imode can be used as an initial guess for the subspace iteration method (method 1).

Notes

  • The default full matrix algorithm (method=2) cleans results of the MATLAB eig function. Computed modes are mass normalized and complex parts, which are known to be spurious for symmetric eigenvalue problems considered here, are eliminated. The alternate algorithm for full matrices (method=0) uses a singular value decomposition to make sure that all frequencies are real. The results are thus wrong, if the matrices are not symmetric and positive definite (semi-positive definite for the stiffness matrix).
  • The Lanczos algorithm (methods 3,4,5) is much faster than the subspace iteration algorithm (method 1). A double orthogonalization scheme and double restart usually detects multiple modes.
  • Method 6 calls eigs (ARPACK) properly and cleans up results. This solver sometimes fails to reach convergence, use method 5 then.
  • The subspace iteration and Lanczos algorithms are rather free interpretation of the standard algorithms (see Ref. [44] for example).
  • For systems with rigid body modes, you must specify a mass-shift. A good value is about one tenth of the first flexible frequency squared, but the Lanczos algorithm tends to be sensitive to this value (you may occasionally need to play around a little). If you do not find the expected number of rigid body modes, this is can be reason. For large frequency bands, consider using a shift at 75% of the largest estimated frequency, using −(0.75*2*pi*f_max)2.
  • DOFs with zero values on the stiffness diagonal are eliminated by default. You can bypass this behavior by giving a shift with unit imaginary value (eigopt(3)=1e3+1i for example).
  • For performance, optimization matters, please refer to section section 4.10.6.

Example#

Here is an example containing a high level call

 model =demosdt('demo gartfe');
 cf=feplot;cf.model=model;
 cf.def=fe_eig(model,[5 20 1e3 11]);
 fecom chc10

and the same example with low level commands

 model =demosdt('demo gartfe');
 [m,k,mdof] = fe_mknl(model);
 cf=feplot;cf.model=model;
 cf.def=fe_eig({m,k,mdof},[6 20 1e3]);fecom chc10

See also

fe_ceig, fe_mk, nor2ss, nor2xf