tmdsimpy.roms.vprnm.constant_h1_displacement¶
- tmdsimpy.roms.vprnm.constant_h1_displacement(epmc_fund_bb, h_fund, epmc_rhi_bb, h_rhi, vprnm_bb, h_vprnm, rhi, control_point_h1, control_amp_h1, control_point_rhi, Flcos, extra_Omega=None, correct_force=True)¶
Reduced order model (ROM) for a superharmonic resonance based VPRNM, EPMC, and on constant first harmonic amplitude.
- Parameters:
- epmc_fund_bb(Mfund, Nhc_fund*Ndof+3) numpy.ndarray
Each row corresponds to an EPMC solution at a given amplitude level for the mode that response at the forcing frequency (fundamental). The first Nhc*Ndof entries of each row are the displacements of the harmonics in h (all of the first harmonic component, then all of next etc.). The last three entries of each row are frequency (rad/s), xi in EPMC formulation (coefficient in front of mass matrix to create a damping matrix), and log10(modal amplitude). Harmonic displacements must be multiplied by the modal amplitude to get the physical displacements except for displacements corresponding to the zeroth harmonic.
- h_fundnumpy.ndarray, sorted
Array of the harmonics used to calculate epmc_fund_bb, must be sorted.
- epmc_rhi_bb(Mrhi, Nhc_rhi*Ndof+3) numpy.ndarray
EPMC solution for the mode that corresponds in superharmonic resonance. Nhc_rhi = utils.harmonic.Nhc(h_rhi).
- h_rhinumpy.ndarray, sorted
Array of the harmonics used to calculate epmc_rhi_bb, must be sorted.
- vprnm_bb(Mvprnm, Nhc_v*Ndof+4) or (Mvprnm, Nhc_v*Ndof+2) numpy.ndarray
VPRNM solution points. The first Nhc_vprnm*Ndof columns correspond to harmonic displacements in physical coordinates. If the second dimension is Nhc_v*Ndof+4, then the last 4 columns are force scaling for cosine, force scaling for sine, frequency (rad/s), controlled amplitude (amplitude and phase were controlled in VPRNM solution). If the second dimension is Nhc_v*Ndof+2, then the last 2 columns are frequency (rad/s) and force magnitude scaling (amplitude and phase not controlled in VPRNM solution). Nhc_v = utils.harmonic.Nhc(h_vprnm).
- h_vprnmnumpy.ndarray
Array of the harmonics used to calculate vprnm_bb, must be sorted.
- rhiint
The resonant superharmonic of interest (multiple that multiplies forcing frequency).
- control_point_h1(Ndof,) numpy.ndarray
Vector for extracting the degree of freedom that should be controlled to constant amplitude as control_point @ U(t) where U(t) is the Ndof displacement vector.
- control_amp_h1float
Desired amplitude that the control point should be controlled to for the first harmonic motion. This is displacement amplitude (0th derivative of motion)
- control_point_rhi(Ndof,) numpy.ndarray
Vector for extracting the degree of freedom that is considered to interpolate/estimate the amplitude of the superharmonic resonance.
- Flcos(Ndof,) numpy.ndarray
External constant excitation to be considered. Excitation is just applied to the first harmonic (at undetermined phase)
- extra_OmegaNone or 1D numpy.ndarray, optional
Extra forcing frequencies (rad/s) to calculate the force scaling magnitude at. Superharmonic responses are linearly interpolated from calculated points to these points. Responses at these frequencies are returned as well as directly calculated frequencies for the superharmonic resonance EPMC ROM. The default is None.
- correct_forcebool, optional
Flag to calculate force correction for VPRNM compared to the EPMC ROM for the first harmonic force. This correction is required to give good normalized responses for some cases. However, if extra_Omega includes values outside of what the superharmonic EPMC ROM identifies, settings this to False will allow for getting the EPMC ROM for the 1st mode at those points (and nan for the superharmonic) The default is True.
- Returns:
- Uw_reconstruct(Mout,Nhc_out*Ndof+1) numpy.ndarray
Reconstructed Frequency Response Curve (FRC) based on the VPRNM and EPMC solutions. Each row corresponds to a different forcing frequency. The first Nhc_out*Ndof columns are harmonic displacements in physical coordinates corresponding to the harmonics in h_reconstruct. The last column is the forcing frequency in rad/s. Nhc_out = utils.harmonic.Nhc(h_reconstruct).
- force_reconstruct(Mout,) numpy.ndarray
External excitation force magnitude (scaling for Flcos) to create the response. Phase information of this force is not calculated.
- h_reconstructnumpy.ndarray
List of sorted harmonics for the output. This list includes harmonics in h_fund and rhi*h_rhi.
See also
tmdsimpy.roms.vprnmModule includes more details and discussion on VPRNM based ROMs
tmdsimpy.roms.epmc.constant_forceConstant force EPMC based ROM for a single mode
tmdsimpy.roms.epmc.constant_displacementConstant amplitude EPMC based ROM for a single mode.
tmdsimpy.VibrationSystem.epmc_resEPMC residual method for finding EPMC backbones.
tmdsimpy.VibrationSystem.vprnm_resVPRNM residual method for finding VPRNM backbone.
tmdsimpy.VibrationSystem.hbm_amp_control_resHBM residual method calculates something similar to a truth solution to compare this ROM to.
tmdsimpy.Continuation.continuationMethod of obtaining solutions to EPMC and VPRNM at multiple points to create inputs to this function
Notes
Formulation of this ROM is derived and explained in detail in [1], [2].
control_point_h1 and control_point_rhi need not be the same, but can be the same. Depending on mode shapes it may make sense to use different locations to approximate the amplitudes of the two modes used.
In general, epmc_fund_bb should be calculated with h_fund not including rhi.
Static (harmonic 0) displacements in Uw_reconstruct are taken solely from the VPRNM solution. If vprnm is calculated without a zeroth harmonic and a zeroth harmonic is included in one or both of the EPMC solutions, then zeros are returned for the static displacements ignoring static components of the EPMC solutions.
In general, specific response frequencies cannot be requested because the resonance of the superharmonic provides the frequency as a direct calculation while requesting a specific frequency may require a nonlinear solution. The parameter extra_Omega allows for requesting higher resolution for the output force scaling (which can be directly calculated at requested frequencies). For these requested frequencies, the superharmonic response is simply linearly interpolated to fill in the points.
For any frequencies in extra_Omega lower than the lowest frequency that this function calculates without extra_Omega, nan values will be returned for parts of the solution.
The superharmonic EPMC backbone is abbreviated to include the highest amplitude equal to the VPRNM amplitude of the superharmonic resonance. However, for a linear system, the peak amplitude is slightly more than the amplitude at phase resonance. This slightly higher amplitude is not currently captured, but could be implemented in the future (at least for the linear case). For the nonlinear case, it may be nontrivial to interpolate the exact resonance amplitude, but extra EPMC points could be allowed in the superharmonic part of the ROM to try to improve resolution rather than the current approach that eliminates any points greater than the VPRNM amplitude.
References
[1]Porter, J. H., and M. R. W. Brake. Under Review. “Efficient Model Reduction and Prediction of Superharmonic Resonances in Frictional and Hysteretic Systems.” Mechanical Systems and Signal Processing. arXiv:2405.15918.
[2]Porter, J. H. 2024. Modal Interactions and Jointed Structures. PhD Thesis. Rice University.