00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066 #include "Riostream.h"
00067 #include "TROOT.h"
00068 #include "TVirtualPad.h"
00069 #include "THelix.h"
00070 #include "TClass.h"
00071 #include "TMath.h"
00072
00073 Int_t THelix::fgMinNSeg=5;
00074
00075 ClassImp(THelix)
00076
00077
00078
00079 void THelix::SetHelix(Double_t *p, Double_t *v, Double_t w,
00080 Double_t *range, EHelixRangeType rType,
00081 Double_t *axis )
00082 {
00083
00084
00085
00086 SetAxis(axis);
00087
00088
00089 fW = w;
00090 Double_t * m = fRotMat->GetMatrix();
00091 Double_t vx0, vy0, vz0;
00092 vx0 = v[0] * m[0] + v[1] * m[1] + v[2] * m[2];
00093 vy0 = v[0] * m[3] + v[1] * m[4] + v[2] * m[5];
00094 vz0 = v[0] * m[6] + v[1] * m[7] + v[2] * m[8];
00095 fVt = TMath::Sqrt(vx0*vx0 + vy0*vy0);
00096 fPhi0 = TMath::ATan2(vy0,vx0);
00097 fVz = vz0;
00098 fX0 = p[0] * m[0] + p[1] * m[1] + p[2] * m[2];
00099 fY0 = p[0] * m[3] + p[1] * m[4] + p[2] * m[5];
00100 fZ0 = p[0] * m[6] + p[1] * m[7] + p[2] * m[8];
00101 if (fW != 0) {
00102 fX0 += fVt / fW * TMath::Sin(fPhi0);
00103 fY0 -= fVt / fW * TMath::Cos(fPhi0);
00104 }
00105
00106
00107 Double_t r1 = 0;
00108 Double_t r2 = 1;
00109 if (range) {r1 = range[0]; r2 = range[1];}
00110 if (rType != kUnchanged) {
00111 fRange[0] = 0.0; fRange[1] = TMath::Pi();
00112 SetRange(r1,r2,rType);
00113 }
00114 }
00115
00116
00117
00118 THelix::THelix()
00119 {
00120
00121
00122 fX0 = fY0 = fZ0 = fVt = fPhi0 = fVz = fAxis[0] = fAxis[1] = 0.0;
00123 fAxis[2] = 1.0;
00124 fW = 1.5E7;
00125 fRange[0] = 0.0;
00126 fRange[1] = 1.0;
00127 fRotMat = 0;
00128 }
00129
00130
00131
00132 THelix::THelix(Double_t x, Double_t y, Double_t z,
00133 Double_t vx, Double_t vy, Double_t vz,
00134 Double_t w)
00135 : TPolyLine3D()
00136 {
00137
00138 Double_t p[3], v[3];
00139 p[0] = x;
00140 p[1] = y;
00141 p[2] = z;
00142 v[0] = vx;
00143 v[1] = vy;
00144 v[2] = vz;
00145 Double_t *range = 0;
00146 fRotMat = 0;
00147
00148 SetHelix(p, v, w, range, kHelixZ);
00149 fOption = "";
00150 }
00151
00152
00153
00154 THelix::THelix(Double_t * p, Double_t * v, Double_t w,
00155 Double_t * range, EHelixRangeType rType, Double_t * axis)
00156 : TPolyLine3D()
00157 {
00158
00159
00160 Double_t r[2];
00161 if ( range ) {
00162 r[0] = range[0]; r[1] = range[1];
00163 } else {
00164 r[0] = 0.0; r[1] = 1.0;
00165 }
00166
00167 fRotMat = 0;
00168 if ( axis ) {
00169 SetHelix(p, v, w, r, rType, axis);
00170 } else {
00171 SetHelix(p, v, w, r, rType);
00172 }
00173
00174 fOption = "";
00175 }
00176
00177
00178 #if 0
00179
00180 THelix::THelix(const THelix &h) : TPolyLine3D()
00181 {
00182
00183
00184 fX0 = h.fX0;
00185 fY0 = h.fY0;
00186 fZ0 = h.fZ0;
00187 fVt = h.fVt;
00188 fPhi0 = h.fPhi0;
00189 fVz = h.fVz;
00190 fW = h.fW;
00191 for (Int_t i=0; i<3; i++) fAxis[i] = h.fAxis[i];
00192 fRotMat = new TRotMatrix(*(h.fRotMat));
00193 fRange[0] = h.fRange[0];
00194 fRange[1] = h.fRange[1];
00195
00196 fOption = h.fOption;
00197 }
00198 #endif
00199
00200
00201 THelix& THelix::operator=(const THelix& hx)
00202 {
00203
00204 if(this!=&hx) {
00205 TPolyLine3D::operator=(hx);
00206 fX0=hx.fX0;
00207 fY0=hx.fY0;
00208 fZ0=hx.fZ0;
00209 fVt=hx.fVt;
00210 fPhi0=hx.fPhi0;
00211 fVz=hx.fVz;
00212 fW=hx.fW;
00213 for(Int_t i=0; i<3; i++)
00214 fAxis[i]=hx.fAxis[i];
00215 fRotMat=hx.fRotMat;
00216 for(Int_t i=0; i<2; i++)
00217 fRange[i]=hx.fRange[i];
00218 }
00219 return *this;
00220 }
00221
00222
00223 THelix::~THelix()
00224 {
00225
00226
00227 if (fRotMat) delete fRotMat;
00228 }
00229
00230
00231
00232 THelix::THelix(const THelix &helix) : TPolyLine3D(helix)
00233 {
00234
00235
00236 fRotMat=0;
00237 ((THelix&)helix).THelix::Copy(*this);
00238 }
00239
00240
00241
00242 void THelix::Copy(TObject &obj) const
00243 {
00244
00245
00246 TObject::Copy(obj);
00247 TAttLine::Copy(((THelix&)obj));
00248
00249 ((THelix&)obj).fX0 = fX0;
00250 ((THelix&)obj).fY0 = fY0;
00251 ((THelix&)obj).fZ0 = fZ0;
00252 ((THelix&)obj).fVt = fVt;
00253 ((THelix&)obj).fPhi0 = fPhi0;
00254 ((THelix&)obj).fVz = fVz;
00255 ((THelix&)obj).fW = fW;
00256 for (Int_t i=0; i<3; i++)
00257 ((THelix&)obj).fAxis[i] = fAxis[i];
00258
00259 if (((THelix&)obj).fRotMat)
00260 delete ((THelix&)obj).fRotMat;
00261 ((THelix&)obj).fRotMat = new TRotMatrix(*fRotMat);
00262
00263 ((THelix&)obj).fRange[0] = fRange[0];
00264 ((THelix&)obj).fRange[1] = fRange[1];
00265
00266 ((THelix&)obj).fOption = fOption;
00267
00268
00269
00270
00271 ((THelix&)obj).SetRange(fRange[0], fRange[1], kHelixT);
00272 }
00273
00274
00275
00276 void THelix::Draw(Option_t *option)
00277 {
00278
00279
00280 AppendPad(option);
00281 }
00282
00283
00284
00285 void THelix::Print(Option_t *option) const
00286 {
00287
00288
00289 cout <<" THelix Printing N=" <<fN<<" Option="<<option<<endl;
00290 }
00291
00292
00293
00294 void THelix::SavePrimitive(ostream &out, Option_t * )
00295 {
00296
00297
00298 char quote = '"';
00299 out<<" "<<endl;
00300 if (gROOT->ClassSaved(THelix::Class())) {
00301 out<<" ";
00302 } else {
00303 out<<" THelix *";
00304 }
00305 out<<"helix = new THelix("<<fX0<<","<<fY0<<","<<fZ0<<","
00306 <<fVt*TMath::Cos(fPhi0)<<","<<fVt*TMath::Sin(fPhi0)<<","<<fVz<<","
00307 <<fW<<","<<fRange[0]<<","<<fRange[1]<<","<<(Int_t)kHelixT<<","
00308 <<fAxis[0]<<","<<fAxis[1]<<","<<fAxis[2]<<","
00309 <<quote<<fOption<<quote<<");"<<endl;
00310
00311 SaveLineAttributes(out,"helix",1,1,1);
00312
00313 out<<" helix->Draw();"<<endl;
00314 }
00315
00316
00317
00318 void THelix::SetAxis(Double_t * axis)
00319 {
00320
00321
00322 if (axis) {
00323 Double_t len = TMath::Sqrt(axis[0]*axis[0] + axis[1]*axis[1] + axis[2]*axis[2]);
00324 if (len <= 0) {
00325 Error("SetAxis()", "Impossible! axis length %lf <= 0!", len);
00326 return;
00327 }
00328 fAxis[0] = axis[0]/len;
00329 fAxis[1] = axis[1]/len;
00330 fAxis[2] = axis[2]/len;
00331 } else {
00332 fAxis[0] = 0;
00333 fAxis[1] = 0;
00334 fAxis[2] = 1;
00335 }
00336
00337
00338 SetRotMatrix();
00339 }
00340
00341
00342
00343 void THelix::SetAxis(Double_t x, Double_t y, Double_t z)
00344 {
00345
00346
00347 Double_t axis[3]; axis[0] = x; axis[1] = y; axis[2] = z;
00348 SetAxis(axis);
00349 }
00350
00351
00352
00353 void THelix::SetRange(Double_t * range, EHelixRangeType rType)
00354 {
00355
00356
00357 Double_t a[2];
00358 Double_t halfpi = TMath::Pi()/2.0;
00359 Int_t i;
00360 Double_t vx = fVt * TMath::Cos(fPhi0);
00361 Double_t vy = fVt * TMath::Sin(fPhi0);
00362 Double_t phase;
00363
00364 if ( fW != 0 && fVz != 0 ) {
00365 switch ( rType ) {
00366 case kHelixT :
00367 fRange[0] = range[0]; fRange[1] = range[1]; break;
00368
00369 case kHelixX :
00370 for (i=0; i<2; i++ ) {
00371 a[i] = fW / fVt * (range[i] - fX0);
00372 if ( a[i] < -1 || a[i] > 1 ) {
00373 Error("SetRange()",
00374 "range out of bound (%lf:%lf): %lf. Default used: %lf",
00375 fX0-fVt/fW, fX0+fVt/fW, range[i], fRange[i]);
00376 return;
00377 }
00378 phase = FindClosestPhase(fPhi0+halfpi, a[i]);
00379 fRange[i] = ( fPhi0 + halfpi - phase ) / fW;
00380 }
00381 break;
00382
00383 case kHelixY :
00384 for (i=0; i<2; i++ ) {
00385 a[i] = fW / fVt * (range[i] - fY0);
00386 if ( a[i] < -1 || a[i] > 1 ) {
00387 Error("SetRange()",
00388 "range out of bound (%lf:%lf): %lf. Default used: %lf",
00389 fY0-fVt/fW, fY0+fVt/fW, range[i], fRange[i]);
00390 return;
00391 }
00392 phase = FindClosestPhase(fPhi0, a[i]);
00393 fRange[i] = ( fPhi0 - phase ) / fW;
00394 }
00395 break;
00396
00397 case kHelixZ :
00398 if ( fVz != 0 ) {
00399 for (i=0; i<2; i++ ) {
00400 fRange[i] = (range[i] - fZ0) / fVz;
00401 }
00402 } else {
00403 Error("SetRange()",
00404 "Vz = 0 and attempts to set range along helix axis!");
00405 return;
00406 }
00407 break;
00408
00409 case kLabX :
00410 case kLabY :
00411 case kLabZ :
00412 printf("setting range in lab axes is not implemented yet\n");
00413 break;
00414 default:
00415 Error("SetRange()","unknown range type %d", rType);
00416 break;
00417 }
00418 } else if ( fW == 0 ) {
00419 switch ( rType ) {
00420 case kHelixT :
00421 fRange[0] = range[0]; fRange[1] = range[1];
00422 break;
00423 case kHelixX :
00424 if ( vx != 0 ) {
00425 fRange[0] = (range[0] - fX0) / vx;
00426 fRange[1] = (range[1] - fX0) / vx;
00427 } else {
00428 Error("SetRange()",
00429 "Vx = 0 and attempts to set range on helix x axis!");
00430 return;
00431 }
00432 break;
00433 case kHelixY :
00434 if ( vy != 0 ) {
00435 fRange[0] = (range[0] - fY0) / vy;
00436 fRange[1] = (range[1] - fY0) / vy;
00437 } else {
00438 Error("SetRange()",
00439 "Vy = 0 and attempts to set range on helix y axis!");
00440 return;
00441 }
00442 break;
00443 case kHelixZ :
00444 if ( fVz != 0 ) {
00445 fRange[0] = (range[0] - fZ0) / fVz;
00446 fRange[1] = (range[1] - fZ0) / fVz;
00447 } else {
00448 Error("SetRange()",
00449 "Vz = 0 and attempts to set range on helix z axis!");
00450 return;
00451 }
00452 break;
00453 case kLabX :
00454 case kLabY :
00455 case kLabZ :
00456 printf("setting range in lab axes is not implemented yet\n");
00457 break;
00458 default :
00459 Error("SetRange()","unknown range type %d", rType);
00460 break;
00461 }
00462 } else if ( fVz == 0 ) {
00463 switch ( rType ) {
00464 case kHelixT :
00465 fRange[0] = range[0]; fRange[1] = range[1]; break;
00466 case kHelixX :
00467 if ( vx != 0 ) {
00468 fRange[0] = (range[0] - fX0) / vx;
00469 fRange[1] = (range[1] - fX0) / vx;
00470 } else {
00471 Error("SetRange()",
00472 "Vx = 0 and attempts to set range on helix x axis!");
00473 return;
00474 }
00475 break;
00476 case kHelixY :
00477 if ( vy != 0 ) {
00478 fRange[0] = (range[0] - fY0) / vy;
00479 fRange[1] = (range[1] - fY0) / vy;
00480 } else {
00481 Error("SetRange()",
00482 "Vy = 0 and attempts to set range on helix y axis!");
00483 return;
00484 }
00485 break;
00486 case kHelixZ :
00487 Error("SetRange()",
00488 "Vz = 0 and attempts to set range on helix z axis!");
00489 return;
00490 case kLabX :
00491 case kLabY :
00492 case kLabZ :
00493 printf("setting range in lab axes is not implemented yet\n");
00494 break;
00495 default :
00496 Error("SetRange()","unknown range type %d", rType);
00497 break;
00498 }
00499 }
00500
00501 if ( fRange[0] > fRange[1] ) {
00502 Double_t temp = fRange[1]; fRange[1] = fRange[0]; fRange[0] = temp;
00503 }
00504
00505
00506 Double_t degrad = TMath::Pi() / 180.0;
00507 Double_t segment = 5.0 * degrad;
00508 Double_t dt = segment / TMath::Abs(fW);
00509
00510 Int_t nSeg = Int_t((fRange[1]-fRange[0]) / dt) + 1;
00511 if (nSeg < THelix::fgMinNSeg) {
00512 nSeg = THelix::fgMinNSeg;
00513 dt = (fRange[1]-fRange[0])/nSeg;
00514 }
00515
00516 Double_t * xl = new Double_t[nSeg+1];
00517 Double_t * yl = new Double_t[nSeg+1];
00518 Double_t * zl = new Double_t[nSeg+1];
00519
00520 for (i=0; i<=nSeg; i++) {
00521 Double_t t, phase2;
00522 if (i==nSeg) t = fRange[1];
00523 else t = fRange[0] + dt * i;
00524 phase2 = -fW * t + fPhi0;
00525 xl[i] = fX0 - fVt/fW * TMath::Sin(phase2);
00526 yl[i] = fY0 + fVt/fW * TMath::Cos(phase2);
00527 zl[i] = fZ0 + fVz * t;
00528 }
00529
00530 Float_t xg, yg,zg;
00531
00532 Double_t * m = fRotMat->GetMatrix();
00533 TPolyLine3D::SetPolyLine(nSeg+1);
00534 for (i=0; i<=nSeg; i++) {
00535 xg = xl[i] * m[0] + yl[i] * m[3] + zl[i] * m[6] ;
00536 yg = xl[i] * m[1] + yl[i] * m[4] + zl[i] * m[7] ;
00537 zg = xl[i] * m[2] + yl[i] * m[5] + zl[i] * m[8] ;
00538 TPolyLine3D::SetPoint(i,xg,yg,zg);
00539 }
00540
00541 delete[] xl; delete[] yl; delete[] zl;
00542 }
00543
00544
00545
00546 void THelix::SetRange(Double_t r1, Double_t r2, EHelixRangeType rType)
00547 {
00548
00549
00550 Double_t range[2];
00551 range[0] = r1; range[1] = r2;
00552 SetRange(range, rType);
00553 }
00554
00555
00556
00557
00558
00559
00560
00561
00562
00563
00564 void THelix::SetRotMatrix()
00565 {
00566
00567
00568
00569
00570 Double_t raddeg = 180.0 / TMath::Pi();
00571 Double_t halfpi = TMath::Pi()/2.0 * raddeg;
00572
00573 Double_t theta3 = TMath::ACos(fAxis[2]) * raddeg;
00574 Double_t phi3 = TMath::ATan2(fAxis[1], fAxis[0]) * raddeg;
00575
00576 Double_t theta1 = theta3 + halfpi;
00577 Double_t phi1 = phi3;
00578
00579 Double_t theta2 = halfpi;
00580 Double_t phi2 = phi1 + halfpi;
00581
00582
00583 if (fRotMat) delete fRotMat;
00584
00585
00586 fRotMat = new TRotMatrix("HelixRotMat", "Master frame -> Helix frame",
00587 theta1, phi1, theta2, phi2, theta3, phi3 );
00588 return;
00589 }
00590
00591
00592
00593 Double_t THelix::FindClosestPhase(Double_t phi0, Double_t cosine)
00594 {
00595
00596
00597 const Double_t pi = TMath::Pi();
00598 const Double_t twopi = TMath::Pi() * 2.0;
00599 Double_t phi1 = TMath::ACos(cosine);
00600 Double_t phi2 = - phi1;
00601
00602 while ( phi1 - phi0 > pi ) phi1 -= twopi;
00603 while ( phi1 - phi0 < -pi ) phi1 += twopi;
00604
00605 while ( phi2 - phi0 > pi ) phi2 -= twopi;
00606 while ( phi2 - phi0 < -pi ) phi2 += twopi;
00607
00608
00609
00610 if ( TMath::Abs(phi1-phi0) < TMath::Abs(phi2-phi0) ) return phi1;
00611 else return phi2;
00612 }
00613
00614
00615
00616 void THelix::Streamer(TBuffer &R__b)
00617 {
00618
00619
00620 if (R__b.IsReading()) {
00621 UInt_t R__s, R__c;
00622 Version_t R__v = R__b.ReadVersion(&R__s, &R__c); if (R__v) { }
00623 if (R__v > 1) {
00624 R__b.ReadClassBuffer(THelix::Class(), this, R__v, R__s, R__c);
00625 return;
00626 }
00627
00628 TPolyLine3D::Streamer(R__b);
00629 R__b >> fX0;
00630 R__b >> fY0;
00631 R__b >> fZ0;
00632 R__b >> fVt;
00633 R__b >> fPhi0;
00634 R__b >> fVz;
00635 R__b >> fW;
00636 R__b.ReadStaticArray(fAxis);
00637 R__b >> fRotMat;
00638 R__b.ReadStaticArray(fRange);
00639 R__b.CheckByteCount(R__s, R__c, THelix::IsA());
00640
00641
00642 } else {
00643 R__b.WriteClassBuffer(THelix::Class(),this);
00644 }
00645 }