Automatic Control Knowledge Repository

You currently have javascript disabled. Some features will be unavailable. Please consider enabling javascript.

Details for: "triple pendulum"

Name: triple pendulum (Key: HBWN5)
Path: ackrep_data/system_models/triple_pendulum_system View on GitHub
Type: system_model
Short Description: triple pendulum on a cart, chaotic system
Created: 2022-11-11
Compatible Environment: default_conda_environment (Key: CDAMA)
Source Code [ / ] simulation.py
import numpy as np
import system_model
from scipy.integrate import odeint

from ackrep_core import ResultContainer
from ackrep_core.system_model_management import save_plot_in_dir
import matplotlib.pyplot as plt
import os

# link to documentation with examples: https://ackrep-doc.readthedocs.io/en/latest/devdoc/contributing_data.html


def simulate():
    """
    simulate the system model with scipy.integrate.solve_ivp

    :return: result of solve_ivp, might contains input function
    """

    model = system_model.Model()

    eqns, ff, xx, gg, _ = model.get_symbolic_model()
    print("Lagrange Equations:\n")
    for i, eq in enumerate(eqns):
        print(f"0 =", eq)

    rhs = model.get_rhs_odeint_fnc(ff, xx, gg)

    # initial state values
    xx0 = np.array([2, 2, 2, 0, 0, 0, 0, 0])

    t_end = 15
    tt = np.linspace(0, t_end, 10000)
    simulation_data = odeint(rhs, xx0, tt)

    save_plot(simulation_data, tt)

    return simulation_data


def save_plot(simulation_data, tt):
    """
    plot your data and save the plot
    access to data via: simulation_data.t   array of time values
                        simulation_data.y   array of data components
                        simulation_data.uu  array of input values

    :param simulation_data: simulation_data of system_model
    :return: None
    """
    # ---------start of edit section--------------------------------------
    # plot of your data
    for i in range(len(simulation_data[0, :])):
        plt.plot(tt, simulation_data[:, i], label="$x_{}$".format(i + 1))
    plt.grid()
    plt.legend()
    plt.xlabel("Time [s]")

    # ---------end of edit section----------------------------------------

    plt.tight_layout()

    save_plot_in_dir()


def evaluate_simulation(simulation_data):
    """
    assert that the simulation results are as expected

    :param simulation_data: simulation_data of system_model
    :return:
    """
    # ---------start of edit section--------------------------------------
    # fill in final states of simulation to check your model
    # simulation_data.y[i][-1]
    expected_final_state = [
        3.3544161195899775,
        3.6128436581303522,
        4.184663788192004,
        0.0,
        -2.0778538809498333,
        -5.473145553078993,
        -2.3347177255025864,
        0.0,
    ]

    # ---------end of edit section----------------------------------------

    rc = ResultContainer(score=1.0)
    simulated_final_state = []
    for i in range(simulation_data.shape[1]):
        simulated_final_state.append(simulation_data[:, i][-1])

    rc.final_state_errors = [
        simulated_final_state[i] - expected_final_state[i] for i in np.arange(0, len(simulated_final_state))
    ]
    rc.success = np.allclose(expected_final_state, simulated_final_state, rtol=0, atol=1e-2)

    return rc
system_model.py
import sympy as sp
import symbtools as st
import importlib
import sys, os
import pickle
from symbtools import modeltools as mt

# from ipydex import IPS, activate_ips_on_exception

from ackrep_core.system_model_management import GenericModel, import_parameters

# Import parameter_file
params = import_parameters()


# link to documentation with examples: https://ackrep-doc.readthedocs.io/en/latest/devdoc/contributing_data.html


class Model(GenericModel):
    def initialize(self):
        """
        this function is called by the constructor of GenericModel

        :return: None
        """

        # ---------start of edit section--------------------------------------
        # Define number of inputs -- MODEL DEPENDENT
        self.u_dim = 1

        # Set "sys_dim" to constant value, if system dimension is constant
        self.sys_dim = 8

        # ---------end of edit section----------------------------------------

        # check existence of params file
        self.has_params = True
        self.params = params

    # ----------- SET DEFAULT INPUT FUNCTION ---------- #
    # --------------- Only for non-autonomous Systems
    def uu_default_func(self):
        """
        define input function

        :return:(function with 2 args - t, xx_nv) default input function
        """

        # ---------start of edit section--------------------------------------
        def uu_rhs(t, xx_nv):
            """
            sequence of numerical input values

            :param t:(scalar or vector) time
            :param xx_nv:(vector or array of vectors) numeric state vector
            :return:(list) numeric inputs
            """
            u1 = 1

            return [u1]

        # ---------end of edit section----------------------------------------

        return uu_rhs

    # ----------- SYMBOLIC RHS FUNCTION ---------- #

    def get_symbolic_model(self):
        """
        define model

        :return: model
        """

        x1, x2, x3, x4, x5, x6, x7, x8 = self.xx_symb  # state components
        m0, m1, m2, m3, J1, J2, J3, l1, l2, l3, a1, a2, a3, g = self.pp_symb  # parameters

        u1 = self.uu_symb[0]  # inputs

        ttheta = sp.Matrix([[x1], [x2], [x3], [x4]])
        xdot1, xdot2, xdot3, xdot4 = sp.symbols("xdot1, xdot2, xdot3, xdot4")

        ex = sp.Matrix([1, 0])
        ey = sp.Matrix([0, 1])

        Rz = mt.Rz

        S0 = G0 = ex * x4

        G1 = G0 + Rz(x1) * ey * l1
        S1 = G0 + Rz(x1) * ey * a1

        G2 = G1 + Rz(x2) * ey * l2
        S2 = G1 + Rz(x2) * ey * a2

        G3 = G2 + Rz(x3) * ey * l3
        S3 = G2 + Rz(x3) * ey * a2

        # Time derivatives of coordinates of the centers of masses
        S0dt, S1dt, S2dt, S3dt = st.col_split(st.time_deriv(st.col_stack(S0, S1, S2, S3), ttheta))

        # kinetic energy of the cart
        T0 = 0.5 * m0 * xdot4**2
        # kinetic energy of pendulum1
        T1 = 0.5 * m1 * (S1dt.T * S1dt)[0] + 0.5 * J1 * xdot1**2
        # kinetic energy of pendulum2
        T2 = 0.5 * m2 * (S2dt.T * S2dt)[0] + 0.5 * J2 * xdot2**2
        # kinetic energy of pendulum3
        T3 = 0.5 * m3 * (S3dt.T * S3dt)[0] + 0.5 * J3 * xdot3**2

        # total kinetic energy
        T = T0 + T1 + T2 + T3

        # total potential energy
        V = g * (m1 * S1[1] + m2 * S2[1] + m3 * S3[1])

        external_forces = [0, 0, 0, u1]

        dir_of_this_file = os.path.dirname(os.path.abspath(sys.modules.get(__name__).__file__))
        fpath = os.path.join(dir_of_this_file, "pendulum.pcl")

        if not os.path.isfile(fpath):
            # Calculate the model based on lagrange equation (about 30s)
            mod = mt.generate_symbolic_model(T, V, ttheta, external_forces)

            # perform patial linearization such that system input is acceleration and not force (about 9min)
            mod.calc_coll_part_lin_state_eq()

            # write the model to disk to save time in the next run of the notebook
            with open(fpath, "wb") as pfile:
                pickle.dump(mod, pfile)
        else:
            with open(fpath, "rb") as pfile:
                mod = pickle.load(pfile)

        eqns = mod.eqns.subs([(xdot1, x5), (xdot2, x6), (xdot3, x7), (xdot4, x8)])

        ff = mod.ff.subs([(xdot1, x5), (xdot2, x6), (xdot3, x7), (xdot4, x8)])
        xx = mod.xx.subs([(xdot1, x5), (xdot2, x6), (xdot3, x7), (xdot4, x8)])
        gg = mod.gg.subs([(xdot1, x5), (xdot2, x6), (xdot3, x7), (xdot4, x8)])

        return [eqns, ff, xx, gg, mod]

    def get_rhs_odeint_fnc(self, ff, xx, gg):
        """
        Creates an executable function of the model for the odeint solver

        :return: rhs function
        """
        parameter_values = list(self.pp_dict.items())

        ff = ff.subs(parameter_values)
        xx = xx.subs(parameter_values)
        gg = gg.subs(parameter_values)
        simmod = st.SimulationModel(ff, gg, xx)
        rhs = simmod.create_simfunction()

        return rhs

    def get_rhs_symbolic(self):
        """
        define symbolic rhs function

        :return: matrix of symbolic rhs-functions
        """

        _, _, _, _, mod = self.get_symbolic_model()
        mod.calc_state_eq(simplify=False)
        x5, x6, x7, x8, xdot1, xdot2, xdot3, xdot4 = sp.symbols("x5, x6, x7, x8, xdot1, xdot2, xdot3, xdot4")
        state_eq = mod.state_eq.subs([(xdot1, x5), (xdot2, x6), (xdot3, x7), (xdot4, x8)])

        return state_eq
parameters.py
import sys
import os
import numpy as np
import sympy as sp

import tabulate as tab


# link to documentation with examples: https://ackrep-doc.readthedocs.io/en/latest/devdoc/contributing_data.html


# set model name
model_name = "triple pendulum"


# ---------- create symbolic parameters
pp_symb = [m0, m1, m2, m3, J1, J2, J3, l1, l2, l3, a1, a2, a3, g] = sp.symbols(
    "m0, m1, m2, m3, J1, J2, J3, l1, l2, l3, a1, a2, a3, g", real=True
)


# ---------- create symbolic parameter functions
# parameter values can be constant/fixed values OR set in relation to other parameters (for example: a = 2*b)
m0_sf = 3.34
m1_sf = 0.8512
m2_sf = 0.8973
m3_sf = 0.5519
J1_sf = 0.01980194
J2_sf = 0.02105375
J3_sf = 0.01818537
l1_sf = 0.32
l2_sf = 0.419
l3_sf = 0.485
a1_sf = 0.20001517
a2_sf = 0.26890449
a3_sf = 0.21666087
g_sf = 9.81

# list of symbolic parameter functions
# tailing "_sf" stands for "symbolic parameter function"
pp_sf = [m0_sf, m1_sf, m2_sf, m3_sf, J1_sf, J2_sf, J3_sf, l1_sf, l2_sf, l3_sf, a1_sf, a2_sf, a3_sf, g_sf]


#  ---------- list for substitution
# -- entries are tuples like: (independent symbolic parameter, numerical value)
pp_subs_list = []


# OPTONAL: Dictionary which defines how certain variables shall be written
# in the table - key: Symbolic Variable, Value: LaTeX Representation/Code
# useful for example for complex variables: {Z: r"\underline{Z}"}
latex_names = {}


# ---------- Define LaTeX table

# Define table header
# DON'T CHANGE FOLLOWING ENTRIES: "Symbol", "Value"
tabular_header = ["Parameter Name", "Symbol", "Value", "Unit"]

# Define column text alignments
col_alignment = ["left", "center", "left", "center"]


# Define Entries of all columns before the Symbol-Column
# --- Entries need to be latex code
col_1 = [
    "mass of the cart",
    "mass of link 1",
    "mass of link 2",
    "mass of link 3",
    "moment of inertia of link 1",
    "moment of inertia of link 2",
    "moment of inertia of link 3",
    "length of link 1",
    "length of link 2",
    "length of link 3",
    "distance from the joint to the center of gravity of link 1",
    "distance from the joint to the center of gravity of link 2",
    "distance from the joint to the center of gravity of link 3",
    "acceleration due to gravity",
]

# contains all lists of the columns before the "Symbol" Column
# --- Empty list, if there are no columns before the "Symbol" Column
start_columns_list = [col_1]


# Define Entries of the columns after the Value-Column
# --- Entries need to be latex code
col_4 = [
    "kg",
    "kg",
    "kg",
    "kg",
    r"$kg \cdot m^2$",
    r"$kg \cdot m^2$",
    r"$kg \cdot m^2$",
    "m",
    "m",
    "m",
    "m",
    "m",
    "m",
    r"$\frac{m}{s^2}$",
]

# contains all lists of columns after the FIX ENTRIES
# --- Empty list, if there are no columns after the "Value" column
end_columns_list = [col_4]

Related Problems:
Extensive Material:
Download pdf
Result: Success.
Last Build: Checkout CI Build
Runtime: 5.5 (estimated: 10s)