00001 #include <cmath>
00002 #include <iostream>
00003 #include "TopQuarkAnalysis/TopTools/interface/TopologyWorker.h"
00004
00005 Int_t TopologyWorker::m_maxpart = 1000;
00006
00007 TopologyWorker::TopologyWorker(bool boost):
00008 m_dSphMomPower(2.0),m_dDeltaThPower(0),
00009 m_iFast(4),m_dConv(0.0001),m_iGood(2)
00010 {
00011 m_dAxes.ResizeTo(4,4);
00012 m_mom.ResizeTo(m_maxpart,6);
00013 m_mom2.ResizeTo(m_maxpart,6);
00014 m_np=-1;
00015 m_np2=-1;
00016 m_sanda_called=false;
00017 m_fowo_called=false;
00018 m_sumangles_called=false;
00019 m_verbose=false;
00020 m_boost=boost;
00021 m_sph=-1;
00022 m_apl=-1;
00023 m_h10=-1;
00024 m_h20=-1;
00025 m_h30=-1;
00026 m_h40=-1;
00027 m_sqrts=0;
00028 m_ht=0;
00029 m_ht3=0;
00030 m_et56=0;
00031 m_njetsweighed=0;
00032 m_et0=0;
00033 }
00034
00035
00036
00037
00038
00039
00040
00041
00042 void TopologyWorker::setPartList(TObjArray* e1, TObjArray* e2)
00043 {
00044
00045
00046
00047 TMatrix mom(m_maxpart,6);
00048 TMatrix mom2(m_maxpart,6);
00049 double tmax = 0;
00050 double phi = 0.;
00051 double the = 0.;
00052 double sgn;
00053 TMatrix fast(m_iFast + 1,6);
00054 TMatrix work(11,6);
00055 double tdi[4] = {0.,0.,0.,0.};
00056 double tds;
00057 double tpr[4] = {0.,0.,0.,0.};
00058 double thp;
00059 double thps;
00060 double pxtot=0;
00061 double pytot=0;
00062 double pztot=0;
00063 double etot=0;
00064
00065 TMatrix temp(3,5);
00066 Int_t np = 0;
00067 Int_t numElements = e1->GetEntries();
00068 Int_t numElements2 = e2->GetEntries();
00069
00070
00071
00072
00073
00074 m_np=0;
00075 for(Int_t elem=0;elem<numElements;elem++) {
00076 if(m_verbose){
00077 std::cout << "looping over array, element " << elem << std::endl;
00078 }
00079 TObject* o = (TObject*) e1->At(elem);
00080 if(m_verbose){
00081 std::cerr << "TopologyWorker:SetPartList(): adding jet " << elem << "." << std::endl;
00082 }
00083 if (np >= m_maxpart) {
00084 printf("Too many particles input to TopologyWorker");
00085 return;
00086 }
00087 if(m_verbose){
00088 std::cout << "getting name of object..." << std::endl;
00089 }
00090 TString nam(o->IsA()->GetName());
00091 if(m_verbose){
00092 std::cout << "TopologyWorker::setPartList(): object is of type " << nam << std::endl;
00093 }
00094 if (nam.Contains("TVector3")) {
00095 TVector3 v(((TVector3 *) o)->X(),
00096 ((TVector3 *) o)->Y(),
00097 ((TVector3 *) o)->Z());
00098 mom(np,1) = v.X();
00099 mom(np,2) = v.Y();
00100 mom(np,3) = v.Z();
00101 mom(np,4) = v.Mag();
00102 }
00103 else if (nam.Contains("TLorentzVector")) {
00104 TVector3 v(((TLorentzVector *) o)->X(),
00105 ((TLorentzVector *) o)->Y(),
00106 ((TLorentzVector *) o)->Z());
00107 mom(np,1) = v.X();
00108 mom(np,2) = v.Y();
00109 mom(np,3) = v.Z();
00110 mom(np,4) = ((TLorentzVector *) o)->T();
00111 }
00112 else {
00113 printf("TopologyWorker::setEvent input is not a TVector3 or a TLorentzVector\n");
00114 continue;
00115 }
00116
00117
00118 if ( TMath::Abs( m_dDeltaThPower ) <= 0.001 ) {
00119 mom(np,5) = 1.0;
00120 }
00121 else {
00122 mom(np,5) = TMath::Power(mom(np,4),m_dDeltaThPower);
00123 }
00124 tmax = tmax + mom(np,4)*mom(np,5);
00125 pxtot+=mom(np,1);
00126 pytot+=mom(np,2);
00127 pztot+=mom(np,3);
00128 etot+=mom(np,4);
00129 np++;
00130 m_np=np;
00131 }
00132 Int_t np2=0;
00133
00134 for(Int_t elem=0;elem<numElements2;elem++) {
00135
00136 TObject* o = (TObject*) e2->At(elem);
00137 if (np2 >= m_maxpart) {
00138 printf("Too many particles input to TopologyWorker");
00139 return;
00140 }
00141
00142 TString nam(o->IsA()->GetName());
00143 if (nam.Contains("TVector3")) {
00144 TVector3 v(((TVector3 *) o)->X(),
00145 ((TVector3 *) o)->Y(),
00146 ((TVector3 *) o)->Z());
00147 mom2(np2,1) = v.X();
00148 mom2(np2,2) = v.Y();
00149 mom2(np2,3) = v.Z();
00150 mom2(np2,4) = v.Mag();
00151 }
00152 else if (nam.Contains("TLorentzVector")) {
00153 TVector3 v(((TLorentzVector *) o)->X(),
00154 ((TLorentzVector *) o)->Y(),
00155 ((TLorentzVector *) o)->Z());
00156 mom2(np2,1) = v.X();
00157 mom2(np2,2) = v.Y();
00158 mom2(np2,3) = v.Z();
00159 mom2(np2,4) = ((TLorentzVector *) o)->T();
00160
00161 }
00162 else {
00163 printf("TopologyWorker::setEvent Second vector input is not a TVector3 or a TLorentzVector\n");
00164 continue;
00165 }
00166 np2++;
00167 m_np2=np2;
00168 }
00169 m_mom2=mom2;
00170
00171 if (m_boost && m_np>1) {
00172 printf("TopologyWorker::setEvent Only boosting first vector so watch out when you do this!!!");
00173 TVector3 booze(-pxtot/etot,-pytot/etot,-pztot/etot);
00174 TLorentzVector l1;
00175 for (int k=0; k<m_np ; k++) {
00176 l1.SetPx(mom(k,1));
00177 l1.SetPy(mom(k,2));
00178 l1.SetPz(mom(k,3));
00179 l1.SetE(mom(k,4));
00180 l1.Boost(booze);
00181 mom(k,1)=l1.Px();
00182 mom(k,2)=l1.Py();
00183 mom(k,3)=l1.Pz();
00184 mom(k,4)=l1.E();
00185 }
00186 }
00187
00188 m_sanda_called=false;
00189 m_fowo_called=false;
00190 for (int ip=0;ip<m_maxpart;ip++) {
00191 for (int id=0;id<6;id++) {
00192 m_mom(ip,id)=mom(ip,id);
00193 }
00194 }
00195
00196 if ( np < 2 ) {
00197 m_dThrust[1] = -1.0;
00198 m_dOblateness = -1.0;
00199 return;
00200 }
00201
00202
00203 for ( Int_t pass=1; pass < 3; pass++) {
00204 if ( pass == 2 ) {
00205 phi = ulAngle(m_dAxes(1,1), m_dAxes(1,2));
00206 ludbrb( &mom, 0, -phi, 0., 0., 0. );
00207 for ( Int_t i = 0; i < 3; i++ ) {
00208 for ( Int_t j = 1; j < 4; j++ ) {
00209 temp(i,j) = m_dAxes(i+1,j);
00210 }
00211 temp(i,4) = 0;
00212 }
00213 ludbrb(&temp,0.,-phi,0.,0.,0.);
00214 for ( Int_t ib = 0; ib < 3; ib++ ) {
00215 for ( Int_t j = 1; j < 4; j++ ) {
00216 m_dAxes(ib+1,j) = temp(ib,j);
00217 }
00218 }
00219 the = ulAngle( m_dAxes(1,3), m_dAxes(1,1) );
00220 ludbrb( &mom, -the, 0., 0., 0., 0. );
00221 for ( Int_t ic = 0; ic < 3; ic++ ) {
00222 for ( Int_t j = 1; j < 4; j++ ) {
00223 temp(ic,j) = m_dAxes(ic+1,j);
00224 }
00225 temp(ic,4) = 0;
00226 }
00227 ludbrb(&temp,-the,0.,0.,0.,0.);
00228 for ( Int_t id = 0; id < 3; id++ ) {
00229 for ( Int_t j = 1; j < 4; j++ ) {
00230 m_dAxes(id+1,j) = temp(id,j);
00231 }
00232 }
00233 }
00234 for ( Int_t ifas = 0; ifas < m_iFast + 1 ; ifas++ ) {
00235 fast(ifas,4) = 0.;
00236 }
00237
00238
00239
00240 for ( Int_t i = 0; i < np; i++ ) {
00241 if ( pass == 2 ) {
00242 mom(i,4) = TMath::Sqrt( mom(i,1)*mom(i,1)
00243 + mom(i,2)*mom(i,2) );
00244 }
00245 for ( Int_t ifas = m_iFast - 1; ifas > -1; ifas-- ) {
00246 if ( mom(i,4) > fast(ifas,4) ) {
00247 for ( Int_t j = 1; j < 6; j++ ) {
00248 fast(ifas+1,j) = fast(ifas,j);
00249 if ( ifas == 0 ) fast(ifas,j) = mom(i,j);
00250 }
00251 }
00252 else {
00253 for ( Int_t j = 1; j < 6; j++ ) {
00254 fast(ifas+1,j) = mom(i,j);
00255 }
00256 break;
00257 }
00258 }
00259 }
00260
00261 for ( Int_t ie = 0; ie < work.GetNrows(); ie++ ) {
00262 work(ie,4) = 0.;
00263 }
00264 Int_t p = TMath::Min( m_iFast, np ) - 1;
00265
00266
00267 Int_t nc = iPow(2,p);
00268 for ( Int_t n = 0; n < nc; n++ ) {
00269 for ( Int_t j = 1; j < 4; j++ ) {
00270 tdi[j] = 0.;
00271 }
00272 for ( Int_t i = 0; i < TMath::Min(m_iFast,n); i++ ) {
00273 sgn = fast(i,5);
00274 if (iPow(2,(i+1))*((n+iPow(2,i))/iPow(2,(i+1))) >= i+1){
00275 sgn = -sgn;
00276 }
00277 for ( Int_t j = 1; j < 5-pass; j++ ) {
00278 tdi[j] = tdi[j] + sgn*fast(i,j);
00279 }
00280 }
00281 tds = tdi[1]*tdi[1] + tdi[2]*tdi[2] + tdi[3]*tdi[3];
00282 for ( Int_t iw = TMath::Min(n,9); iw > -1; iw--) {
00283 if( tds > work(iw,4) ) {
00284 for ( Int_t j = 1; j < 5; j++ ) {
00285 work(iw+1,j) = work(iw,j);
00286 if ( iw == 0 ) {
00287 if ( j < 4 ) {
00288 work(iw,j) = tdi[j];
00289 }
00290 else {
00291 work(iw,j) = tds;
00292 }
00293 }
00294 }
00295 }
00296 else {
00297 for ( Int_t j = 1; j < 4; j++ ) {
00298 work(iw+1,j) = tdi[j];
00299 }
00300 work(iw+1,4) = tds;
00301 }
00302 }
00303 }
00304
00305 m_dThrust[pass] = 0;
00306 thp = -99999.;
00307 Int_t nagree = 0;
00308 for ( Int_t iw = 0;
00309 iw < TMath::Min(nc,10) && nagree < m_iGood; iw++ ){
00310 thp = 0.;
00311 thps = -99999.;
00312 while ( thp > thps + m_dConv ) {
00313 thps = thp;
00314 for ( Int_t j = 1; j < 4; j++ ) {
00315 if ( thp <= 1E-10 ) {
00316 tdi[j] = work(iw,j);
00317 }
00318 else {
00319 tdi[j] = tpr[j];
00320 tpr[j] = 0;
00321 }
00322 }
00323 for ( Int_t i = 0; i < np; i++ ) {
00324 sgn = sign(mom(i,5),
00325 tdi[1]*mom(i,1) +
00326 tdi[2]*mom(i,2) +
00327 tdi[3]*mom(i,3));
00328 for ( Int_t j = 1; j < 5 - pass; j++ ){
00329 tpr[j] = tpr[j] + sgn*mom(i,j);
00330 }
00331 }
00332 thp = TMath::Sqrt(tpr[1]*tpr[1]
00333 + tpr[2]*tpr[2]
00334 + tpr[3]*tpr[3])/tmax;
00335 }
00336
00337
00338 if ( thp < m_dThrust[pass] - m_dConv ) {
00339 break;
00340 }
00341 if ( thp > m_dThrust[pass] + m_dConv ) {
00342 nagree = 0;
00343 sgn = iPow( -1, (Int_t)TMath::Nint(m_random.Rndm()) );
00344 for ( Int_t j = 1; j < 4; j++ ) {
00345 m_dAxes(pass,j) = sgn*tpr[j]/(tmax*thp);
00346 }
00347 m_dThrust[pass] = thp;
00348 }
00349 nagree = nagree + 1;
00350 }
00351 }
00352
00353 sgn = iPow( -1, (Int_t)TMath::Nint(m_random.Rndm()));
00354 m_dAxes(3,1) = -sgn*m_dAxes(2,2);
00355 m_dAxes(3,2) = sgn*m_dAxes(2,1);
00356 m_dAxes(3,3) = 0.;
00357 thp = 0.;
00358 for ( Int_t i = 0; i < np; i++ ) {
00359 thp += mom(i,5)*TMath::Abs(m_dAxes(3,1)*mom(i,1) +
00360 m_dAxes(3,2)*mom(i,2));
00361 }
00362 m_dThrust[3] = thp/tmax;
00363
00364 for ( Int_t i6 = 0; i6 < 3; i6++ ) {
00365 for ( Int_t j = 1; j < 4; j++ ) {
00366 temp(i6,j) = m_dAxes(i6+1,j);
00367 }
00368 temp(i6,4) = 0;
00369 }
00370 ludbrb(&temp,the,phi,0.,0.,0.);
00371 for ( Int_t i7 = 0; i7 < 3; i7++ ) {
00372 for ( Int_t j = 1; j < 4; j++ ) {
00373 m_dAxes(i7+1,j) = temp(i7,j);
00374 }
00375 }
00376 m_dOblateness = m_dThrust[2] - m_dThrust[3];
00377
00378
00379
00380
00381 CalcWmul();
00382 CalcHTstuff();
00383 CalcSqrts();
00384
00385 }
00386
00387
00388
00389
00390 void TopologyWorker::setThMomPower(double tp)
00391 {
00392
00393 if ( tp > 0. ) m_dDeltaThPower = tp - 1.0;
00394 return;
00395 }
00396
00397
00398 double TopologyWorker::getThMomPower()
00399 {
00400 return 1.0 + m_dDeltaThPower;
00401 }
00402
00403
00404 void TopologyWorker::setFast(Int_t nf)
00405 {
00406
00407 if ( nf > 3 ) m_iFast = nf;
00408 return;
00409 }
00410
00411
00412 Int_t TopologyWorker::getFast()
00413 {
00414 return m_iFast;
00415 }
00416
00417
00418
00419
00420 TVector3 TopologyWorker::thrustAxis() {
00421 TVector3 m_ThrustAxis(m_dAxes(1,1),m_dAxes(1,2),m_dAxes(1,3));
00422 return m_ThrustAxis;
00423 }
00424
00425
00426 TVector3 TopologyWorker::majorAxis() {
00427 TVector3 m_MajorAxis(m_dAxes(2,1),m_dAxes(2,2),m_dAxes(2,3));
00428 return m_MajorAxis;
00429 }
00430
00431
00432 TVector3 TopologyWorker::minorAxis() {
00433 TVector3 m_MinorAxis(m_dAxes(3,1),m_dAxes(3,2),m_dAxes(3,3));
00434 return m_MinorAxis;
00435 }
00436
00437
00438 TVector3 TopologyWorker::thrust() {
00439 TVector3 m_Thrust(m_dThrust[1],m_dThrust[2],m_dThrust[3]);
00440 return m_Thrust;
00441 }
00442
00443
00444 double TopologyWorker::oblateness() {
00445 return m_dOblateness;
00446 }
00447
00448
00449
00450 double TopologyWorker::ulAngle(double x, double y)
00451 {
00452 double ulangl = 0;
00453 double r = TMath::Sqrt(x*x + y*y);
00454 if ( r < 1.0E-20 ) {
00455 return ulangl;
00456 }
00457 if ( TMath::Abs(x)/r < 0.8 ) {
00458 ulangl = sign(TMath::ACos(x/r),y);
00459 }
00460 else {
00461 ulangl = TMath::ASin(y/r);
00462 if ( x < 0. && ulangl >= 0. ) {
00463 ulangl = TMath::Pi() - ulangl;
00464 }
00465 else if ( x < 0. ) {
00466 ulangl = -TMath::Pi() - ulangl;
00467 }
00468 }
00469 return ulangl;
00470 }
00471
00472
00473 double TopologyWorker::sign(double a, double b) {
00474 if ( b < 0 ) {
00475 return -TMath::Abs(a);
00476 }
00477 else {
00478 return TMath::Abs(a);
00479 }
00480 }
00481
00482
00483 void TopologyWorker::ludbrb(TMatrix* mom,
00484 double the,
00485 double phi,
00486 double bx,
00487 double by,
00488 double bz)
00489 {
00490
00491
00492 TMatrix rot(4,4);
00493 double pr[4];
00494 double dp[5];
00495 Int_t np = mom->GetNrows();
00496 if ( the*the + phi*phi > 1.0E-20 )
00497 {
00498 rot(1,1) = TMath::Cos(the)*TMath::Cos(phi);
00499 rot(1,2) = -TMath::Sin(phi);
00500 rot(1,3) = TMath::Sin(the)*TMath::Cos(phi);
00501 rot(2,1) = TMath::Cos(the)*TMath::Sin(phi);
00502 rot(2,2) = TMath::Cos(phi);
00503 rot(2,3) = TMath::Sin(the)*TMath::Sin(phi);
00504 rot(3,1) = -TMath::Sin(the);
00505 rot(3,2) = 0.0;
00506 rot(3,3) = TMath::Cos(the);
00507 for ( Int_t i = 0; i < np; i++ )
00508 {
00509 for ( Int_t j = 1; j < 4; j++ )
00510 {
00511 pr[j] = (*mom)(i,j);
00512 (*mom)(i,j) = 0;
00513 }
00514 for ( Int_t jb = 1; jb < 4; jb++)
00515 {
00516 for ( Int_t k = 1; k < 4; k++)
00517 {
00518 (*mom)(i,jb) = (*mom)(i,jb) + rot(jb,k)*pr[k];
00519 }
00520 }
00521 }
00522 double beta = TMath::Sqrt( bx*bx + by*by + bz*bz );
00523 if ( beta*beta > 1.0E-20 )
00524 {
00525 if ( beta > 0.99999999 )
00526 {
00527
00528 bx = bx*(0.99999999/beta);
00529 by = by*(0.99999999/beta);
00530 bz = bz*(0.99999999/beta);
00531 beta = 0.99999999;
00532 }
00533 double gamma = 1.0/TMath::Sqrt(1.0 - beta*beta);
00534 for ( Int_t i = 0; i < np; i++ )
00535 {
00536 for ( Int_t j = 1; j < 5; j++ )
00537 {
00538 dp[j] = (*mom)(i,j);
00539 }
00540 double bp = bx*dp[1] + by*dp[2] + bz*dp[3];
00541 double gbp = gamma*(gamma*bp/(1.0 + gamma) + dp[4]);
00542 (*mom)(i,1) = dp[1] + gbp*bx;
00543 (*mom)(i,2) = dp[2] + gbp*by;
00544 (*mom)(i,3) = dp[3] + gbp*bz;
00545 (*mom)(i,4) = gamma*(dp[4] + bp);
00546 }
00547 }
00548 }
00549 return;
00550 }
00551
00552
00553
00554
00555
00556 void TopologyWorker::sanda() {
00557 float SPH=-1;
00558 float APL=-1;
00559 m_sanda_called=true;
00560
00561
00562
00563
00564
00565
00566
00567
00568
00569
00570
00571
00572
00573
00574
00575 float P[1000][6];
00576 double SM[4][4],SV[4][4];
00577 double PA,PWT,PS,SQ,SR,SP,SMAX,SGN;
00578 int NP;
00579 int J, J1, J2, I, N, JA, JB, J3, JC, JB1, JB2;
00580 JA=JB=JC=0;
00581 double RL;
00582 float rlu,rlu1;
00583
00584
00585
00586 N=m_np;
00587 NP=0;
00588 for (I=1;I<N+1;I++){
00589 NP++;
00590 P[NP][1]=m_mom(I-1,1) ;
00591 P[NP][2]=m_mom(I-1,2) ;
00592 P[NP][3]=m_mom(I-1,3) ;
00593 P[NP][4]=m_mom(I-1,4) ;
00594 P[NP][5]=0;
00595 }
00596
00597
00598
00599 N=NP;
00600
00601 for (J1=1;J1<4;J1++) {
00602 for (J2=J1;J2<4;J2++) {
00603 SM[J1][J2]=0.;
00604 }
00605 }
00606 PS=0.;
00607 for (I=1;I<N+1;I++) {
00608 PA=sqrt(pow(P[I][1],2)+pow(P[I][2],2)+pow(P[I][3],2));
00609 PWT=1.;
00610 for (J1=1;J1<4;J1++) {
00611 for (J2=J1;J2<4;J2++) {
00612 SM[J1][J2]=SM[J1][J2]+PWT*P[I][J1]*P[I][J2];
00613 }
00614 }
00615 PS=PS+PWT*PA*PA;
00616 }
00617
00618 if(NP<2) {
00619 SPH=-1.;
00620 APL=-1.;
00621 return;
00622 }
00623 for (J1=1;J1<4;J1++) {
00624 for (J2=J1;J2<4;J2++) {
00625 SM[J1][J2]=SM[J1][J2]/PS;
00626 }
00627 }
00628
00629 SQ=(SM[1][1]*SM[2][2]+SM[1][1]*SM[3][3]+SM[2][2]*SM[3][3]
00630 -pow(SM[1][2],2)
00631 -pow(SM[1][3],2)-pow(SM[2][3],2))/3.-1./9.;
00632 SR=-0.5*(SQ+1./9.+SM[1][1]*pow(SM[2][3],2)+SM[2][2]*pow(SM[1][3],2)+SM[3][3]*
00633 pow(SM[1][2],2)-SM[1][1]*SM[2][2]*SM[3][3])+SM[1][2]*SM[1][3]*SM[2][3]+1./27.;
00634
00635 SP=TMath::Cos(TMath::ACos(TMath::Max(TMath::Min(SR/TMath::Sqrt(-pow(SQ,3)),1.),-1.))/3.);
00636
00637 P[N+1][4]=1./3.+TMath::Sqrt(-SQ)*TMath::Max(2.*SP,TMath::Sqrt(3.*(1.-SP*SP))-SP);
00638 P[N+3][4]=1./3.+TMath::Sqrt(-SQ)*TMath::Min(2.*SP,-TMath::Sqrt(3.*(1.-SP*SP))-SP);
00639 P[N+2][4]=1.-P[N+1][4]-P[N+3][4];
00640 if (P[N+2][4]> 1E-5) {
00641
00642
00643 for (I=1;I<4;I=I+2) {
00644 for (J1=1;J1<4;J1++) {
00645 SV[J1][J1]=SM[J1][J1]-P[N+I][4];
00646 for (J2=J1+1;J2<4;J2++) {
00647 SV[J1][J2]=SM[J1][J2];
00648 SV[J2][J1]=SM[J1][J2];
00649 }
00650 }
00651 SMAX=0.;
00652 for (J1=1;J1<4;J1++) {
00653 for (J2=1;J2<4;J2++) {
00654 if(std::fabs(SV[J1][J2])>SMAX) {
00655 JA=J1;
00656 JB=J2;
00657 SMAX=std::fabs(SV[J1][J2]);
00658 }
00659 }
00660 }
00661 SMAX=0.;
00662 for (J3=JA+1;J3<JA+3;J3++) {
00663 J1=J3-3*((J3-1)/3);
00664 RL=SV[J1][JB]/SV[JA][JB];
00665 for (J2=1;J2<4;J2++) {
00666 SV[J1][J2]=SV[J1][J2]-RL*SV[JA][J2];
00667 if (std::fabs(SV[J1][J2])>SMAX) {
00668 JC=J1;
00669 SMAX=std::fabs(SV[J1][J2]);
00670 }
00671 }
00672 }
00673 JB1=JB+1-3*(JB/3);
00674 JB2=JB+2-3*((JB+1)/3);
00675 P[N+I][JB1]=-SV[JC][JB2];
00676 P[N+I][JB2]=SV[JC][JB1];
00677 P[N+I][JB]=-(SV[JA][JB1]*P[N+I][JB1]+SV[JA][JB2]*P[N+I][JB2])/
00678 SV[JA][JB];
00679 PA=TMath::Sqrt(pow(P[N+I][1],2)+pow(P[N+I][2],2)+pow(P[N+I][3],2));
00680
00681 float pa=P[N-1][I];
00682 rlu=std::fabs(pa)-std::fabs(int(pa)*1.);
00683 rlu1=std::fabs(pa*pa)-std::fabs(int(pa*pa)*1.);
00684 SGN=pow((-1.),1.*int(rlu+0.5));
00685 for (J=1;J<4;J++) {
00686 P[N+I][J]=SGN*P[N+I][J]/PA;
00687 }
00688 }
00689
00690
00691 SGN=pow((-1.),1.*int(rlu1+0.5));
00692 P[N+2][1]=SGN*(P[N+1][2]*P[N+3][3]-P[N+1][3]*P[N+3][2]);
00693 P[N+2][2]=SGN*(P[N+1][3]*P[N+3][1]-P[N+1][1]*P[N+3][3]);
00694 P[N+2][3]=SGN*(P[N+1][1]*P[N+3][2]-P[N+1][2]*P[N+3][1]);
00695
00696
00697 SPH=1.5*(P[N+2][4]+P[N+3][4]);
00698 APL=1.5*P[N+3][4];
00699
00700 }
00701
00702 m_sph=SPH;
00703 m_apl=APL;
00704 return;
00705 }
00706
00707
00708
00709
00710 void TopologyWorker::planes_sphe(double& pnorm, double& p2, double& p3) {
00711 float SPH=-1;
00712 float APL=-1;
00713
00714 float P[1000][6];
00715 double SM[4][4],SV[4][4];
00716 double PA,PWT,PS,SQ,SR,SP,SMAX,SGN;
00717 int NP;
00718 int J, J1, J2, I, N, JA, JB, J3, JC, JB1, JB2;
00719 JA=JB=JC=0;
00720 double RL;
00721 float rlu,rlu1;
00722
00723
00724
00725 N=m_np;
00726 NP=0;
00727 for (I=1;I<N+1;I++){
00728 NP++;
00729 P[NP][1]=m_mom(I-1,1) ;
00730 P[NP][2]=m_mom(I-1,2) ;
00731 P[NP][3]=m_mom(I-1,3) ;
00732 P[NP][4]=m_mom(I-1,4) ;
00733 P[NP][5]=0;
00734 }
00735
00736
00737
00738 N=NP;
00739
00740 for (J1=1;J1<4;J1++) {
00741 for (J2=J1;J2<4;J2++) {
00742 SM[J1][J2]=0.;
00743 }
00744 }
00745 PS=0.;
00746 for (I=1;I<N+1;I++) {
00747 PA=sqrt(pow(P[I][1],2)+pow(P[I][2],2)+pow(P[I][3],2));
00748 double eta,phi;
00749 getetaphi(P[I][1],P[I][2],P[I][3],eta,phi);
00750 PWT=exp(-std::fabs(eta));
00751 PWT=1.;
00752 for (J1=1;J1<4;J1++) {
00753 for (J2=J1;J2<4;J2++) {
00754 SM[J1][J2]=SM[J1][J2]+PWT*P[I][J1]*P[I][J2];
00755 }
00756 }
00757 PS=PS+PWT*PA*PA;
00758 }
00759
00760 if(NP<2) {
00761 SPH=-1.;
00762 APL=-1.;
00763 return;
00764 }
00765 for (J1=1;J1<4;J1++) {
00766 for (J2=J1;J2<4;J2++) {
00767 SM[J1][J2]=SM[J1][J2]/PS;
00768 }
00769 }
00770
00771 SQ=(SM[1][1]*SM[2][2]+SM[1][1]*SM[3][3]+SM[2][2]*SM[3][3]
00772 -pow(SM[1][2],2)
00773 -pow(SM[1][3],2)-pow(SM[2][3],2))/3.-1./9.;
00774 SR=-0.5*(SQ+1./9.+SM[1][1]*pow(SM[2][3],2)+SM[2][2]*pow(SM[1][3],2)+SM[3][3]*
00775 pow(SM[1][2],2)-SM[1][1]*SM[2][2]*SM[3][3])+SM[1][2]*SM[1][3]*SM[2][3]+1./27.;
00776
00777 SP=TMath::Cos(TMath::ACos(TMath::Max(TMath::Min(SR/TMath::Sqrt(-pow(SQ,3)),1.),-1.))/3.);
00778
00779 P[N+1][4]=1./3.+TMath::Sqrt(-SQ)*TMath::Max(2.*SP,TMath::Sqrt(3.*(1.-SP*SP))-SP);
00780 P[N+3][4]=1./3.+TMath::Sqrt(-SQ)*TMath::Min(2.*SP,-TMath::Sqrt(3.*(1.-SP*SP))-SP);
00781 P[N+2][4]=1.-P[N+1][4]-P[N+3][4];
00782 if (P[N+2][4]> 1E-5) {
00783
00784
00785 for (I=1;I<4;I=I+2) {
00786 for (J1=1;J1<4;J1++) {
00787 SV[J1][J1]=SM[J1][J1]-P[N+I][4];
00788 for (J2=J1+1;J2<4;J2++) {
00789 SV[J1][J2]=SM[J1][J2];
00790 SV[J2][J1]=SM[J1][J2];
00791 }
00792 }
00793 SMAX=0.;
00794 for (J1=1;J1<4;J1++) {
00795 for (J2=1;J2<4;J2++) {
00796 if(std::fabs(SV[J1][J2])>SMAX) {
00797 JA=J1;
00798 JB=J2;
00799 SMAX=std::fabs(SV[J1][J2]);
00800 }
00801 }
00802 }
00803 SMAX=0.;
00804 for (J3=JA+1;J3<JA+3;J3++) {
00805 J1=J3-3*((J3-1)/3);
00806 RL=SV[J1][JB]/SV[JA][JB];
00807 for (J2=1;J2<4;J2++) {
00808 SV[J1][J2]=SV[J1][J2]-RL*SV[JA][J2];
00809 if (std::fabs(SV[J1][J2])>SMAX) {
00810 JC=J1;
00811 SMAX=std::fabs(SV[J1][J2]);
00812 }
00813 }
00814 }
00815 JB1=JB+1-3*(JB/3);
00816 JB2=JB+2-3*((JB+1)/3);
00817 P[N+I][JB1]=-SV[JC][JB2];
00818 P[N+I][JB2]=SV[JC][JB1];
00819 P[N+I][JB]=-(SV[JA][JB1]*P[N+I][JB1]+SV[JA][JB2]*P[N+I][JB2])/
00820 SV[JA][JB];
00821 PA=TMath::Sqrt(pow(P[N+I][1],2)+pow(P[N+I][2],2)+pow(P[N+I][3],2));
00822
00823 float pa=P[N-1][I];
00824 rlu=std::fabs(pa)-std::fabs(int(pa)*1.);
00825 rlu1=std::fabs(pa*pa)-std::fabs(int(pa*pa)*1.);
00826 SGN=pow((-1.),1.*int(rlu+0.5));
00827 for (J=1;J<4;J++) {
00828 P[N+I][J]=SGN*P[N+I][J]/PA;
00829 }
00830 }
00831
00832
00833 SGN=pow((-1.),1.*int(rlu1+0.5));
00834 P[N+2][1]=SGN*(P[N+1][2]*P[N+3][3]-P[N+1][3]*P[N+3][2]);
00835 P[N+2][2]=SGN*(P[N+1][3]*P[N+3][1]-P[N+1][1]*P[N+3][3]);
00836 P[N+2][3]=SGN*(P[N+1][1]*P[N+3][2]-P[N+1][2]*P[N+3][1]);
00837
00838
00839 SPH=1.5*(P[N+2][4]+P[N+3][4]);
00840 APL=1.5*P[N+3][4];
00841
00842 }
00843
00844
00845 double etstot[4];
00846 double eltot[4];
00847 double sum23=0;
00848 double sum22=0;
00849 double sum33=0;
00850 double pina[4];
00851 double ax[4], ay[4], az[4];
00852 for (int ia=1;ia<4;ia++) {
00853 etstot[ia]=0;
00854 eltot[ia]=0;
00855 pina[ia]=0;
00856 ax[ia]=P[N+ia][1];
00857 ay[ia]=P[N+ia][2];
00858 az[ia]=P[N+ia][3];
00859 ax[ia]/=sqrt(ax[ia]*ax[ia]+ay[ia]*ay[ia]+az[ia]*az[ia]);
00860 ay[ia]/=sqrt(ax[ia]*ax[ia]+ay[ia]*ay[ia]+az[ia]*az[ia]);
00861 az[ia]/=sqrt(ax[ia]*ax[ia]+ay[ia]*ay[ia]+az[ia]*az[ia]);
00862 }
00863
00864
00865 for (int k =0 ; k<m_np ; k++) {
00866
00867
00868 double W=1.0;
00869 for (int ia=1;ia<4;ia++) {
00870 double e=sqrt(m_mom(k,1)*m_mom(k,1) +
00871 m_mom(k,2)*m_mom(k,2) +
00872 m_mom(k,3)*m_mom(k,3));
00873 double el=ax[ia]*m_mom(k,1) + ay[ia]*m_mom(k,2) + az[ia]*m_mom(k,3);
00874 pina[ia]=el;
00875 double ets=(e*e-el*el);
00876 etstot[ia]+=ets*W;
00877 eltot[ia]+=el*el;
00878 }
00879 double a2=pina[2];
00880 double a3=pina[3];
00881
00882
00883
00884 sum22+=a2*a2*W;
00885 sum23+=a2*a3*W;
00886 sum33+=a3*a3*W;
00887 }
00888
00889
00890
00891 double pi=3.1415927;
00892 double phi=pi/2.0;
00893 double phip=pi/2.0;
00894 double a=sum23;
00895 double c=-a;
00896 double b=sum22-sum33;
00897 double disc=b*b-4*a*c;
00898
00899 if (disc>=0) {
00900 double x1=(sqrt(disc)-b)/2/a;
00901 double x2=(-sqrt(disc)-b)/2/a;
00902 phi=atan(x1);
00903 phip=atan(x2);
00904 if (phi<0) phi=2.*pi+phi;
00905 if (phip<0) phip=2.*pi+phip;
00906 }
00907 double p21=sum22*cos(phi)*cos(phi)+sum33*sin(phi)*sin(phi)+2*sum23*cos(phi)*sin(phi);
00908 double p31=sum22*sin(phi)*sin(phi)+sum33*cos(phi)*cos(phi)-2*sum23*cos(phi)*sin(phi);
00909
00910 double p22=sum22*cos(phip)*cos(phip)+sum33*sin(phip)*sin(phip)+2*sum23*cos(phip)*sin(phip);
00911 double p32=sum22*sin(phip)*sin(phip)+sum33*cos(phip)*cos(phip)-2*sum23*cos(phip)*sin(phip);
00912
00913
00914 double d1=std::fabs(p31*p31 - p21*p21);
00915 double d2=std::fabs(p32*p32 - p22*p22);
00916
00917
00918
00919 p2=p21;
00920 p3=p31;
00921 if (d2>d1) {
00922 p2=p22;
00923 p3=p32;
00924 }
00925 pnorm=sqrt(eltot[1]);
00926 if (p2>p3) {
00927 p3=sqrt(p3);
00928 p2=sqrt(p2);
00929 }else {
00930 double p4=p3;
00931 p3=sqrt(p2);
00932 p2=sqrt(p4);
00933 }
00934
00935
00936
00937
00938 return;
00939 }
00940
00941
00942 void TopologyWorker::planes_sphe_wei(double& pnorm, double& p2, double& p3) {
00943 float SPH=-1;
00944 float APL=-1;
00945
00946 float P[1000][6];
00947 double SM[4][4],SV[4][4];
00948 double PA,PWT,PS,SQ,SR,SP,SMAX,SGN;
00949 int NP;
00950 int J, J1, J2, I, N, JA, JB, J3, JC, JB1, JB2;
00951 JA=JB=JC=0;
00952 double RL;
00953 float rlu,rlu1;
00954
00955
00956
00957 N=m_np;
00958 NP=0;
00959 for (I=1;I<N+1;I++){
00960 NP++;
00961 P[NP][1]=m_mom(I-1,1) ;
00962 P[NP][2]=m_mom(I-1,2) ;
00963 P[NP][3]=m_mom(I-1,3) ;
00964 P[NP][4]=m_mom(I-1,4) ;
00965 P[NP][5]=0;
00966 }
00967
00968
00969
00970 N=NP;
00971
00972 for (J1=1;J1<4;J1++) {
00973 for (J2=J1;J2<4;J2++) {
00974 SM[J1][J2]=0.;
00975 }
00976 }
00977 PS=0.;
00978 for (I=1;I<N+1;I++) {
00979 PA=sqrt(pow(P[I][1],2)+pow(P[I][2],2)+pow(P[I][3],2));
00980
00981
00982
00983 PWT=1.;
00984 for (J1=1;J1<4;J1++) {
00985 for (J2=J1;J2<4;J2++) {
00986 SM[J1][J2]=SM[J1][J2]+PWT*P[I][J1]*P[I][J2];
00987 }
00988 }
00989 PS=PS+PWT*PA*PA;
00990 }
00991
00992 if(NP<2) {
00993 SPH=-1.;
00994 APL=-1.;
00995 return;
00996 }
00997 for (J1=1;J1<4;J1++) {
00998 for (J2=J1;J2<4;J2++) {
00999 SM[J1][J2]=SM[J1][J2]/PS;
01000 }
01001 }
01002
01003 SQ=(SM[1][1]*SM[2][2]+SM[1][1]*SM[3][3]+SM[2][2]*SM[3][3]
01004 -pow(SM[1][2],2)
01005 -pow(SM[1][3],2)-pow(SM[2][3],2))/3.-1./9.;
01006 SR=-0.5*(SQ+1./9.+SM[1][1]*pow(SM[2][3],2)+SM[2][2]*pow(SM[1][3],2)+SM[3][3]*
01007 pow(SM[1][2],2)-SM[1][1]*SM[2][2]*SM[3][3])+SM[1][2]*SM[1][3]*SM[2][3]+1./27.;
01008
01009 SP=TMath::Cos(TMath::ACos(TMath::Max(TMath::Min(SR/TMath::Sqrt(-pow(SQ,3)),1.),-1.))/3.);
01010
01011 P[N+1][4]=1./3.+TMath::Sqrt(-SQ)*TMath::Max(2.*SP,TMath::Sqrt(3.*(1.-SP*SP))-SP);
01012 P[N+3][4]=1./3.+TMath::Sqrt(-SQ)*TMath::Min(2.*SP,-TMath::Sqrt(3.*(1.-SP*SP))-SP);
01013 P[N+2][4]=1.-P[N+1][4]-P[N+3][4];
01014 if (P[N+2][4]> 1E-5) {
01015
01016
01017 for (I=1;I<4;I=I+2) {
01018 for (J1=1;J1<4;J1++) {
01019 SV[J1][J1]=SM[J1][J1]-P[N+I][4];
01020 for (J2=J1+1;J2<4;J2++) {
01021 SV[J1][J2]=SM[J1][J2];
01022 SV[J2][J1]=SM[J1][J2];
01023 }
01024 }
01025 SMAX=0.;
01026 for (J1=1;J1<4;J1++) {
01027 for (J2=1;J2<4;J2++) {
01028 if(std::fabs(SV[J1][J2])>SMAX) {
01029 JA=J1;
01030 JB=J2;
01031 SMAX=std::fabs(SV[J1][J2]);
01032 }
01033 }
01034 }
01035 SMAX=0.;
01036 for (J3=JA+1;J3<JA+3;J3++) {
01037 J1=J3-3*((J3-1)/3);
01038 RL=SV[J1][JB]/SV[JA][JB];
01039 for (J2=1;J2<4;J2++) {
01040 SV[J1][J2]=SV[J1][J2]-RL*SV[JA][J2];
01041 if (std::fabs(SV[J1][J2])>SMAX) {
01042 JC=J1;
01043 SMAX=std::fabs(SV[J1][J2]);
01044 }
01045 }
01046 }
01047 JB1=JB+1-3*(JB/3);
01048 JB2=JB+2-3*((JB+1)/3);
01049 P[N+I][JB1]=-SV[JC][JB2];
01050 P[N+I][JB2]=SV[JC][JB1];
01051 P[N+I][JB]=-(SV[JA][JB1]*P[N+I][JB1]+SV[JA][JB2]*P[N+I][JB2])/
01052 SV[JA][JB];
01053 PA=TMath::Sqrt(pow(P[N+I][1],2)+pow(P[N+I][2],2)+pow(P[N+I][3],2));
01054
01055 float pa=P[N-1][I];
01056 rlu=std::fabs(pa)-std::fabs(int(pa)*1.);
01057 rlu1=std::fabs(pa*pa)-std::fabs(int(pa*pa)*1.);
01058 SGN=pow((-1.),1.*int(rlu+0.5));
01059 for (J=1;J<4;J++) {
01060 P[N+I][J]=SGN*P[N+I][J]/PA;
01061 }
01062 }
01063
01064
01065 SGN=pow((-1.),1.*int(rlu1+0.5));
01066 P[N+2][1]=SGN*(P[N+1][2]*P[N+3][3]-P[N+1][3]*P[N+3][2]);
01067 P[N+2][2]=SGN*(P[N+1][3]*P[N+3][1]-P[N+1][1]*P[N+3][3]);
01068 P[N+2][3]=SGN*(P[N+1][1]*P[N+3][2]-P[N+1][2]*P[N+3][1]);
01069
01070
01071 SPH=1.5*(P[N+2][4]+P[N+3][4]);
01072 APL=1.5*P[N+3][4];
01073
01074 }
01075
01076
01077 double etstot[4];
01078 double eltot[4];
01079 double sum23=0;
01080 double sum22=0;
01081 double sum33=0;
01082 double pina[4];
01083 double ax[4], ay[4], az[4];
01084 for (int ia=1;ia<4;ia++) {
01085 etstot[ia]=0;
01086 eltot[ia]=0;
01087 pina[ia]=0;
01088 ax[ia]=P[N+ia][1];
01089 ay[ia]=P[N+ia][2];
01090 az[ia]=P[N+ia][3];
01091 ax[ia]/=sqrt(ax[ia]*ax[ia]+ay[ia]*ay[ia]+az[ia]*az[ia]);
01092 ay[ia]/=sqrt(ax[ia]*ax[ia]+ay[ia]*ay[ia]+az[ia]*az[ia]);
01093 az[ia]/=sqrt(ax[ia]*ax[ia]+ay[ia]*ay[ia]+az[ia]*az[ia]);
01094 }
01095
01096 for (int k =0 ; k<m_np ; k++) {
01097 double eta,phi;
01098 getetaphi(m_mom(k,1),m_mom(k,2),m_mom(k,3),eta,phi);
01099 double W=exp(-std::fabs(eta*1.0));
01100 for (int ia=1;ia<4;ia++) {
01101 double e=sqrt(m_mom(k,1)*m_mom(k,1) +
01102 m_mom(k,2)*m_mom(k,2) +
01103 m_mom(k,3)*m_mom(k,3));
01104 double el=ax[ia]*m_mom(k,1) + ay[ia]*m_mom(k,2) + az[ia]*m_mom(k,3);
01105 pina[ia]=el;
01106 double ets=(e*e-el*el);
01107 etstot[ia]+=ets*W;
01108 eltot[ia]+=el*el*W;
01109 }
01110 double a2=pina[2];
01111 double a3=pina[3];
01112
01113
01114
01115 sum22+=a2*a2*W;
01116 sum23+=a2*a3*W;
01117 sum33+=a3*a3*W;
01118 }
01119
01120
01121 double pi=3.1415927;
01122 double phi=pi/2.0;
01123 double phip=pi/2.0;
01124 double a=sum23;
01125 double c=-a;
01126 double b=sum22-sum33;
01127 double disc=b*b-4*a*c;
01128
01129 if (disc>=0) {
01130 double x1=(sqrt(disc)-b)/2/a;
01131 double x2=(-sqrt(disc)-b)/2/a;
01132 phi=atan(x1);
01133 phip=atan(x2);
01134 if (phi<0) phi=2.*pi+phi;
01135 if (phip<0) phip=2.*pi+phip;
01136 }
01137 double p21=sum22*cos(phi)*cos(phi)+sum33*sin(phi)*sin(phi)+2*sum23*cos(phi)*sin(phi);
01138 double p31=sum22*sin(phi)*sin(phi)+sum33*cos(phi)*cos(phi)-2*sum23*cos(phi)*sin(phi);
01139
01140 double p22=sum22*cos(phip)*cos(phip)+sum33*sin(phip)*sin(phip)+2*sum23*cos(phip)*sin(phip);
01141 double p32=sum22*sin(phip)*sin(phip)+sum33*cos(phip)*cos(phip)-2*sum23*cos(phip)*sin(phip);
01142
01143
01144 double d1=std::fabs(p31*p31 - p21*p21);
01145 double d2=std::fabs(p32*p32 - p22*p22);
01146
01147
01148
01149 p2=p21;
01150 p3=p31;
01151 if (d2>d1) {
01152 p2=p22;
01153 p3=p32;
01154 }
01155 pnorm=sqrt(eltot[1]);
01156 if (p2>p3) {
01157 p3=sqrt(p3);
01158 p2=sqrt(p2);
01159 }else {
01160 double p4=p3;
01161 p3=sqrt(p2);
01162 p2=sqrt(p4);
01163 }
01164
01165
01166
01167
01168 return;
01169 }
01170
01171
01172
01173 void TopologyWorker::planes_thrust(double& pnorm, double& p2, double& p3) {
01174 TVector3 thrustaxis=thrustAxis();
01175 TVector3 majoraxis=majorAxis();
01176 TVector3 minoraxis=minorAxis();
01177
01178 double etstot[4];
01179 double eltot[4];
01180 double sum23=0;
01181 double sum22=0;
01182 double sum33=0;
01183 double pina[4];
01184 double ax[4], ay[4], az[4];
01185 ax[1]=thrustaxis(0); ay[1]=thrustaxis(1); az[1]=thrustaxis(2);
01186 ax[2]=minoraxis(0); ay[2]=minoraxis(1); az[2]=minoraxis(2);
01187 ax[3]=majoraxis(0); ay[3]=majoraxis(1); az[3]=majoraxis(2);
01188 for (int ia=1;ia<4;ia++) {
01189 etstot[ia]=0;
01190 eltot[ia]=0;
01191 pina[ia]=0;
01192 ax[ia]/=sqrt(ax[ia]*ax[ia]+ay[ia]*ay[ia]+az[ia]*az[ia]);
01193 ay[ia]/=sqrt(ax[ia]*ax[ia]+ay[ia]*ay[ia]+az[ia]*az[ia]);
01194 az[ia]/=sqrt(ax[ia]*ax[ia]+ay[ia]*ay[ia]+az[ia]*az[ia]);
01195 }
01196
01197 for (int k =0 ; k<m_np ; k++) {
01198 for (int ia=1;ia<4;ia++) {
01199 double e=sqrt(m_mom(k,1)*m_mom(k,1) +
01200 m_mom(k,2)*m_mom(k,2) +
01201 m_mom(k,3)*m_mom(k,3));
01202 double el=ax[ia]*m_mom(k,1) + ay[ia]*m_mom(k,2) + az[ia]*m_mom(k,3);
01203 pina[ia]=el;
01204 double ets=(e*e-el*el);
01205 etstot[ia]+=ets;
01206 eltot[ia]+=std::fabs(el);
01207 }
01208 double a2=pina[2];
01209 double a3=pina[3];
01210
01211
01212
01213 sum22+=a2*a2;
01214 sum23+=a2*a3;
01215 sum33+=a3*a3;
01216 }
01217
01218
01219 double pi=3.1415927;
01220 double phi=pi/2.0;
01221 double phip=pi/2.0;
01222 double a=sum23;
01223 double c=-a;
01224 double b=sum22-sum33;
01225 double disc=b*b-4*a*c;
01226
01227 if (disc>=0) {
01228 double x1=(sqrt(disc)-b)/2/a;
01229 double x2=(-sqrt(disc)-b)/2/a;
01230 phi=atan(x1);
01231 phip=atan(x2);
01232 if (phi<0) phi=2.*pi+phi;
01233 if (phip<0) phip=2.*pi+phip;
01234 }
01235 double p21=sum22*cos(phi)*cos(phi)+sum33*sin(phi)*sin(phi)+2*sum23*cos(phi)*sin(phi);
01236 double p31=sum22*sin(phi)*sin(phi)+sum33*cos(phi)*cos(phi)-2*sum23*cos(phi)*sin(phi);
01237
01238 double p22=sum22*cos(phip)*cos(phip)+sum33*sin(phip)*sin(phip)+2*sum23*cos(phip)*sin(phip);
01239 double p32=sum22*sin(phip)*sin(phip)+sum33*cos(phip)*cos(phip)-2*sum23*cos(phip)*sin(phip);
01240
01241
01242 double d1=std::fabs(p31*p31 - p21*p21);
01243 double d2=std::fabs(p32*p32 - p22*p22);
01244
01245
01246
01247 p2=p21;
01248 p3=p31;
01249 if (d2>d1) {
01250 p2=p22;
01251 p3=p32;
01252 }
01253 pnorm=sqrt(etstot[1]);
01254 if (p2>p3) {
01255 p3=sqrt(p3);
01256 p2=sqrt(p2);
01257 }else {
01258 double p4=p3;
01259 p3=sqrt(p2);
01260 p2=sqrt(p4);
01261 }
01262
01263
01264
01265
01266 return;
01267 }
01268
01269
01270 void TopologyWorker::fowo() {
01271
01272 m_fowo_called=true;
01273
01274 float H10=-1;
01275 float H20=-1;
01276 float H30=-1;
01277 float H40=-1;
01278 float H50=-1;
01279 float H60=-1;
01280 if (1==1) {
01281 float P[1000][6],H0,HD,CTHE;
01282 int N,NP,I,J,I1,I2;
01283 H0=HD=0.;
01284 N=m_np;
01285 NP=0;
01286 for (I=1;I<N+1;I++){
01287 NP++;
01288 P[NP][1]=m_mom(I-1,1) ;
01289 P[NP][2]=m_mom(I-1,2) ;
01290 P[NP][3]=m_mom(I-1,3) ;
01291 P[NP][4]=m_mom(I-1,4) ;
01292 P[NP][5]=m_mom(I-1,5) ;
01293 }
01294
01295 N=NP;
01296 NP=0;
01297
01298 for (I=1;I<N+1;I++) {
01299 NP=NP+1;
01300 for (J=1;J<5;J++) {
01301 P[N+NP][J]=P[I][J];
01302 }
01303
01304 P[N+NP][4]=sqrt(pow(P[I][1],2)+pow(P[I][2],2)+pow(P[I][3],2));
01305
01306 P[N+NP][5]=sqrt(pow(P[I][1],2)+pow(P[I][2],2));
01307 H0=H0+P[N+NP][5];
01308 HD=HD+pow(P[N+NP][5],2);
01309 }
01310 H0=H0*H0;
01311
01312
01313
01314
01315 if (NP<2) {
01316 H10=-1.;
01317 H20=-1.;
01318 H30=-1.;
01319 H40=-1.;
01320 H50=-1.;
01321 H60=-1.;
01322 return;
01323 }
01324
01325
01326 H10=0.;
01327 H20=0.;
01328 H30=0.;
01329 H40=0.;
01330 H50=0.;
01331 H60=0.;
01332 for (I1=N+1;I1<N+NP+1;I1++) {
01333 for (I2=I1+1;I2<N+NP+1;I2++) {
01334 CTHE=(P[I1][1]*P[I2][1]+P[I1][2]*P[I2][2]+P[I1][3]*P[I2][3])/
01335 (P[I1][4]*P[I2][4]);
01336 H10=H10+P[I1][5]*P[I2][5]*CTHE;
01337 double C2=(1.5*CTHE*CTHE-0.5);
01338 H20=H20+P[I1][5]*P[I2][5]*C2;
01339 double C3=(2.5*CTHE*CTHE*CTHE-1.5*CTHE);
01340 H30=H30+P[I1][5]*P[I2][5]*C3;
01341
01342 double C4=(7*CTHE*C3 - 3*C2)/4.;
01343 double C5=(9*CTHE*C4 - 4*C3)/5.;
01344 double C6=(11*CTHE*C5 - 5*C4)/6.;
01345
01346
01347
01348
01349
01350 H40=H40+P[I1][5]*P[I2][5]*C4;
01351 H50=H50+P[I1][5]*P[I2][5]*C5;
01352 H60=H60+P[I1][5]*P[I2][5]*C6;
01353 }
01354 }
01355
01356 H10=(HD+2.*H10)/H0;
01357 H20=(HD+2.*H20)/H0;
01358 H30=(HD+2.*H30)/H0;
01359 H40=(HD+2.*H40)/H0;
01360 H50=(HD+2.*H50)/H0;
01361 H60=(HD+2.*H60)/H0;
01362 m_h10=H10;
01363 m_h20=H20;
01364 m_h30=H30;
01365 m_h40=H40;
01366 m_h50=H50;
01367 m_h60=H60;
01368 }
01369
01370 }
01371
01372 double TopologyWorker::get_h10() {
01373 if (!m_fowo_called) fowo();
01374 return m_h10;
01375 }
01376 double TopologyWorker::get_h20() {
01377 if (!m_fowo_called) fowo();
01378 return m_h20;
01379 }
01380 double TopologyWorker::get_h30() {
01381 if (!m_fowo_called) fowo();
01382 return m_h30;
01383 }
01384 double TopologyWorker::get_h40() {
01385 if (!m_fowo_called) fowo();
01386 return m_h40;
01387 }
01388
01389 double TopologyWorker::get_h50() {
01390 if (!m_fowo_called) fowo();
01391 return m_h50;
01392 }
01393 double TopologyWorker::get_h60() {
01394 if (!m_fowo_called) fowo();
01395 return m_h60;
01396 }
01397
01398
01399 double TopologyWorker::get_sphericity() {
01400 if (!m_sanda_called) sanda();
01401 return m_sph;
01402 }
01403 double TopologyWorker::get_aplanarity() {
01404 if (!m_sanda_called) sanda();
01405 return m_apl;
01406 }
01407
01408 void TopologyWorker::getetaphi(double px, double py, double pz, double& eta, double& phi) {
01409
01410 double pi=3.1415927;
01411
01412 double p=sqrt(px*px+py*py+pz*pz);
01413
01414 double th=pi/2.;
01415 if (p!=0) {
01416 th=acos(pz/p);
01417 }
01418 float thg=th;
01419 if (th<=0) {
01420 thg = pi + th;
01421 }
01422 eta=-9.;
01423 if (tan( thg/2.0 )>0.000001) {
01424 eta = -log( tan( thg/2.0 ) );
01425 }
01426 phi = atan2(py,px);
01427 if(phi<=0) phi += 2.0*pi;
01428 return;
01429 }
01430
01431
01432
01433 void TopologyWorker::sumangles(float& sdeta, float& sdr) {
01434 double eta1,eta2,phi1,phi2,deta,dphi,dr;
01435 m_sumangles_called=true;
01436 sdeta=0;
01437 sdr=0;
01438 for (int k=0;k<m_np;k++){
01439 for (int kp=k;kp<m_np;kp++){
01440 getetaphi(m_mom(k,1) , m_mom(k,2), m_mom(k,3), eta1,phi1);
01441 getetaphi(m_mom(kp,1) , m_mom(kp,2), m_mom(kp,3), eta2,phi2);
01442 dphi=std::fabs(phi1-phi2);
01443 if (dphi>3.1415) dphi=2*3.1415-dphi;
01444 deta=std::fabs(eta1-eta2);
01445 dr=sqrt(dphi*dphi+deta*deta);
01446 sdeta+=deta;
01447 sdr+=dr;
01448 }
01449 }
01450 return;
01451 }
01452
01453
01454
01455 Int_t TopologyWorker::iPow(Int_t man, Int_t exp)
01456 {
01457 Int_t ans = 1;
01458 for( Int_t k = 0; k < exp; k++)
01459 {
01460 ans = ans*man;
01461 }
01462 return ans;
01463 }
01464
01465
01466
01467 void TopologyWorker::CalcWmul(){
01468
01469 Int_t njets = m_np;
01470 double result=0;
01471 for(Int_t ijet=0; ijet<njets-1; ijet++){
01472 double emin=55;
01473 double emax=55;
01474 if(CalcPt(ijet)<55)
01475 emax=CalcPt(ijet);
01476 if(CalcPt(ijet+1)<55)
01477 emin=CalcPt(ijet+1);
01478 result+=0.5 * (emax*emax-emin*emin)*(ijet+1);
01479 }
01480 double elo=15;
01481 if(CalcPt(njets-1)>elo){
01482 elo=CalcPt(njets-1);
01483 }
01484
01485 result+=0.5 * (elo*elo-(15*15))*(njets);
01486 result/=((55*55)-100)/2.0;
01487 m_njetsweighed=result;
01488 }
01489
01490
01491
01492 void TopologyWorker::CalcSqrts(){
01493 TLorentzVector event(0,0,0,0);
01494 TLorentzVector worker(0,0,0,0);
01495
01496 for(int i=0; i< m_np; i++){
01497 double energy=m_mom(i,4);
01498 if(m_mom(i,4)<0.00001)
01499 energy=sqrt(pow(m_mom(i,1),2)+pow(m_mom(i,2),2)+pow(m_mom(i,3),2));
01500
01501 worker.SetXYZT(m_mom(i,1),m_mom(i,2),m_mom(i,3),energy);
01502 event+=worker;
01503 }
01504 m_sqrts=event.M();
01505 }
01506
01507
01508 void TopologyWorker::CalcHTstuff(){
01509 m_ht=0;
01510 m_ht3=0;
01511 m_et56=0;
01512 m_et0=0;
01513 double ptlead=0;
01514 double h=0;
01515 for(int i=0; i< m_np; i++){
01516
01517 m_ht+=CalcPt(i);
01518 h+=m_mom(i,4);
01519 if(i>1)
01520 m_ht3+=CalcPt(i);
01521 if(i==5)
01522 m_et56=sqrt(CalcPt(i)*CalcPt(i-1));
01523 }
01524
01525 for(int j=0; j< m_np2; j++){
01526
01527 if(ptlead<CalcPt2(j))
01528 ptlead=CalcPt2(j);
01529
01530 }
01531 if(m_ht>0.0001){
01532 m_et0=ptlead/m_ht;
01533
01534 }
01535 if(h>0.00001)
01536 m_centrality=m_ht/h;
01537 }
01538
01540