Example C3: Response Constrained Forced Response of a 2DoF Hysteretic Oscillator
It is sometimes the case that it is necessary to conduct forced response simulations where the response amplitude is kept fixed. In fact, experimentally it is significantly easier to achieve response-controlled forced response measurements (in the context of, say, stepped sine testing) than force-controlled measurements. While there are several reasons for this (large displacement effects on the exciter at resonance, low SNR at anti-resonance points, etc.), we will not go into this here.
MDOFUTILS provides HBRESFUN_A!, a routine that returns the harmonic balance residue along with an amplitude constraint. The usage of this is only slightly different from HBRESFUN! (which we have already encountered in Examples C1 and C2 for forced response analysis).
The same system as before will be studied here:
\[\underbrace{\begin{bmatrix} 1&0\\0&1 \end{bmatrix}}_{\mx{M}} \begin{bmatrix} \ddot{x_1}\\ \ddot{x_2} \end{bmatrix} + \left( 0.01 \mx{M} + 0.001 \mx{K}\right) \begin{bmatrix} \dot{x_1}\\ \dot{x_2} \end{bmatrix} + \underbrace{\begin{bmatrix} 2&-1\\-1&2 \end{bmatrix}}_{\mx{K}} \begin{bmatrix} x_1\\x_2 \end{bmatrix} + \begin{bmatrix} 0\\ f_{fr}(x_2) \end{bmatrix} = \begin{bmatrix} 1\\0 \end{bmatrix}\cos\Omega t.\]
Preamble: Load Packages
using GLMakie
using LinearAlgebra
using SparseArrays
using ForwardDiff
using NonlinearSolve
using DSP
using juliajim.HARMONIC
using juliajim.CONTINUATION
using juliajim.MDOFUTILSSystem Setup
MCK matrices and nonlinearity setup: identical to Example C2.
M = [1. 0.;0. 1.];
K = [2. -1.;-1. 2.];
C = 0.01*M+0.001*K;
mdl = MDOFGEN(M, C, K);
# Nonlinearity
kt = 1.0;
fs = 1.0;
fnl = (t,u,up,fp)-> if all(abs.(fp+kt*(u-up)).<fs)
return fp+kt*(u-up), kt*ones(size(u)), -kt*ones(size(u)), ones(size(u)); else
return fs*sign.(fp+kt*(u-up)), zeros(size(u)), zeros(size(u)), zeros(size(u));
end
L = [0. 1.];
mdl = ADDNL(mdl, :Hyst, fnl, L);Setup HB
HB setup: identical to Example C2.
h = 1:2:5;
# h = 0:5;
N = 256;
t = (0:N-1)*2π/N;
Nhc = sum(all(h.==0, dims=2) + 2any(h.!=0, dims=2));
_, _, zinds, rinds, iinds = HINDS(mdl.Ndofs, h)
Fl = zeros(Nhc*mdl.Ndofs);
Fl[rinds[1]] = 1.0;
Wst = 0.6;
E, dEdw = HARMONICSTIFFNESS(mdl.M, mdl.C, mdl.K, Wst, h);
U0 = E\Fl;Test HB Residue
The residue function HBRESFUN_A! expects a vector of the form (\begin{bmatrix} \vc{u}& F& \Omega \end{bmatrix}^T), where (\vc{u}) is the vector of harmonics, (F) is the (scalar) excitation amplitude (also an unknown), and (\Omega) is the (scalar) excitation frequency. Apart from the same arguments as HBRESFUN!, the required response amplitude Amp must also be provided here. By default, this amplitude is interpreted as the first harmonic amplitude of the first DoF. This can be customized by supplying optional keywords atype (one of :H1 or :RMS, specifying type of amplitude measure) and ashape (either DoF index or vector of weights to DoFs); see the documentation. Here we set this up using NonlinearSolve.jl and solve for a single point.
R = zeros(mdl.Ndofs*Nhc+1);
dRdUf = zeros(mdl.Ndofs*Nhc+1, mdl.Ndofs*Nhc+1);
dRdw = zeros(mdl.Ndofs*Nhc+1);
Amp = 1e0;
dof = 1;
HBRESFUN_A!([U0; 1.0; Wst], mdl, Amp, Fl, h, N; R=R, dRdUf=dRdUf, dRdw=dRdw)
fun = NonlinearFunction((r,uf,p)->HBRESFUN_A!([uf;p], mdl, Amp, Fl, h, N; R=r),
jac=(J,uf,p)->HBRESFUN_A!([uf;p], mdl, Amp, Fl, h, N; dRdUf=J),
paramjac=(Jp,uf,p)->HBRESFUN_A!([uf;p], mdl, Amp, Fl, h, N; dRdw=Jp));
prob = NonlinearProblem(fun, [U0;1.0], Wst);
sol = solve(prob, show_trace=Val(true));
Algorithm: NewtonRaphson(
descent = NewtonDescent(),
autodiff = AutoForwardDiff(),
vjp_autodiff = AutoFiniteDiff(
fdtype = Val{:forward}(),
fdjtype = Val{:forward}(),
fdhtype = Val{:hcentral}(),
dir = true
),
jvp_autodiff = AutoForwardDiff(),
concrete_jac = Val{false}()
)
---- ------------- -----------
Iter f(u) inf-norm Step 2-norm
---- ------------- -----------
0 5.91774639e-01 0.00000000e+00
1 8.98242782e-04 3.38756045e-01
2 2.01529002e-07 7.42307628e-04
3 1.02140518e-14 1.66618285e-07
Final 1.02140518e-14
----------------------Continuation
Now we continue as before with continuation. In fact this piece of code is exactly the same as in Example C2.
Om0 = 3.0;
Om1 = 0.1;
dOm = 0.5;
# Use the unscaled version if scaling doesn't work.
# cpars = (parm=:arclength, nmax=2000, Dsc=:none, itopt=4); # Autoscaling is not working
cpars = (parm=:arclength, nmax=2000, Dsc=:auto, angopt=deg2rad(30), itopt=4); # Autoscaling is not working
sols, its, dss, xis, Dsc = CONTINUATE([U0;1.0], fun, [Om0, Om1], dOm; cpars...);
uh = zeros(Complex, maximum(h)+1, length(sols), 2);
for i in 1:2
uh[h.+1, :, i] = hcat([[s.up[zinds[i:2:end]];
s.up[rinds[i:2:end]]+1im*s.up[iinds[i:2:end]]]
for s in sols]...);
end
Oms = sols.p;
Fs = [s.up[end-1] for s in sols];
Algorithm: NewtonRaphson(
descent = NewtonDescent(),
autodiff = AutoForwardDiff(),
vjp_autodiff = AutoFiniteDiff(
fdtype = Val{:forward}(),
fdjtype = Val{:forward}(),
fdhtype = Val{:hcentral}(),
dir = true
),
jvp_autodiff = AutoForwardDiff(),
concrete_jac = Val{false}()
)
---- ------------- -----------
Iter f(u) inf-norm Step 2-norm
---- ------------- -----------
0 8.38540145e+00 4.90676743e-02
1 1.09954902e-03 7.87390001e+00
2 3.01920033e-07 3.79589072e-03
Final 3.01920033e-07
----------------------
2. 2.97 with step -0.0287 (0.0399) converged in 1 iterations. 0.306865 Δα.
3. 2.91 with step -0.0632 (0.0660) converged in 2 iterations. 0.423962 Δα.
4. 2.83 with step -0.0768 (0.0833) converged in 2 iterations. 0.425465 Δα.
5. 2.74 with step -0.0908 (0.1051) converged in 2 iterations. 0.419631 Δα.
6. 2.64 with step -0.1040 (0.1330) converged in 2 iterations. 0.403370 Δα.
7. 2.53 with step -0.1146 (0.1699) converged in 2 iterations. 0.372251 Δα.
8. 2.43 with step -0.1088 (0.1996) converged in 2 iterations. 0.326021 Δα.
9. 2.35 with step -0.0901 (0.1996) converged in 2 iterations. 0.298637 Δα.
10. 2.28 with step -0.0773 (0.1996) converged in 2 iterations. 0.282363 Δα.
11. 2.22 with step -0.0673 (0.1996) converged in 2 iterations. 0.270570 Δα.
12. 2.16 with step -0.0588 (0.1996) converged in 2 iterations. 0.259629 Δα.
13. 2.11 with step -0.0511 (0.1996) converged in 2 iterations. 0.247266 Δα.
14. 2.07 with step -0.0436 (0.1996) converged in 2 iterations. 0.231256 Δα.
15. 2.04 with step -0.0363 (0.1996) converged in 2 iterations. 0.208826 Δα.
16. 2.02 with step -0.0289 (0.1996) converged in 2 iterations. 0.178034 Δα.
17. 2.00 with step -0.0220 (0.1996) converged in 3 iterations. 0.086990 Δα.
18. 2.00 with step -0.0026 (0.1996) converged in 1 iterations. 0.084397 Δα.
19. 1.99 with step -0.0029 (0.1996) converged in 2 iterations. 0.082151 Δα.
20. 1.99 with step -0.0029 (0.1996) converged in 1 iterations. 0.059838 Δα.
21. 1.99 with step -0.0032 (0.1996) converged in 1 iterations. 0.107583 Δα.
22. 1.98 with step -0.0033 (0.1996) converged in 1 iterations. 0.059009 Δα.
23. 1.98 with step -0.0036 (0.1996) converged in 1 iterations. 0.090222 Δα.
24. 1.98 with step -0.0037 (0.1996) converged in 1 iterations. 0.013927 Δα.
25. 1.97 with step -0.0041 (0.1996) converged in 1 iterations. 0.058424 Δα.
26. 1.97 with step -0.0043 (0.1996) converged in 1 iterations. 0.036204 Δα.
27. 1.96 with step -0.0046 (0.1996) converged in 1 iterations. 0.028656 Δα.
28. 1.96 with step -0.0050 (0.1996) converged in 2 iterations. 0.036686 Δα.
29. 1.95 with step -0.0053 (0.1996) converged in 1 iterations. 0.027079 Δα.
30. 1.95 with step -0.0058 (0.1996) converged in 2 iterations. 0.059840 Δα.
31. 1.94 with step -0.0061 (0.1996) converged in 2 iterations. 0.027133 Δα.
32. 1.94 with step -0.0066 (0.1996) converged in 2 iterations. 0.056641 Δα.
33. 1.93 with step -0.0069 (0.1996) converged in 2 iterations. 0.031351 Δα.
34. 1.92 with step -0.0074 (0.1996) converged in 2 iterations. 0.063410 Δα.
35. 1.91 with step -0.0075 (0.1996) converged in 2 iterations. 0.232716 Δα.
36. 1.91 with step -0.0079 (0.1996) converged in 2 iterations. 0.249988 Δα.
37. 1.90 with step -0.0077 (0.1996) converged in 2 iterations. 0.214593 Δα.
38. 1.89 with step -0.0070 (0.1996) converged in 2 iterations. 0.202705 Δα.
39. 1.89 with step -0.0061 (0.1996) converged in 2 iterations. 0.072448 Δα.
40. 1.88 with step -0.0052 (0.1996) converged in 2 iterations. 0.094164 Δα.
41. 1.88 with step -0.0041 (0.1996) converged in 2 iterations. 0.025058 Δα.
42. 1.87 with step -0.0033 (0.1996) converged in 2 iterations. 0.012881 Δα.
43. 1.87 with step -0.0025 (0.1996) converged in 1 iterations. 0.026998 Δα.
44. 1.87 with step -0.0020 (0.1996) converged in 1 iterations. 0.033316 Δα.
45. 1.87 with step -0.0017 (0.1996) converged in 1 iterations. 0.037251 Δα.
46. 1.87 with step -0.0017 (0.1996) converged in 1 iterations. 0.034179 Δα.
47. 1.86 with step -0.0017 (0.1996) converged in 1 iterations. 0.037312 Δα.
48. 1.86 with step -0.0017 (0.1996) converged in 1 iterations. 0.035107 Δα.
49. 1.86 with step -0.0017 (0.1996) converged in 1 iterations. 0.034950 Δα.
50. 1.86 with step -0.0017 (0.1996) converged in 1 iterations. 0.035042 Δα.
51. 1.86 with step -0.0017 (0.1996) converged in 1 iterations. 0.024433 Δα.
52. 1.86 with step -0.0017 (0.1996) converged in 1 iterations. 0.012884 Δα.
53. 1.85 with step -0.0017 (0.1996) converged in 1 iterations. 0.001505 Δα.
54. 1.85 with step -0.0018 (0.1996) converged in 1 iterations. 0.013461 Δα.
55. 1.85 with step -0.0020 (0.1996) converged in 1 iterations. 0.023035 Δα.
56. 1.85 with step -0.0022 (0.1996) converged in 1 iterations. 0.042910 Δα.
57. 1.85 with step -0.0024 (0.1996) converged in 1 iterations. 0.002973 Δα.
58. 1.84 with step -0.0026 (0.1996) converged in 1 iterations. 0.027035 Δα.
59. 1.84 with step -0.0027 (0.1996) converged in 1 iterations. 0.041029 Δα.
60. 1.84 with step -0.0027 (0.1996) converged in 1 iterations. 0.097557 Δα.
61. 1.83 with step -0.0030 (0.1996) converged in 1 iterations. 0.022513 Δα.
62. 1.83 with step -0.0031 (0.1996) converged in 1 iterations. 0.038264 Δα.
63. 1.83 with step -0.0032 (0.1996) converged in 1 iterations. 0.105127 Δα.
64. 1.82 with step -0.0036 (0.1996) converged in 1 iterations. 0.039754 Δα.
65. 1.82 with step -0.0037 (0.1996) converged in 1 iterations. 0.042279 Δα.
66. 1.82 with step -0.0040 (0.1996) converged in 1 iterations. 0.027536 Δα.
67. 1.81 with step -0.0042 (0.1996) converged in 1 iterations. 0.041250 Δα.
68. 1.81 with step -0.0046 (0.1996) converged in 1 iterations. 0.025314 Δα.
69. 1.80 with step -0.0047 (0.1996) converged in 1 iterations. 0.046109 Δα.
70. 1.80 with step -0.0050 (0.1996) converged in 1 iterations. 0.089411 Δα.
71. 1.79 with step -0.0054 (0.1996) converged in 2 iterations. 0.142186 Δα.
72. 1.79 with step -0.0067 (0.1996) converged in 1 iterations. 0.078338 Δα.
73. 1.78 with step -0.0081 (0.1996) converged in 2 iterations. 0.117076 Δα.
74. 1.77 with step -0.0101 (0.1996) converged in 3 iterations. 0.128194 Δα.
75. 1.75 with step -0.0128 (0.1996) converged in 2 iterations. 0.099065 Δα.
76. 1.74 with step -0.0155 (0.1996) converged in 2 iterations. 0.148778 Δα.
77. 1.72 with step -0.0193 (0.1996) converged in 2 iterations. 0.074061 Δα.
78. 1.69 with step -0.0216 (0.1996) converged in 2 iterations. 0.112074 Δα.
79. 1.67 with step -0.0240 (0.1996) converged in 2 iterations. 0.087439 Δα.
80. 1.65 with step -0.0250 (0.1996) converged in 2 iterations. 0.080135 Δα.
81. 1.62 with step -0.0254 (0.1996) converged in 2 iterations. 0.108032 Δα.
82. 1.60 with step -0.0245 (0.1996) converged in 2 iterations. 0.145266 Δα.
83. 1.58 with step -0.0225 (0.1996) converged in 2 iterations. 0.136929 Δα.
84. 1.56 with step -0.0199 (0.1996) converged in 2 iterations. 0.116508 Δα.
85. 1.54 with step -0.0167 (0.1996) converged in 2 iterations. 0.086887 Δα.
86. 1.53 with step -0.0135 (0.1996) converged in 2 iterations. 0.004448 Δα.
87. 1.52 with step -0.0099 (0.1996) converged in 2 iterations. 0.023424 Δα.
88. 1.51 with step -0.0088 (0.1996) converged in 2 iterations. 0.083395 Δα.
89. 1.51 with step -0.0077 (0.1996) converged in 2 iterations. 0.162303 Δα.
90. 1.50 with step -0.0067 (0.1996) converged in 2 iterations. 0.280972 Δα.
91. 1.49 with step -0.0058 (0.1996) converged in 2 iterations. 0.470452 Δα.
92. 1.49 with step -0.0051 (0.1996) converged in 2 iterations. 1.044048 Δα.
93. 1.48 with step -0.0044 (0.1996) converged in 2 iterations. 0.151697 Δα.
94. 1.48 with step -0.0038 (0.1996) converged in 2 iterations. 0.952525 Δα.
95. 1.48 with step -0.0033 (0.1996) converged in 2 iterations. 0.598313 Δα.
96. 1.48 with step -0.0028 (0.1996) converged in 2 iterations. 0.424231 Δα.
97. 1.47 with step -0.0026 (0.1996) converged in 2 iterations. 0.493307 Δα.
98. 1.47 with step -0.0024 (0.1996) converged in 2 iterations. 0.263467 Δα.
99. 1.47 with step -0.0023 (0.1996) converged in 2 iterations. 0.318345 Δα.
100. 1.47 with step -0.0020 (0.1996) converged in 2 iterations. 0.214118 Δα.
101. 1.46 with step -0.0020 (0.1996) converged in 2 iterations. 0.405864 Δα.
102. 1.46 with step -0.0017 (0.1996) converged in 2 iterations. 0.380054 Δα.
103. 1.46 with step -0.0012 (0.1996) converged in 2 iterations. 0.204723 Δα.
104. 1.46 with step -0.0014 (0.1996) converged in 2 iterations. 0.312877 Δα.
105. 1.46 with step -0.0011 (0.1996) converged in 2 iterations. 0.510007 Δα.
106. 1.46 with step -0.0012 (0.1996) converged in 2 iterations. 0.813096 Δα.
107. 1.46 with step -0.0018 (0.1996) converged in 2 iterations. 0.503413 Δα.
108. 1.45 with step -0.0024 (0.1996) converged in 2 iterations. 0.416089 Δα.
109. 1.45 with step -0.0030 (0.1996) converged in 2 iterations. 0.379177 Δα.
110. 1.45 with step -0.0036 (0.1996) converged in 2 iterations. 0.347959 Δα.
111. 1.44 with step -0.0042 (0.1996) converged in 2 iterations. 0.303667 Δα.
112. 1.44 with step -0.0045 (0.1996) converged in 2 iterations. 0.247640 Δα.
113. 1.43 with step -0.0044 (0.1996) converged in 2 iterations. 0.306971 Δα.
114. 1.43 with step -0.0043 (0.1996) converged in 1 iterations. 0.161446 Δα.
115. 1.42 with step -0.0041 (0.1996) converged in 2 iterations. 0.160979 Δα.
116. 1.42 with step -0.0046 (0.1996) converged in 2 iterations. 0.227466 Δα.
117. 1.38 with step -0.0056 (0.1996) converged in 5 iterations. 0.582534 Δα.
118. 1.35 with step -0.0300 (0.1838) converged in 2 iterations. 0.192633 Δα.
119. 1.32 with step -0.0307 (0.1996) converged in 2 iterations. 0.163625 Δα.
120. 1.29 with step -0.0271 (0.1996) converged in 2 iterations. 0.129304 Δα.
121. 1.27 with step -0.0225 (0.1996) converged in 2 iterations. 0.094986 Δα.
122. 1.26 with step -0.0180 (0.1996) converged in 2 iterations. 0.063740 Δα.
123. 1.24 with step -0.0140 (0.1996) converged in 2 iterations. 0.037061 Δα.
124. 1.23 with step -0.0108 (0.1996) converged in 2 iterations. 0.380982 Δα.
125. 1.22 with step -0.0084 (0.1996) converged in 2 iterations. 0.342173 Δα.
126. 1.22 with step -0.0067 (0.1996) converged in 2 iterations. 0.313021 Δα.
127. 1.21 with step -0.0056 (0.1996) converged in 1 iterations. 0.295920 Δα.
128. 1.21 with step -0.0047 (0.1996) converged in 1 iterations. 0.280490 Δα.
129. 1.21 with step -0.0040 (0.1996) converged in 1 iterations. 0.268027 Δα.
130. 1.20 with step -0.0035 (0.1996) converged in 1 iterations. 0.262847 Δα.
131. 1.20 with step -0.0032 (0.1996) converged in 1 iterations. 0.261282 Δα.
132. 1.20 with step -0.0029 (0.1996) converged in 1 iterations. 0.262146 Δα.
133. 1.19 with step -0.0026 (0.1996) converged in 1 iterations. 0.262192 Δα.
134. 1.19 with step -0.0024 (0.1996) converged in 1 iterations. 0.261009 Δα.
135. 1.19 with step -0.0022 (0.1996) converged in 2 iterations. 0.257698 Δα.
136. 1.19 with step -0.0020 (0.1996) converged in 2 iterations. 0.253153 Δα.
137. 1.19 with step -0.0018 (0.1996) converged in 2 iterations. 0.245231 Δα.
138. 1.18 with step -0.0016 (0.1996) converged in 2 iterations. 0.233052 Δα.
139. 1.18 with step -0.0015 (0.1996) converged in 2 iterations. 0.347640 Δα.
140. 1.18 with step -0.0014 (0.1996) converged in 2 iterations. 0.334354 Δα.
141. 1.18 with step -0.0012 (0.1996) converged in 2 iterations. 0.312927 Δα.
142. 1.18 with step -0.0011 (0.1996) converged in 2 iterations. 0.282482 Δα.
143. 1.18 with step -0.0009 (0.1996) converged in 2 iterations. 0.248128 Δα.
144. 1.18 with step -0.0007 (0.1996) converged in 2 iterations. 0.203437 Δα.
145. 1.18 with step -0.0006 (0.1996) converged in 2 iterations. 0.146438 Δα.
146. 1.18 with step -0.0005 (0.1996) converged in 2 iterations. 0.102133 Δα.
147. 1.18 with step -0.0004 (0.1996) converged in 2 iterations. 0.070232 Δα.
148. 1.18 with step -0.0003 (0.1996) converged in 2 iterations. 0.048045 Δα.
149. 1.18 with step -0.0002 (0.1996) converged in 1 iterations. 0.032795 Δα.
150. 1.18 with step -0.0002 (0.1996) converged in 1 iterations. 0.022456 Δα.
151. 1.18 with step -0.0002 (0.1996) converged in 1 iterations. 0.015351 Δα.
152. 1.18 with step -0.0001 (0.1996) converged in 1 iterations. 0.010685 Δα.
153. 1.18 with step -0.0001 (0.1996) converged in 1 iterations. 0.007790 Δα.
154. 1.18 with step -0.0001 (0.1996) converged in 1 iterations. 0.005462 Δα.
155. 1.18 with step -0.0001 (0.1996) converged in 1 iterations. 0.004250 Δα.
156. 1.18 with step -0.0001 (0.1996) converged in 1 iterations. 0.003045 Δα.
157. 1.18 with step -0.0001 (0.1996) converged in 1 iterations. 0.001845 Δα.
158. 1.18 with step -0.0001 (0.1996) converged in 1 iterations. 0.000657 Δα.
159. 1.18 with step -0.0001 (0.1996) converged in 1 iterations. 0.000586 Δα.
160. 1.18 with step -0.0001 (0.1996) converged in 1 iterations. 0.001772 Δα.
161. 1.18 with step -0.0001 (0.1996) converged in 1 iterations. 0.002974 Δα.
162. 1.18 with step -0.0001 (0.1996) converged in 1 iterations. 0.004181 Δα.
163. 1.18 with step -0.0001 (0.1996) converged in 1 iterations. 0.005396 Δα.
164. 1.17 with step -0.0001 (0.1996) converged in 1 iterations. 0.006618 Δα.
165. 1.17 with step -0.0001 (0.1996) converged in 1 iterations. 0.008021 Δα.
166. 1.17 with step -0.0001 (0.1996) converged in 1 iterations. 0.009952 Δα.
167. 1.17 with step -0.0001 (0.1996) converged in 1 iterations. 0.012374 Δα.
168. 1.17 with step -0.0001 (0.1996) converged in 1 iterations. 0.015255 Δα.
169. 1.17 with step -0.0001 (0.1996) converged in 1 iterations. 0.018999 Δα.
170. 1.17 with step -0.0002 (0.1996) converged in 1 iterations. 0.024024 Δα.
171. 1.17 with step -0.0002 (0.1996) converged in 1 iterations. 0.030665 Δα.
172. 1.17 with step -0.0002 (0.1996) converged in 1 iterations. 0.039408 Δα.
173. 1.17 with step -0.0003 (0.1996) converged in 1 iterations. 0.050950 Δα.
174. 1.17 with step -0.0003 (0.1996) converged in 2 iterations. 0.066375 Δα.
175. 1.17 with step -0.0004 (0.1996) converged in 2 iterations. 0.086897 Δα.
176. 1.17 with step -0.0004 (0.1996) converged in 2 iterations. 0.114688 Δα.
177. 1.17 with step -0.0005 (0.1996) converged in 2 iterations. 0.175311 Δα.
178. 1.17 with step -0.0006 (0.1996) converged in 2 iterations. 0.205700 Δα.
179. 1.17 with step -0.0008 (0.1996) converged in 2 iterations. 0.251774 Δα.
180. 1.17 with step -0.0010 (0.1996) converged in 2 iterations. 0.311653 Δα.
181. 1.17 with step -0.0013 (0.1996) converged in 2 iterations. 0.375100 Δα.
182. 1.17 with step -0.0015 (0.1996) converged in 2 iterations. 0.427289 Δα.
183. 1.16 with step -0.0018 (0.1996) converged in 2 iterations. 0.302351 Δα.
184. 1.16 with step -0.0021 (0.1996) converged in 2 iterations. 0.328154 Δα.
185. 1.16 with step -0.0024 (0.1996) converged in 2 iterations. 0.339438 Δα.
186. 1.16 with step -0.0028 (0.1996) converged in 2 iterations. 0.340478 Δα.
187. 1.15 with step -0.0031 (0.1996) converged in 2 iterations. 0.334512 Δα.
188. 1.15 with step -0.0035 (0.1996) converged in 2 iterations. 0.323880 Δα.
189. 1.14 with step -0.0038 (0.1996) converged in 2 iterations. 0.310315 Δα.
190. 1.14 with step -0.0042 (0.1996) converged in 2 iterations. 0.306161 Δα.
191. 1.13 with step -0.0046 (0.1996) converged in 1 iterations. 0.296338 Δα.
192. 1.13 with step -0.0052 (0.1996) converged in 1 iterations. 0.298366 Δα.
193. 1.12 with step -0.0060 (0.1996) converged in 1 iterations. 0.302518 Δα.
194. 1.12 with step -0.0070 (0.1996) converged in 1 iterations. 0.306823 Δα.
195. 1.11 with step -0.0082 (0.1996) converged in 1 iterations. 0.315994 Δα.
196. 1.10 with step -0.0099 (0.1996) converged in 2 iterations. 0.329508 Δα.
197. 1.08 with step -0.0121 (0.1996) converged in 2 iterations. 0.344969 Δα.
198. 1.07 with step -0.0146 (0.1996) converged in 2 iterations. 0.037209 Δα.
199. 1.05 with step -0.0177 (0.1996) converged in 2 iterations. 0.051065 Δα.
200. 1.03 with step -0.0215 (0.1996) converged in 2 iterations. 0.065191 Δα.
201. 1.00 with step -0.0261 (0.1996) converged in 2 iterations. 0.080138 Δα.
202. 0.97 with step -0.0317 (0.1996) converged in 2 iterations. 0.096537 Δα.
203. 0.93 with step -0.0386 (0.1996) converged in 2 iterations. 0.115052 Δα.
204. 0.88 with step -0.0470 (0.1996) converged in 2 iterations. 0.136395 Δα.
205. 0.82 with step -0.0573 (0.1996) converged in 2 iterations. 0.161341 Δα.
206. 0.74 with step -0.0694 (0.1996) converged in 2 iterations. 0.190578 Δα.
207. 0.66 with step -0.0830 (0.1996) converged in 2 iterations. 0.224090 Δα.
208. 0.56 with step -0.0954 (0.1996) converged in 2 iterations. 0.259269 Δα.
209. 0.45 with step -0.1016 (0.1996) converged in 2 iterations. 0.289244 Δα.
210. 0.35 with step -0.0971 (0.1996) converged in 2 iterations. 0.308357 Δα.
211. 0.27 with step -0.0837 (0.1996) converged in 2 iterations. 0.318260 Δα.
212. 0.20 with step -0.0673 (0.1996) converged in 1 iterations. 0.323792 Δα.
213. 0.15 with step -0.0520 (0.1996) converged in 1 iterations. 0.327549 Δα.
214. 0.11 with step -0.0393 (0.1996) converged in 1 iterations. 0.330553 Δα.
215. 0.08 with step -0.0293 (0.1996) converged in 1 iterations. 0.332990 Δα.Plotting
We plot first just the responses, to demonstrate how the amplitude constraint has been effective.
Amplitude Responses
his = [1, 3, 5];
set_theme!(theme_latexfonts())
fsz = 24;
fig1 = Figure(fontsize=fsz, size=(1000, 600));
for i in eachindex(his[his.<=maximum(h)])
ax = Axis(fig1[1, i],
ylabel=L"$H_%$(his[i])$ Response (m)", yscale=log10);
scatterlines!(ax, Oms, abs.(uh[his[i].+1, :, 1]), label="x1");
scatterlines!(ax, Oms, abs.(uh[his[i].+1, :, 2]), label="x2");
if i==1
axislegend(ax)
end
ax = Axis(fig1[2, i], xlabel=L"Excitation Frequency $\Omega$",
ylabel=L"$H_%$(his[i])$ Phase (rad)");
scatterlines!(ax, Oms, unwrap(angle.(uh[his[i].+1, :, 1])), label="x1");
scatterlines!(ax, Oms, unwrap(angle.(uh[his[i].+1, :, 2])), label="x2");
end
fig1
It can be observed that the first harmonic of the first DoF is always fixed to the value of Amp ($1.0$ here), while the other DoF's amplitude changes, showing that we have a fixed amplitude forced response here.
Forced Responses
In order to visualize the nonlinear forced response function, we plot the response over force next. Here the resonances are clearly visible, albeit very different from what was seen in Example C2.
his = [1, 3, 5];
set_theme!(theme_latexfonts())
fsz = 24;
fig2 = Figure(fontsize=fsz, size=(1000, 600));
for i in eachindex(his[his.<=maximum(h)])
ax = Axis(fig2[1, i],
ylabel=L"$H_%$(his[i])$ Response (m/N)", yscale=log10);
scatterlines!(ax, Oms, abs.(uh[his[i].+1, :, 1]./Fs), label="x1");
scatterlines!(ax, Oms, abs.(uh[his[i].+1, :, 2]./Fs), label="x2");
ax = Axis(fig2[2, i], xlabel=L"Excitation Frequency $\Omega$",
ylabel=L"$H_%$(his[i])$ Phase (rad)");
scatterlines!(ax, Oms, unwrap(angle.(uh[his[i].+1, :, 1])), label="x1");
scatterlines!(ax, Oms, unwrap(angle.(uh[his[i].+1, :, 2])), label="x2");
end
fig2
This page was generated using Literate.jl.