|
| 1 | +function choreo |
| 2 | +%CHOREO Compute planar choreographies of the n-body problem. |
| 3 | +% CHEB.CHOREO computes a planar choreography using hand-drawn initial |
| 4 | +% guesses. At the end of the computation, the user is asked to press <1> to |
| 5 | +% display the motion of the planets. To stop the program, press <CTRL>-<C>. |
| 6 | +% |
| 7 | +% Planar choreographies are periodic solutions of the n-body problem in which |
| 8 | +% the bodies share a single orbit. This orbit can be fixed or rotating with some |
| 9 | +% angular velocity relative to an inertial reference frame. |
| 10 | +% |
| 11 | +% The algorithm uses trigonometric interpolation and quasi-Newton methods. |
| 12 | +% See [1] for details. |
| 13 | +% |
| 14 | +% Example: at the prompt, specify 5 bodies, angular rotation 0 or 1.2. |
| 15 | +% Then either draw a curve (if your machine has imfreehand) or click |
| 16 | +% in 10 or 15 points (if it doesn't) roughly along a figure-8. Type |
| 17 | +% <Enter>, and in a few seconds you will see a choreograpy. |
| 18 | +% |
| 19 | +% [1] H. Montanelli and N. I. Gushterov, Computing planar and spherical |
| 20 | +% choreographies, SIAM Journal on Applied Dynamical Systems 15 (2016), |
| 21 | + |
| 22 | +% Copyright 2017 by The University of Oxford and The Chebfun Developers. |
| 23 | +% See http://www.chebfun.org/ for Chebfun information. |
| 24 | + |
| 25 | +close all |
| 26 | +LW = 'linewidth'; MS = 'markersize'; FS = 'fontsize'; |
| 27 | +lw = 2; ms = 30; fs = 18; |
| 28 | +format long, format compact |
| 29 | +n = input('How many bodies? (e.g. 5) '); |
| 30 | +w = input('Angular velocity? (either 0 or a nonzero noninteger) '); |
| 31 | +k = 15; |
| 32 | +N = k*n; |
| 33 | +dom = [0 2*pi]; |
| 34 | + |
| 35 | +while 1 |
| 36 | + |
| 37 | + % Hand-drawn initial guess: |
| 38 | + xmax = 2.5; ymax = 2.5; |
| 39 | + clf, hold off, axis equal, axis([-xmax xmax -ymax ymax]), box on |
| 40 | + set(gca,FS,fs), title('Draw a curve!') |
| 41 | + if ( exist('imfreehand') == 2 ) |
| 42 | + % Try to use IMFREEHAND (need the image processing toolbox): |
| 43 | + h = imfreehand(); |
| 44 | + z = getPosition(h); |
| 45 | + delete(h) |
| 46 | + z = z(:,1) + 1i*z(:,2); |
| 47 | + q0 = chebfun(z,dom,'trig'); |
| 48 | + q0 = chebfun(q0,dom,N,'trig'); |
| 49 | + c0 = trigcoeffs(q0); |
| 50 | + c0(1+floor(N/2)) = 0; |
| 51 | + q0 = chebfun(c0,dom,'coeffs','trig'); |
| 52 | + else |
| 53 | + % Otherwise, use GINPUT: |
| 54 | + x = []; y = []; button = 1; |
| 55 | + disp('Input points with mouse, press <enter> for final point.') |
| 56 | + while ( button == 1 ) |
| 57 | + [xx,yy,button] = ginput(1); |
| 58 | + x = [x; xx]; y = [y; yy]; %#ok<*AGROW> |
| 59 | + hold on, plot(xx,yy,'xb',MS,ms/2), drawnow |
| 60 | + end |
| 61 | + z = x + 1i*y; |
| 62 | + q0 = chebfun(z,dom,'trig'); |
| 63 | + q0 = chebfun(q0,dom,N,'trig'); |
| 64 | + end |
| 65 | + |
| 66 | + % Solve the problem: |
| 67 | + hold on, plot(q0,'.-b',LW,lw), drawnow |
| 68 | + c0 = trigcoeffs(q0); |
| 69 | + c0 = [real(c0);imag(c0)]; |
| 70 | + [A0,G0] = actiongradeval(c0,n,w); |
| 71 | + fprintf('\nInitial action: %.6f\n',A0) |
| 72 | + options = optimoptions('fminunc'); |
| 73 | + options.Algorithm = 'quasi-newton'; |
| 74 | + options.HessUpdate = 'bfgs'; |
| 75 | + options.GradObj = 'on'; |
| 76 | + options.Display = 'off'; |
| 77 | + [c,A,~, ~,G] = fminunc(@(x)actiongradeval(x,n,w),c0,options); |
| 78 | + fprintf('Action after optimization: %.6f\n',A) |
| 79 | + fprintf('Norm of the gradient: %.3e\n',norm(G)/norm(G0)) |
| 80 | + |
| 81 | + % Plot the result: |
| 82 | + c = c(1:N) + 1i*c(N+1:2*N); |
| 83 | + c(1+floor(N/2)) = 0; |
| 84 | + q = chebfun(c,dom,'coeffs','trig'); |
| 85 | + hold on, plot(q,'r',LW,lw) |
| 86 | + hold on, plot(q(2*pi*(0:n-1)/n),'.k',MS,ms), title('') |
| 87 | + s = input('Want to see the planets dancing? (Yes=<1>, No=<0>) '); |
| 88 | + if ( s == 1 ) |
| 89 | + hold off, plotonplane(q,n,w), pause |
| 90 | + end |
| 91 | + |
| 92 | +end |
| 93 | + |
| 94 | +end |
| 95 | + |
| 96 | +%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
| 97 | +function [A,G] = actiongradeval(c,n,w) |
| 98 | +%ACTIONGRADEVAL Compute the action and its gradient in the plane. |
| 99 | + |
| 100 | +% Set up: |
| 101 | +if ( nargin < 3 ) |
| 102 | + w = 0; |
| 103 | +end |
| 104 | +N = length(c)/2; |
| 105 | +u = c(1:N); |
| 106 | +v = c(N+1:end); |
| 107 | +c = u + 1i*v; |
| 108 | +if ( mod(N, 2) == 0 ) |
| 109 | + k = (-N/2:N/2-1)'; |
| 110 | +else |
| 111 | + k = (-(N-1)/2:(N-1)/2)'; |
| 112 | +end |
| 113 | +[t,s] = trigpts(N,[0 2*pi]); |
| 114 | + |
| 115 | +% Evaluate action A: |
| 116 | +c = bsxfun(@times,exp(1i*k*(0:n-1)*2*pi/n),c(:,1)); |
| 117 | +vals = ifft(N*c); |
| 118 | +dc = 1i*k.*c(:,1); |
| 119 | +if ( mod(N, 2) == 0 ) |
| 120 | + dc(1,1) = 0; |
| 121 | +end |
| 122 | +dvals = ifft(N*dc); |
| 123 | +K = abs(dvals+1i*w*vals(:,1)).^2; |
| 124 | +U = 0; |
| 125 | +for i = 2:n |
| 126 | + U = U - 1./abs(vals(:,1)-vals(:,i)); |
| 127 | +end |
| 128 | +A = n/2*s*(K-U); |
| 129 | + |
| 130 | +if ( nargout > 1 ) |
| 131 | + |
| 132 | + % Initialize gradient G: |
| 133 | + G = zeros(2*N,1); |
| 134 | + |
| 135 | + % Loop over the bodies for the AU contribution: |
| 136 | + for j = 1:n-1 |
| 137 | + a = bsxfun(@times,1-cos(k'*j*2*pi/n),cos(t*k')) + ... |
| 138 | + bsxfun(@times,sin(k'*j*2*pi/n),sin(t*k')); |
| 139 | + b = bsxfun(@times,-1+cos(k'*j*2*pi/n),sin(t*k')) + ... |
| 140 | + bsxfun(@times,sin(k'*j*2*pi/n),cos(t*k')); |
| 141 | + f = (a*u + b*v).^2 + (-b*u + a*v).^2; |
| 142 | + df = 2*(bsxfun(@times,a*u + b*v,a) + bsxfun(@times,-b*u + a*v,-b)); |
| 143 | + G(1:N) = G(1:N) - n/4*bsxfun(@rdivide,df,f.^(3/2))'*s'; |
| 144 | + df = 2*(bsxfun(@times,a*u + b*v, b) + bsxfun(@times,-b*u + a*v,a)); |
| 145 | + G(N+1:2*N) = G(N+1:2*N) - n/4*bsxfun(@rdivide,df,f.^(3/2))'*s'; |
| 146 | + end |
| 147 | + |
| 148 | + % Add the AK contribution: |
| 149 | + if ( mod(N, 2) == 0 ) |
| 150 | + G = G + 2*pi*n*[w;(k(2:end)+w);w;(k(2:end)+w)].^2.*[u;v]; |
| 151 | + else |
| 152 | + G = G + 2*pi*n*[(k+w);(k+w)].^2.*[u;v]; |
| 153 | + end |
| 154 | +end |
| 155 | + |
| 156 | +end |
| 157 | + |
| 158 | +%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
| 159 | +function plotonplane(q,n,w) |
| 160 | +%PLOTONPLANE Plot a planar choreography. |
| 161 | + |
| 162 | +vscale = max(max(real(q)),max(imag(q))); |
| 163 | +xmax = 1.3*vscale; ymax = xmax; |
| 164 | +LW = 'linewidth'; MS = 'markersize'; FS = 'fontsize'; |
| 165 | +lw = 2; ms = 30; fs = 18; |
| 166 | +dom = [0 2*pi]; T = 2*pi; dt = .1; |
| 167 | +xStars = 2*xmax*rand(250,1)-xmax; |
| 168 | +yStars = 2*ymax*rand(250,1)-ymax; |
| 169 | +for t = dt:dt:10*T |
| 170 | + clf, hold off |
| 171 | + fill(xmax*[-1 1 1 -1 -1],ymax*[-1 -1 1 1 -1],'k') |
| 172 | + axis equal, axis([-xmax xmax -ymax ymax]), box on |
| 173 | + set(gca,FS,fs) |
| 174 | + hold on, plot(xStars,yStars,'.w',MS,4) |
| 175 | + if ( w ~= 0 ) |
| 176 | + q = chebfun(@(t)exp(1i*w*dt).*q(t),dom,length(q),'trig'); |
| 177 | + end |
| 178 | + if ( t < T ) |
| 179 | + hold on, plot(q,'r',LW,lw) |
| 180 | + end |
| 181 | + for j = 1:n |
| 182 | + hold on, plot(q(t+2*pi*j/n),'.',MS,ms) |
| 183 | + end |
| 184 | + drawnow |
| 185 | +end |
| 186 | + |
| 187 | +end |
0 commit comments