1 | <!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN" |
---|
2 | "http://www.w3.org/TR/REC-html40/loose.dtd"> |
---|
3 | <html> |
---|
4 | <head> |
---|
5 | <title>Description of gettwiss</title> |
---|
6 | <meta name="keywords" content="gettwiss"> |
---|
7 | <meta name="description" content="GETWISS - Calculate the twiss parameters"> |
---|
8 | <meta http-equiv="Content-Type" content="text/html; charset=iso-8859-1"> |
---|
9 | <meta name="generator" content="m2html © 2003 Guillaume Flandin"> |
---|
10 | <meta name="robots" content="index, follow"> |
---|
11 | <link type="text/css" rel="stylesheet" href="../m2html.css"> |
---|
12 | </head> |
---|
13 | <body> |
---|
14 | <a name="_top"></a> |
---|
15 | <div><a href="../index.html">Home</a> > <a href="index.html">at</a> > gettwiss.m</div> |
---|
16 | |
---|
17 | <!--<table width="100%"><tr><td align="left"><a href="../index.html"><img alt="<" border="0" src="../left.png"> Master index</a></td> |
---|
18 | <td align="right"><a href="index.html">Index for at <img alt=">" border="0" src="../right.png"></a></td></tr></table>--> |
---|
19 | |
---|
20 | <h1>gettwiss |
---|
21 | </h1> |
---|
22 | |
---|
23 | <h2><a name="_name"></a>PURPOSE <a href="#_top"><img alt="^" border="0" src="../up.png"></a></h2> |
---|
24 | <div class="box"><strong>GETWISS - Calculate the twiss parameters</strong></div> |
---|
25 | |
---|
26 | <h2><a name="_synopsis"></a>SYNOPSIS <a href="#_top"><img alt="^" border="0" src="../up.png"></a></h2> |
---|
27 | <div class="box"><strong>function [Optics] = gettwiss(THERING, DP) </strong></div> |
---|
28 | |
---|
29 | <h2><a name="_description"></a>DESCRIPTION <a href="#_top"><img alt="^" border="0" src="../up.png"></a></h2> |
---|
30 | <div class="fragment"><pre class="comment">GETWISS - Calculate the twiss parameters |
---|
31 | [Optics] = gettwiss(THERING, DP) |
---|
32 | |
---|
33 | GETTWISS calls LINOPT2 for entire ring and returns Twiss parameters |
---|
34 | See LinOpt2 for description of nomenclarture and mathematics |
---|
35 | |
---|
36 | Phase calculation from transfer matrices from '0' to location '1' (AIP 184 p. 50) |
---|
37 | Only need initial alfa, beta |
---|
38 | M=M44 |
---|
39 | M11: sqrt(beta1/beta0)*(cos(phi) + alfa0*sin(phi)) |
---|
40 | M12: sqrt(beta0*beta1)*sin(phi) |
---|
41 | M21: (1/sqrt(beta0*beta1))*[(alfa0-alfa1)*cos(phi)-(1+alfa0*alfa1)*sin(phi)) |
---|
42 | M22: sqrt(beta0/beta1)*(cos(phi)-alfa1*sin(phi)) |
---|
43 | phi=atan2(M12,beta0*M11-alfa0*M12) |
---|
44 | |
---|
45 | Written by Jeff Corbett</pre></div> |
---|
46 | |
---|
47 | <!-- crossreference --> |
---|
48 | <h2><a name="_cross"></a>CROSS-REFERENCE INFORMATION <a href="#_top"><img alt="^" border="0" src="../up.png"></a></h2> |
---|
49 | This function calls: |
---|
50 | <ul style="list-style-image:url(../matlabicon.gif)"> |
---|
51 | <li><a href="linopt2.html" class="code" title="function [LinData, varargout] = linopt2(RING,DP,varargin);">linopt2</a> LINOPT2 performs linear analysis of the COUPLED lattice RING</li></ul> |
---|
52 | This function is called by: |
---|
53 | <ul style="list-style-image:url(../matlabicon.gif)"> |
---|
54 | <li><a href="modelDA.html" class="code" title="function [DA, Data] = modelDA( r0, nsteps, nturns, dp, res)">modelDA</a> modelDA( r0, nsteps, nturns, dp, res)</li></ul> |
---|
55 | <!-- crossreference --> |
---|
56 | |
---|
57 | |
---|
58 | <h2><a name="_source"></a>SOURCE CODE <a href="#_top"><img alt="^" border="0" src="../up.png"></a></h2> |
---|
59 | <div class="fragment"><pre>0001 <a name="_sub0" href="#_subfunctions" class="code">function [Optics] = gettwiss(THERING, DP)</a> |
---|
60 | 0002 <span class="comment">%GETWISS - Calculate the twiss parameters</span> |
---|
61 | 0003 <span class="comment">% [Optics] = gettwiss(THERING, DP)</span> |
---|
62 | 0004 <span class="comment">%</span> |
---|
63 | 0005 <span class="comment">% GETTWISS calls LINOPT2 for entire ring and returns Twiss parameters</span> |
---|
64 | 0006 <span class="comment">% See LinOpt2 for description of nomenclarture and mathematics</span> |
---|
65 | 0007 <span class="comment">%</span> |
---|
66 | 0008 <span class="comment">% Phase calculation from transfer matrices from '0' to location '1' (AIP 184 p. 50)</span> |
---|
67 | 0009 <span class="comment">% Only need initial alfa, beta</span> |
---|
68 | 0010 <span class="comment">% M=M44</span> |
---|
69 | 0011 <span class="comment">% M11: sqrt(beta1/beta0)*(cos(phi) + alfa0*sin(phi))</span> |
---|
70 | 0012 <span class="comment">% M12: sqrt(beta0*beta1)*sin(phi)</span> |
---|
71 | 0013 <span class="comment">% M21: (1/sqrt(beta0*beta1))*[(alfa0-alfa1)*cos(phi)-(1+alfa0*alfa1)*sin(phi))</span> |
---|
72 | 0014 <span class="comment">% M22: sqrt(beta0/beta1)*(cos(phi)-alfa1*sin(phi))</span> |
---|
73 | 0015 <span class="comment">% phi=atan2(M12,beta0*M11-alfa0*M12)</span> |
---|
74 | 0016 <span class="comment">%</span> |
---|
75 | 0017 <span class="comment">% Written by Jeff Corbett</span> |
---|
76 | 0018 |
---|
77 | 0019 |
---|
78 | 0020 <span class="keyword">if</span> nargin < 1 |
---|
79 | 0021 <span class="keyword">global</span> THERING |
---|
80 | 0022 <span class="keyword">end</span> |
---|
81 | 0023 <span class="keyword">if</span> nargin < 2 |
---|
82 | 0024 DP = 1e-8; |
---|
83 | 0025 <span class="keyword">end</span> |
---|
84 | 0026 |
---|
85 | 0027 |
---|
86 | 0028 <span class="comment">%load linear optics structure for entire ring</span> |
---|
87 | 0029 <span class="comment">%disp(' Computing Coupled Lattice Parameters ...');</span> |
---|
88 | 0030 LinData = <a href="linopt2.html" class="code" title="function [LinData, varargout] = linopt2(RING,DP,varargin);">linopt2</a>(THERING,DP,1:length(THERING)+1); |
---|
89 | 0031 |
---|
90 | 0032 <span class="comment">%compute fractional tune</span> |
---|
91 | 0033 nux=acos(trace(LinData(1).A/2))/(2*pi); |
---|
92 | 0034 nuy=acos(trace(LinData(1).B/2))/(2*pi); |
---|
93 | 0035 |
---|
94 | 0036 <span class="comment">%pre-define arrays</span> |
---|
95 | 0037 NR=length(THERING)+1; |
---|
96 | 0038 name=cell(NR,1); <span class="comment">%initialize cell for names</span> |
---|
97 | 0039 s=zeros(NR,1); |
---|
98 | 0040 len=s; |
---|
99 | 0041 betax=s; betay=s; |
---|
100 | 0042 alfax=s; alfay=s; |
---|
101 | 0043 gammax=s; gammay=s; |
---|
102 | 0044 etax=s; etay=s; |
---|
103 | 0045 etapx=s; etapy=s; |
---|
104 | 0046 phix=s; phiy=s; |
---|
105 | 0047 curlhx=s; curlhy=s; |
---|
106 | 0048 ch1=s; ch2=s; ch3=s; |
---|
107 | 0049 |
---|
108 | 0050 <span class="comment">%compute optics for all elements</span> |
---|
109 | 0051 <span class="keyword">for</span> ii=1:NR |
---|
110 | 0052 |
---|
111 | 0053 <span class="comment">%load position</span> |
---|
112 | 0054 s(ii)=LinData(ii).SPos(1); |
---|
113 | 0055 |
---|
114 | 0056 <span class="comment">%horizontal Twiss parameters</span> |
---|
115 | 0057 m22=LinData(ii).A; |
---|
116 | 0058 a=m22(1,1); b=m22(1,2); c=m22(2,1); d=m22(2,2); |
---|
117 | 0059 betax(ii)=b/sin(2*pi*nux); |
---|
118 | 0060 alfax(ii)=(a-d)/(2*sin(2*pi*nux)); |
---|
119 | 0061 gammax(ii)=-c/sin(2*pi*nux); |
---|
120 | 0062 |
---|
121 | 0063 <span class="comment">%vertical Twiss parameters</span> |
---|
122 | 0064 m22=LinData(ii).B; |
---|
123 | 0065 a=m22(1,1); b=m22(1,2); c=m22(2,1); d=m22(2,2); |
---|
124 | 0066 betay(ii)=b/sin(2*pi*nuy); |
---|
125 | 0067 alfay(ii)=(a-d)/(2*sin(2*pi*nuy)); |
---|
126 | 0068 gammay(ii)=-c/sin(2*pi*nuy); |
---|
127 | 0069 |
---|
128 | 0070 <span class="comment">%horizontal phase</span> |
---|
129 | 0071 num=LinData(ii).M44(1,2); |
---|
130 | 0072 den=LinData(1).beta(1)*LinData(ii).M44(1,1)-LinData(1).alfa(1)*num; |
---|
131 | 0073 phix(ii)=atan2(num,den); |
---|
132 | 0074 |
---|
133 | 0075 <span class="comment">%vertical phase</span> |
---|
134 | 0076 num=LinData(ii).M44(3,4); |
---|
135 | 0077 den=LinData(1).beta(2)*LinData(ii).M44(3,3)-LinData(1).alfa(2)*num; |
---|
136 | 0078 phiy(ii)=atan2(num,den); |
---|
137 | 0079 |
---|
138 | 0080 <span class="comment">%load element names</span> |
---|
139 | 0081 <span class="keyword">if</span> ii<NR |
---|
140 | 0082 <span class="keyword">if</span> isfield(THERING{ii},<span class="string">'Name'</span>) |
---|
141 | 0083 name{ii}=THERING{ii}.Name; |
---|
142 | 0084 <span class="keyword">else</span> |
---|
143 | 0085 name{ii}=THERING{ii}.FamName; |
---|
144 | 0086 <span class="keyword">end</span> |
---|
145 | 0087 |
---|
146 | 0088 <span class="comment">%load element lengths</span> |
---|
147 | 0089 len(ii)=THERING{ii}.Length; |
---|
148 | 0090 <span class="keyword">end</span> |
---|
149 | 0091 <span class="keyword">end</span> |
---|
150 | 0092 |
---|
151 | 0093 <span class="keyword">if</span> ii==NR |
---|
152 | 0094 name{NR}=<span class="string">'End'</span>; <span class="comment">%last element has name 'End' and zero length</span> |
---|
153 | 0095 len(NR)=0; <span class="comment">%THERING only has NR-1 elements</span> |
---|
154 | 0096 <span class="keyword">end</span> |
---|
155 | 0097 |
---|
156 | 0098 <span class="comment">%compute dispersion, derivative</span> |
---|
157 | 0099 <span class="comment">%disp(' Computing Dispersion ...');</span> |
---|
158 | 0100 dp=1e-6; |
---|
159 | 0101 orb4=findorbit4(THERING,dp,1:NR); |
---|
160 | 0102 etax =orb4(1,:)'/dp; |
---|
161 | 0103 etapx=orb4(2,:)'/dp; |
---|
162 | 0104 etay =orb4(3,:)'/dp; |
---|
163 | 0105 etapy=orb4(4,:)'/dp; |
---|
164 | 0106 |
---|
165 | 0107 <span class="comment">%unwrap phase</span> |
---|
166 | 0108 phix=unwrap(phix)/(2*pi); |
---|
167 | 0109 phiy=unwrap(phiy)/(2*pi); |
---|
168 | 0110 |
---|
169 | 0111 <span class="comment">%compute integer tunes</span> |
---|
170 | 0112 qx=round(phix(NR)); |
---|
171 | 0113 <span class="keyword">if</span> phix(NR)-qx <=0 |
---|
172 | 0114 qx=qx-1; <span class="comment">% Below half integer</span> |
---|
173 | 0115 <span class="keyword">end</span> |
---|
174 | 0116 qx=qx+nux; |
---|
175 | 0117 qy=round(phiy(NR)); |
---|
176 | 0118 <span class="keyword">if</span> phiy(NR)-qy <=0 |
---|
177 | 0119 qy=qy-1; <span class="comment">% Below half integer</span> |
---|
178 | 0120 <span class="keyword">end</span> |
---|
179 | 0121 qy=qy+nuy; |
---|
180 | 0122 |
---|
181 | 0123 <span class="comment">%compute curly-H</span> |
---|
182 | 0124 <span class="keyword">for</span> ii=1:NR; |
---|
183 | 0125 <span class="comment">%horizontal</span> |
---|
184 | 0126 ch1(ii)=gammax(ii)*etax(ii)*etax(ii); |
---|
185 | 0127 ch2(ii)=2*alfax(ii)*etax(ii)*etapx(ii); |
---|
186 | 0128 ch3(ii)=betax(ii)*etapx(ii)*etapx(ii); |
---|
187 | 0129 curlhx(ii)= ch1(ii)+ch2(ii)+ch3(ii); |
---|
188 | 0130 <span class="comment">%vertical</span> |
---|
189 | 0131 ch1(ii)=gammay(ii)*etay(ii)*etay(ii); |
---|
190 | 0132 ch2(ii)=2*alfay(ii)*etay(ii)*etapy(ii); |
---|
191 | 0133 ch3(ii)=betay(ii)*etapy(ii)*etapy(ii); |
---|
192 | 0134 curlhy(ii)= ch1(ii)+ch2(ii)+ch3(ii); |
---|
193 | 0135 <span class="keyword">end</span> |
---|
194 | 0136 |
---|
195 | 0137 <span class="comment">%load output</span> |
---|
196 | 0138 Optics.name =char(name); |
---|
197 | 0139 Optics.len =len; |
---|
198 | 0140 Optics.s =s; |
---|
199 | 0141 Optics.betax =betax; |
---|
200 | 0142 Optics.alfax =alfax; |
---|
201 | 0143 Optics.gammax=gammax; |
---|
202 | 0144 Optics.phix =phix; |
---|
203 | 0145 <span class="comment">%Optics.qx =qx;</span> |
---|
204 | 0146 <span class="comment">%Optics.nux =nux;</span> |
---|
205 | 0147 Optics.etax =etax; |
---|
206 | 0148 Optics.etapx =etapx; |
---|
207 | 0149 Optics.curlhx=curlhx; |
---|
208 | 0150 <span class="comment">%load transport matrices 'A'</span> |
---|
209 | 0151 <span class="keyword">for</span> ii=1:NR |
---|
210 | 0152 Optics.r11x(ii,1)=LinData(ii).M44(1,1); |
---|
211 | 0153 Optics.r12x(ii,1)=LinData(ii).M44(1,2); |
---|
212 | 0154 Optics.r21x(ii,1)=LinData(ii).M44(2,1); |
---|
213 | 0155 Optics.r22x(ii,1)=LinData(ii).M44(2,2); |
---|
214 | 0156 <span class="keyword">end</span> |
---|
215 | 0157 |
---|
216 | 0158 Optics.betay =betay; |
---|
217 | 0159 Optics.alfay =alfay; |
---|
218 | 0160 Optics.gammay=gammay; |
---|
219 | 0161 Optics.phiy =phiy; |
---|
220 | 0162 <span class="comment">%Optics.qy =qy;</span> |
---|
221 | 0163 <span class="comment">%Optics.nuy =nuy;</span> |
---|
222 | 0164 Optics.etay =etay; |
---|
223 | 0165 Optics.etapy =etapy; |
---|
224 | 0166 Optics.curlhy=curlhy; |
---|
225 | 0167 <span class="comment">%load transport matrices 'B'</span> |
---|
226 | 0168 <span class="keyword">for</span> ii=1:NR |
---|
227 | 0169 Optics.r11y(ii,1)=LinData(ii).M44(3,3); |
---|
228 | 0170 Optics.r12y(ii,1)=LinData(ii).M44(3,4); |
---|
229 | 0171 Optics.r21y(ii,1)=LinData(ii).M44(4,3); |
---|
230 | 0172 Optics.r22y(ii,1)=LinData(ii).M44(4,4); |
---|
231 | 0173 <span class="keyword">end</span> |
---|
232 | 0174 |
---|
233 | 0175 <span class="comment">% plot(s,phix); hold on; plot(s,phiy,'r');</span> |
---|
234 | 0176 <span class="comment">% figure;</span> |
---|
235 | 0177 <span class="comment">% ip=1:round(NR/4);</span> |
---|
236 | 0178 <span class="comment">% plot(s(ip),betax(ip)); hold on; plot(s(ip),betay(ip),'r'); plot(s(ip),10*Optics.etax(ip),'k');</span> |
---|
237 | 0179 <span class="comment">% figure</span> |
---|
238 | 0180 <span class="comment">% plot(s(ip),curlhx(ip));</span> |
---|
239 | 0181 |
---|
240 | 0182 <span class="comment">%disp([' Horizontal Tune: ', num2str(qx,'%6.3f')]);</span> |
---|
241 | 0183 <span class="comment">%disp([' Vertical Tune: ', num2str(qy,'%6.3f')]);</span></pre></div> |
---|
242 | <hr><address>Generated on Fri 18-May-2007 17:13:39 by <strong><a href="http://www.artefact.tk/software/matlab/m2html/">m2html</a></strong> © 2003</address> |
---|
243 | </body> |
---|
244 | </html> |
---|