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