From b221eceb7d7cadbfef37c624b3b05376b873c9e0 Mon Sep 17 00:00:00 2001 From: fpacaud Date: Fri, 19 Apr 2024 19:58:23 +0200 Subject: [PATCH] add missing COPS instances --- src/PureJuMP/catmix.jl | 76 +++++++++++++++++++++ src/PureJuMP/gasoil.jl | 117 ++++++++++++++++++++++++++++++++ src/PureJuMP/glider.jl | 74 +++++++++++++++++++++ src/PureJuMP/methanol.jl | 140 +++++++++++++++++++++++++++++++++++++++ src/PureJuMP/minsurf.jl | 57 ++++++++++++++++ src/PureJuMP/pinene.jl | 106 +++++++++++++++++++++++++++++ src/PureJuMP/rocket.jl | 58 ++++++++++++++++ src/PureJuMP/steering.jl | 48 ++++++++++++++ src/PureJuMP/torsion.jl | 47 +++++++++++++ 9 files changed, 723 insertions(+) create mode 100644 src/PureJuMP/catmix.jl create mode 100644 src/PureJuMP/gasoil.jl create mode 100644 src/PureJuMP/glider.jl create mode 100644 src/PureJuMP/methanol.jl create mode 100644 src/PureJuMP/minsurf.jl create mode 100644 src/PureJuMP/pinene.jl create mode 100644 src/PureJuMP/rocket.jl create mode 100644 src/PureJuMP/steering.jl create mode 100644 src/PureJuMP/torsion.jl diff --git a/src/PureJuMP/catmix.jl b/src/PureJuMP/catmix.jl new file mode 100644 index 00000000..44f4f986 --- /dev/null +++ b/src/PureJuMP/catmix.jl @@ -0,0 +1,76 @@ +# Catalyst Mixing Problem +# Collocation formulation +# COPS 3.0 - November 2002 +# COPS 3.1 - March 2004 + +function catmix(args...; n::Int = default_nvar, kwargs...) + ne = 2 + nc = 3 + + tf = 1 + h = tf / n # Final time + fact = [factorial(k) for k in 0:nc] + + rho = [ + 0.11270166537926, + 0.50000000000000, + 0.88729833462074, + ] + bc = [1.0, 0.0] # Boundary conditions for x + alpha = 0.0 # Smoothing parameter + + model = Model() + + @variable(model, 0.0 <= u[i=1:n, j=1:nc] <= 1.0, start=0.0) + @variable(model, v[i=1:n, s=1:ne], start=mod(s, ne)) + @variable(model, w[i=1:n, j=1:nc, s=1:ne], start=0.0) + @variable(model, pp[i=1:n, j=1:nc, s=1:ne], start=mod(s, ne)) + @variable(model, Dpp[i=1:n, j=1:nc, s=1:ne], start=0.0) + @variable(model, ppf[s=1:ne], start=mod(s, ne)) + + @objective( + model, + Min, + -1.0 + ppf[1] + ppf[2] + + alpha/h*sum((u[i+1, j] - u[i, j])^2 for i in 1:n-1, j in 1:nc) + ) + + # Collocation model + @constraint( + model, + [i=1:n, k=1:nc, s=1:ne], + pp[i, k, s] == v[i, s] + h*sum(w[i, j, s]*(rho[k]^j/fact[j+1]) for j in 1:nc) + ) + @constraint( + model, + [i=1:n, k=1:nc, s=1:ne], + Dpp[i, k, s] == sum(w[i, j, s]*(rho[k]^(j-1)/fact[j]) for j in 1:nc) + ) + @constraint( + model, + [s=1:ne], + ppf[s] == v[n, s] + h * sum(w[n, j, s] / fact[j+1] for j in 1:nc) + ) + # Continuity + @constraint( + model, + continuity[i=1:n-1, s=1:ne], + v[i, s] + sum(w[i, j, s] * h / fact[j+1] for j in 1:nc) == v[i+1, s] + ) + # Dynamics + @NLconstraint( + model, + de1[i=1:n, j=1:nc], + Dpp[i,j,1] == u[i,j] * (10.0*pp[i,j,2] - pp[i,j,1]), + ) + @NLconstraint( + model, + de2[i=1:n, j=1:nc], + Dpp[i,j,2] == u[i,j] * (pp[i,j,1] - 10.0*pp[i,j,2]) - (1 - u[i,j])*pp[i,j,2] + ) + @constraint(model, b_eqn[s=1:ne], v[1, s] == bc[s]) + + return model +end + + diff --git a/src/PureJuMP/gasoil.jl b/src/PureJuMP/gasoil.jl new file mode 100644 index 00000000..bdba5de5 --- /dev/null +++ b/src/PureJuMP/gasoil.jl @@ -0,0 +1,117 @@ +# Catalytic Cracking of Gas Oil Problem +# Collocation formulation +# Michael Merritt - Summer 2000 +# COPS 2.0 - September 2000 +# COPS 3.0 - November 2002 +# COPS 3.1 - March 2004 + +function gasoil(; n::Int = default_nvar, kwargs...) + nc = 4 # number of collocation points + ne = 2 # number of differential equations + np = 3 # number of ODE parameters + nm = 21 # number of measurements + + # roots of k-th degree Legendre polynomial + rho = [0.06943184420297, 0.33000947820757, 0.66999052179243, 0.93056815579703] + # ODE initial conditions + bc = [1, 1, 2, 0] + # times at which observations made + tau = [0.0, 0.025, 0.05, 0.075, 0.10, 0.125, 0.150, 0.175, 0.20, 0.225, 0.250, 0.30, 0.35, 0.40, 0.45, 0.50, 0.55, 0.65, 0.75, 0.85, 0.95] + # ODEs defined in [0,tf] + tf = tau[nm] + # uniform interval length + h = tf / n + t = [(i-1)*h for i in 1:n+1] + fact = [factorial(k) for k in 0:nc] + + itau = Int[min(n, floor(tau[i]/h)+1) for i in 1:nm] + + # Concentrations + z = [ + 1.0000 0; + 0.8105 0.2000; + 0.6208 0.2886; + 0.5258 0.3010; + 0.4345 0.3215; + 0.3903 0.3123; + 0.3342 0.2716; + 0.3034 0.2551; + 0.2735 0.2258; + 0.2405 0.1959; + 0.2283 0.1789; + 0.2071 0.1457; + 0.1669 0.1198; + 0.1530 0.0909; + 0.1339 0.0719; + 0.1265 0.0561; + 0.1200 0.0460; + 0.0990 0.0280; + 0.0870 0.0190; + 0.0770 0.0140; + 0.0690 0.0100; + ] + + v0 = zeros(n, ne) + # Starting-value + for i in 1:itau[1], s in 1:ne + v0[i, s] = bc[s] + end + for j in 2:nm, i =itau[j-1]+1:itau[j], s in 1:ne + v0[i, s] = z[j, s] + end + for i in itau[nm]+1:n, s in 1:ne + v0[i, s] = z[nm, s] + end + + model = Model() + + # ODE parameters + @variable(model, theta[1:np] >= 0.0, start=0.0) + # The collocation approximation u is defined by the parameters v and w. + # uc and Duc are, respectively, u and u' evaluated at the collocation points. + @variable(model, v[i=1:n, s=1:ne], start=v0[i, s]) + @variable(model, w[1:n, 1:nc, 1:ne], start=0.0) + @variable(model, uc[i=1:n, j=1:nc, s=1:ne], start=v0[i,s]) + @variable(model, Duc[i=1:n, j=1:nc, s=1:ne], start=0.0) + + @NLexpression( + model, + error[j=1:nm, s=1:ne], + v[itau[j],s] + sum(w[itau[j],k,s]*(tau[j]-t[itau[j]])^k/(fact[k+1]*h^(k-1)) for k in 1:nc) - z[j,s], + ) + # L2 error + @NLobjective(model, Min, sum(error[j, s]^2 for s in 1:ne, j in 1:nm)) + + # Collocation model + @constraint( + model, + [i=1:n, j=1:nc, s=1:ne], + uc[i, j, s] == v[i,s] + h*sum(w[i,k,s]*(rho[j]^k/fact[k+1]) for k in 1:nc), + ) + @constraint( + model, + [i=1:n, j=1:nc, s=1:ne], + Duc[i, j, s] == sum(w[i,k,s]*(rho[j]^(k-1)/fact[k]) for k in 1:nc), + ) + + # Boundary + @constraint(model, [s=1:ne], v[1, s] == z[1, s]) #TODO + # Continuity + @constraint( + model, + [i=1:n-1, s=1:ne], + v[i, s] + sum(w[i, j, s]*h/fact[j+1] for j in 1:nc) == v[i+1, s], + ) + @NLconstraint( + model, + [i=1:n, j=1:nc], + Duc[i, j, 1] == -(theta[1]+theta[3])*uc[i, j, 1]^2, + ) + @NLconstraint( + model, + [i=1:n, j=1:nc], + Duc[i, j, 2] == theta[1]*uc[i,j,1]^2 - theta[2]*uc[i,j,2], + ) + return model +end + diff --git a/src/PureJuMP/glider.jl b/src/PureJuMP/glider.jl new file mode 100644 index 00000000..f4ded3c9 --- /dev/null +++ b/src/PureJuMP/glider.jl @@ -0,0 +1,74 @@ +# Hang Glider Problem +# Trapezoidal formulation +# David Bortz - Summer 1998 +# COPS 2.0 - September 2000 +# COPS 3.0 - November 2002 +# COPS 3.1 - March 2004 + +function glider(; n::Int = default_nvar, kwargs...) + # Design parameters + x_0 = 0.0 + y_0 = 1000.0 + y_f = 900.0 + vx_0 = 13.23 + vx_f = 13.23 + vy_0 = -1.288 + vy_f = -1.288 + u_c = 2.5 + r_0 = 100.0 + m = 100.0 + g = 9.81 + c0 = 0.034 + c1 = 0.069662 + S = 14.0 + rho = 1.13 + cL_min = 0.0 + cL_max = 1.4 + + model = Model() + + @variables(model, begin + 0 <= t_f, (start=1.0) + 0.0 <= x[k=0:n], (start=x_0 + vx_0*(k/n)) + y[k=0:n], (start=y_0 + (k/n)*(y_f - y_0)) + 0.0 <= vx[k=0:n], (start=vx_0) + vy[k=0:n], (start=vy_0) + cL_min <= cL[k=0:n] <= cL_max, (start=cL_max/2.0) + end) + + @objective(model, Max, x[n]) + + @NLexpressions(model, begin + step, t_f / n + r[i=0:n], (x[i]/r_0 - 2.5)^2 + u[i=0:n], u_c*(1 - r[i])*exp(-r[i]) + w[i=0:n], vy[i] - u[i] + v[i=0:n], sqrt(vx[i]^2 + w[i]^2) + D[i=0:n], 0.5*(c0+c1*cL[i]^2)*rho*S*v[i]^2 + L[i=0:n], 0.5*cL[i]*rho*S*v[i]^2 + vx_dot[i=0:n], (-L[i]*(w[i]/v[i]) - D[i]*(vx[i]/v[i]))/m + vy_dot[i=0:n], (L[i]*(vx[i]/v[i]) - D[i]*(w[i]/v[i]))/m - g + end) + + # Dynamics + @NLconstraints(model, begin + x_eqn[j=1:n], x[j] == x[j-1] + 0.5 * step * (vx[j] + vx[j-1]) + y_eqn[j=1:n], y[j] == y[j-1] + 0.5 * step * (vy[j] + vy[j-1]) + vx_eqn[j=1:n], vx[j] == vx[j-1] + 0.5 * step * (vx_dot[j] + vx_dot[j-1]) + vy_eqn[j=1:n], vy[j] == vy[j-1] + 0.5 * step * (vy_dot[j] + vy_dot[j-1]) + end) + # Boundary constraints + @constraints(model, begin + x[0] == x_0 + y[0] == y_0 + y[n] == y_f + vx[0] == vx_0 + vx[n] == vx_f + vy[0] == vy_0 + vy[n] == vy_f + end) + + return model +end + + diff --git a/src/PureJuMP/methanol.jl b/src/PureJuMP/methanol.jl new file mode 100644 index 00000000..87524516 --- /dev/null +++ b/src/PureJuMP/methanol.jl @@ -0,0 +1,140 @@ +# Methanol-to-Hydrocarbons Problem +# Collocation formulation +# Michael Merritt - Summer 2000 +# COPS 2.0 - September 2000 +# COPS 3.0 - November 2002 +# COPS 3.1 - March 2004 + +function methanol(args...; n::Int = default_nvar, kwargs...) + ne = 3 + np = 5 + nc = 3 + nm = 17 + + rho = [0.11270166537926, 0.5, 0.88729833462074] + # times at which observations made + tau = [ + 0, + 0.050, + 0.065, + 0.080, + 0.123, + 0.233, + 0.273, + 0.354, + 0.397, + 0.418, + 0.502, + 0.553, + 0.681, + 0.750, + 0.916, + 0.937, + 1.122, + ] + tf = tau[nm] # ODEs defined in [0,tf] + h = tf / n # uniform interval length + t = [(i-1)*h for i in 1:n+1] # partition + fact = [factorial(k) for k in 0:nc] + + # itau[i] is the largest integer k with t[k] <= tau[i] + itau = Int[min(n, floor(tau[i]/h)+1) for i in 1:nm] + + # Concentrations + z = [ + 1.0000 0 0; + 0.7085 0.1621 0.0811; + 0.5971 0.1855 0.0965; + 0.5537 0.1989 0.1198; + 0.3684 0.2845 0.1535; + 0.1712 0.3491 0.2097; + 0.1198 0.3098 0.2628; + 0.0747 0.3576 0.2467; + 0.0529 0.3347 0.2884; + 0.0415 0.3388 0.2757; + 0.0261 0.3557 0.3167; + 0.0208 0.3483 0.2954; + 0.0085 0.3836 0.2950; + 0.0053 0.3611 0.2937; + 0.0019 0.3609 0.2831; + 0.0018 0.3485 0.2846; + 0.0006 0.3698 0.2899; + ] + + bc = [1.0, 0.0, 0.0] + + # Starting-value + v0 = zeros(n, ne) + for i in 1:itau[1], s in 1:ne + v0[i, s] = bc[s] + end + for j in 2:nm, i in itau[j-1]+1:itau[j], s in 1:ne + v0[i, s] = z[j, s] + end + for i in itau[nm]+1:n, s in 1:ne + v0[i, s] = z[nm, s] + end + v0 .= 0.001 + + model = Model() + + @variable(model, theta[1:np] >= 0.0, start=1.0) # ODE parameters + @variable(model, v[i=1:n, s=1:ne], start=v0[i, s]) + @variable(model, w[i=1:n, j=1:nc, s=1:ne], start=0.0) + @variable(model, uc[i=1:n, j=1:nc, s=1:ne], start=v0[i, s]) + @variable(model, Duc[i=1:n, j=1:nc, s=1:ne], start=0.0) + + # error + @NLexpression( + model, + error[j=1:nm, s=1:ne], + v[itau[j],s] + sum(w[itau[j],k,s]*(tau[j]-t[itau[j]])^k/(fact[k+1]*h^(k-1)) for k in 1:nc) - z[j,s] + ) + + # L2 error + @NLobjective(model, Min, sum(error[j, s]^2 for j in 1:nm, s in 1:ne)) + + # Collocation model + @constraint( + model, + [i=1:n, j=1:nc, s=1:ne], + uc[i, j, s] == v[i,s] + h*sum(w[i,k,s]*(rho[j]^k/fact[k+1]) for k in 1:nc), + ) + @constraint( + model, + [i=1:n, j=1:nc, s=1:ne], + Duc[i, j, s] == sum(w[i,k,s]*(rho[j]^(k-1)/fact[k]) for k in 1:nc), + ) + # Boundary + @constraint(model, bc_cond[s=1:ne], v[1, s] == bc[s]) + # Continuity + @constraint( + model, + [i=1:n-1, s=1:ne], + v[i, s] + sum(w[i, j, s]*h/fact[j+1] for j in 1:nc) == v[i+1, s] + ) + # Dynamics + @NLconstraint( + model, + collocation_eqn1[i=1:n, j=1:nc], + Duc[i,j,1] == - (2*theta[2] - (theta[1]*uc[i,j,2])/((theta[2]+theta[5])*uc[i,j,1]+uc[i,j,2]) + + theta[3] + theta[4])*uc[i,j,1] + ) + @NLconstraint( + model, + collocation_eqn2[i=1:n, j=1:nc], + Duc[i,j,2] == (theta[1]*uc[i,j,1]*(theta[2]*uc[i,j,1]-uc[i,j,2]))/ + ((theta[2]+theta[5])*uc[i,j,1]+uc[i,j,2]) + + theta[3]*uc[i,j,1] + ) + @NLconstraint( + model, + collocation_eqn3[i=1:n, j=1:nc], + Duc[i,j,3] == (theta[1]*uc[i,j,1]*(uc[i,j,2]+theta[5]*uc[i,j,1]))/ + ((theta[2]+theta[5])*uc[i,j,1]+uc[i,j,2]) + + theta[4]*uc[i,j,1] + ) + + return model +end + diff --git a/src/PureJuMP/minsurf.jl b/src/PureJuMP/minsurf.jl new file mode 100644 index 00000000..e8a60ae8 --- /dev/null +++ b/src/PureJuMP/minsurf.jl @@ -0,0 +1,57 @@ +# Minimal surface with obstacle problem + +# Find the surface with minimal area, given boundary conditions, +# and above an obstacle. + +# This is problem 17=the COPS (Version 3) collection of +# E. Dolan and J. More' +# see "Benchmarking Optimization Software with COPS" +# Argonne National Labs Technical Report ANL/MCS-246 (2004) +# classification OBR2-AN-V-V + + +function minsurf(args...; n = default_nvar, kwargs...) + nx, ny = n + x_mesh = LinRange(0, 1, nx + 2) # coordinates of the mesh points x + + v0 = zeros(nx + 2, ny + 2) # Surface matrix initialization + for i = 1:(nx + 2), j = 1:(ny + 2) + v0[i, j] = 1 - (2 * x_mesh[i] - 1)^2 + end + + hx = 1 / (nx + 1) + hy = 1 / (ny + 1) + area = 1 // 2 * hx * hy + + nlp = Model() + + @variable(nlp, v[i = 1:(nx + 2), j = 1:(ny + 2)], start = v0[i, j]) + + @objective( + nlp, + Min, + sum( + area * + (1 + ((v[i + 1, j] - v[i, j]) / hx)^2 + ((v[i, j + 1] - v[i, j]) / hy)^2)^(1 / 2) for + i = 1:(nx + 1), j = 1:(ny + 1) + ) + sum( + area * + (1 + ((v[i - 1, j] - v[i, j]) / hx)^2 + ((v[i, j - 1] - v[i, j]) / hy)^2)^(1 / 2) for + i = 2:(nx + 2), j = 2:(ny + 2) + ) + ) + + @constraint(nlp, [j = 0:(ny + 1)], v[1, j + 1] == 0) + @constraint(nlp, [j = 0:(ny + 1)], v[nx + 2, j + 1] == 0) + @constraint(nlp, [i = 0:(nx + 1)], v[i + 1, 1] == 1 - (2 * i * hx - 1)^2) + @constraint(nlp, [i = 0:(nx + 1)], v[i + 1, ny + 2] == 1 - (2 * i * hx - 1)^2) + @constraint(nlp, [i = 0:(nx + 1), j = 0:(ny + 1)], v[i + 1, j + 1] >= 0) + @constraint( + nlp, + [i = Int(floor(0.25 / hx)):Int(ceil(0.75 / hx)), j = Int(floor(0.25 / hy)):Int(ceil(0.75 / hy))], + v[i + 1, j + 1] >= 1 + ) + + return nlp +end + diff --git a/src/PureJuMP/pinene.jl b/src/PureJuMP/pinene.jl new file mode 100644 index 00000000..193041ef --- /dev/null +++ b/src/PureJuMP/pinene.jl @@ -0,0 +1,106 @@ +# This is problem 8 in the COPS (Version 3) collection of +# E. Dolan and J. More +# see "Benchmarking Optimization Software with COPS" +# Argonne National Labs Technical Report ANL/MCS-246 (2004) + +# Isomerization of Alpha-Pinene Problem +# Collocation formulation +# Alexander S. Bondarenko - Summer 1998 +# COPS 2.0 - September 2000 +# COPS 3.0 - November 2002 +# COPS 3.1 - March 2004 + +function pinene(; n::Int = default_nvar, kwargs...) + nc = 3 # number of collocation points + ne = 5 # number of differential equations + np = 5 # number of ODE parameters + nm = 8 # number of measurements + + # roots of k-th degree Legendre polynomial + rho = [0.11270166537926, 0.5, 0.88729833462074] + # boundary conditions + bc = [100.0, 0.0, 0.0, 0.0, 0.0] + # times at which observations made + tau = [1230.0, 3060.0, 4920.0, 7800.0, 10680.0, 15030.0, 22620.0, 36420.0] + tf = tau[nm] # ODEs defined in [0,tf] + h = tf / n # uniform interval length + t = [(i-1)*h for i in 1:n+1] # partition + + # itau[i] is the largest integer k with t[k] <= tau[i] + itau = Int[min(n, floor(tau[i]/h)+1) for i in 1:nm] + + fact = [factorial(k) for k in 0:nc] + + # Observations + z = [ + 88.35 7.3 2.3 0.4 1.75; + 76.4 15.6 4.5 0.7 2.8; + 65.1 23.1 5.3 1.1 5.8; + 50.4 32.9 6.0 1.5 9.3; + 37.5 42.7 6.0 1.9 12.0; + 25.9 49.1 5.9 2.2 17.0; + 14.0 57.4 5.1 2.6 21.0; + 4.5 63.1 3.8 2.9 25.7; + ] + + v0 = zeros(n, ne) + # Starting-value + for i in 1:itau[1], s in 1:ne + v0[i, s] = bc[s] + end + for j in 2:nm, i =itau[j-1]+1:itau[j], s in 1:ne + v0[i, s] = z[j, s] + end + for i in itau[nm]+1:n, s in 1:ne + v0[i, s] = z[nm, s] + end + + model = Model() + + @variable(model, theta[1:np] >= 0.0, start=0.0) + # The collocation approximation u is defined by the parameters v and w. + # uc and Duc are, respectively, u and u' evaluated at the collocation points. + @variable(model, v[i=1:n, s=1:ne], start=v0[i, s]) + @variable(model, w[1:n, 1:nc, 1:ne], start=0.0) + @variable(model, uc[i=1:n, j=1:nc, s=1:ne], start=v0[i,s]) + @variable(model, Duc[i=1:n, j=1:nc, s=1:ne], start=0.0) + + @NLexpression( + model, + error[j=1:nm, s=1:ne], + v[itau[j],s] + sum(w[itau[j],k,s]*(tau[j]-t[itau[j]])^k/(fact[k+1]*h^(k-1)) for k in 1:nc) - z[j,s] + ) + # l2 error + @NLobjective(model, Min, sum(error[j, s]^2 for j in 1:nm, s in 1:ne)) + + # Collocation model + @constraint( + model, + [i=1:n, j=1:nc, s=1:ne], + uc[i, j, s] == v[i,s] + h*sum(w[i,k,s]*(rho[j]^k/fact[k+1]) for k in 1:nc), + ) + @constraint( + model, + [i=1:n, j=1:nc, s=1:ne], + Duc[i, j, s] == sum(w[i,k,s]*(rho[j]^(k-1)/fact[k]) for k in 1:nc), + ) + # Boundary + @constraint(model, [s=1:ne], v[1, s] == bc[s]) + # Continuity + @constraint( + model, + [i=1:n-1, s=1:ne], + v[i, s] + sum(w[i, j, s]*h/fact[j+1] for j in 1:nc) == v[i+1, s], + ) + @NLconstraints( + model, begin + [i=1:n, j=1:nc], Duc[i,j,1] == - (theta[1]+theta[2])*uc[i,j,1] + [i=1:n, j=1:nc], Duc[i,j,2] == theta[1]*uc[i,j,1] + [i=1:n, j=1:nc], Duc[i,j,3] == theta[2]*uc[i,j,1] - (theta[3]+theta[4])*uc[i,j,3] + theta[5]*uc[i,j,5] + [i=1:n, j=1:nc], Duc[i,j,4] == theta[3]*uc[i,j,3] + [i=1:n, j=1:nc], Duc[i,j,5] == theta[4]*uc[i,j,3] - theta[5]*uc[i,j,5] + end + ) + return model +end + diff --git a/src/PureJuMP/rocket.jl b/src/PureJuMP/rocket.jl new file mode 100644 index 00000000..1486e02d --- /dev/null +++ b/src/PureJuMP/rocket.jl @@ -0,0 +1,58 @@ +# Goddard Rocket Problem +# Trapezoidal formulation +# COPS 2.0 - September 2000 +# COPS 3.0 - November 2002 +# COPS 3.1 - March 2004 + +function rocket(; n::Int = default_nvar, kwargs...) + h_0 = 1.0 + v_0 = 0.0 + m_0 = 1.0 + g_0 = 1.0 + T_c = 3.5 + h_c = 500.0 + v_c = 620.0 + m_c = 0.6 + + c = 0.5*sqrt(g_0 * h_0) + m_f = m_c * m_0 + D_c = 0.5 * v_c * (m_0 / g_0) + T_max = T_c * m_0 * g_0 + + model = Model() + + @variables(model, begin + 1.0 <= h[i=0:n], (start=1.0) + 0.0 <= v[i=0:n], (start=i/n*(1.0 - i/n)) + m_f <= m[i=0:n] <= m_0, (start=(m_f - m_0)*(i/n) + m_0) + 0.0 <= T[i=0:n] <= T_max, (start=T_max/2.0) + 0.0 <= step, (start=1/n) + end) + + @NLexpressions(model, begin + D[i=0:n], D_c*v[i]^2*exp(-h_c*(h[i] - h_0))/h_0 + g[i=0:n], g_0 * (h_0 / h[i])^2 + dh[i=0:n], v[i] + dv[i=0:n], (T[i] - D[i] - m[i]*g[i]) / m[i] + dm[i=0:n], -T[i]/c + end) + + @objective(model, Max, h[n]) + + # Dynamics + @NLconstraints(model, begin + con_dh[i=1:n], h[i] == h[i-1] + 0.5 * step * (dh[i] + dh[i-1]) + con_dv[i=1:n], v[i] == v[i-1] + 0.5 * step * (dv[i] + dv[i-1]) + con_dm[i=1:n], m[i] == m[i-1] + 0.5 * step * (dm[i] + dm[i-1]) + end) + # Boundary constraints + @constraints(model, begin + h_ic, h[0] == h_0 + v_ic, v[0] == v_0 + m_ic, m[0] == m_0 + m_fc, m[n] == m_f + end) + + return model +end + diff --git a/src/PureJuMP/steering.jl b/src/PureJuMP/steering.jl new file mode 100644 index 00000000..3618a0b0 --- /dev/null +++ b/src/PureJuMP/steering.jl @@ -0,0 +1,48 @@ +# Rocket Steering Problem +# Trapezoidal formulation +# COPS 2.0 - September 2000 +# COPS 3.0 - November 2002 +# COPS 3.1 - March 2004 + +function steering(; n::Int=default_nvar, kwargs...) + a = 100.0 # Magnitude of force. + # Bounds on the control + u_min, u_max = -pi/2.0, pi/2.0 + xs = zeros(4) + xf = [NaN, 5.0, 45.0, 0.0] + + function gen_x0(k, i) + if i == 1 || i == 4 + return 0.0 + elseif i == 2 + return 5*k/n + elseif i == 3 + return 45.0*k/n + end + end + + model = Model() + + @variable(model, u_min <= u[i=1:n+1] <= u_max, start=0.0) # control + @variable(model, x[i=1:n+1, j=1:4], start=gen_x0(i, j)) # state + @variable(model, tf, start=1.0) # final time + + @expression(model, h, tf / n) # step size + @objective(model, Min, tf) + + @constraint(model, tf >= 0.0) + # Dynamics + @NLconstraints( + model, begin + [i=1:n], x[i+1,1] == x[i,1] + 0.5*h*(x[i,3] + x[i+1,3]) + [i=1:n], x[i+1,2] == x[i,2] + 0.5*h*(x[i,4] + x[i+1,4]) + [i=1:n], x[i+1,3] == x[i,3] + 0.5*h*(a*cos(u[i]) + a*cos(u[i+1])) + [i=1:n], x[i+1,4] == x[i,4] + 0.5*h*(a*sin(u[i]) + a*sin(u[i+1])) + end + ) + # Boundary conditions + @constraint(model , [j=1:4], x[1, j] == xs[j]) + @constraint(model , [j=2:4], x[n+1, j] == xf[j]) + return model +end + diff --git a/src/PureJuMP/torsion.jl b/src/PureJuMP/torsion.jl new file mode 100644 index 00000000..965318ef --- /dev/null +++ b/src/PureJuMP/torsion.jl @@ -0,0 +1,47 @@ +# Torsion problem +# Liz Dolan - Summer 2000 +# Version 2.0 - October 2000 +# COPS 3.1 - March 2004 + +function torsion(args...; n = default_nvar, kwargs...) + nx, ny = n + c = 5.0 + hx = 1.0 / (nx + 1.0) # grid spacing + hy = 1.0 / (ny + 1.0) # grid spacing + area = 0.5 * hx * hy # area of triangle + + # Distance to the boundary. + D = [min(min(i,nx-i+1)*hx, min(j, ny-j+1)*hy) for i in 0:nx+1, j in 0:ny+1] + + model = Model() + # v defines the finite element approximation. + @variable(model, v[i=0:nx+1, j=0:ny+1], start=D[i+1, j+1]) + + @expression( + model, + linLower, + sum(v[i+1, j] + v[i, j] + v[i, j+1] for i in 0:nx, j in 0:ny) + ) + @expression( + model, + linUpper, + sum(v[i, j] + v[i-1, j] + v[i, j-1] for i in 1:nx+1, j in 1:ny+1) + ) + @expression( + model, + quadLower, + sum(((v[i+1,j] - v[i,j])/hx)^2 + ((v[i,j+1] - v[i,j])/hy)^2 for i in 0:nx, j in 0:ny) + ) + @expression( + model, + quadUpper, + sum(((v[i,j] - v[i-1,j])/hx)^2 + ((v[i,j] - v[i,j-1])/hy)^2 for i in 1:nx+1, j in 1:ny+1) + ) + + @objective(model, Min, area * ((quadLower + quadUpper) / 2.0 - c * (linLower + linUpper) / 3.0)) + + @constraint(model, [i=0:nx+1, j=0:ny+1], -D[i+1, j+1] <= v[i, j] <= D[i+1, j+1]) + + return model +end +