Skip to content

Commit

Permalink
Merge pull request #6378 from iago-mendes/damp-control
Browse files Browse the repository at this point in the history
Damp initial iterations of Broyden's method used in ID parameter control
  • Loading branch information
nilsvu authored Nov 22, 2024
2 parents b4232b3 + 958ca95 commit eeafc0d
Showing 1 changed file with 26 additions and 0 deletions.
26 changes: 26 additions & 0 deletions support/Pipelines/Bbh/ControlId.py
Original file line number Diff line number Diff line change
Expand Up @@ -288,14 +288,40 @@ def Residual(u):
# Initialize Jacobian as an identity matrix
J = np.identity(len(u))

# Indices of parameters for which the control is delayed in the first
# iterations to avoid going off-bounds
#
# Note: We have experimented with other modifications to Broyden's
# method, including damping the initial updates of the free data / Jacobian
# and enforcing a diagonal Jacobian. None of them converged as fast as the
# delay approach used here. When doing a more complete study in parameter
# space, we should try to find a more robust approach that works for
# multiple configurations.
delayed_indices = np.array([], dtype=bool)
delayed_params = ["center_of_mass", "linear_momentum"]
max_delayed_iteration = 1
for key in control_params.keys():
if key in ScalarQuantities:
delayed_indices = np.append(
delayed_indices, [key in delayed_params]
)
else:
delayed_indices = np.append(
delayed_indices, [key in delayed_params] * 3
)

while iteration < max_iterations and np.max(np.abs(F)) > residual_tolerance:
iteration += 1

# Update the free parameters using a quasi-Newton-Raphson method
Delta_u = -np.dot(np.linalg.inv(J), F)
if iteration <= max_delayed_iteration:
Delta_u[delayed_indices] = 0.0
u += Delta_u

F = Residual(u)
if iteration <= max_delayed_iteration:
F[delayed_indices] = 0.0

# Update the Jacobian using Broyden's method
J += np.outer(F, Delta_u) / np.dot(Delta_u, Delta_u)
Expand Down

0 comments on commit eeafc0d

Please sign in to comment.