ZeeAnalysis.C

00001 //#include "UserTreeAnalysis.H"   // remove if copied to user working directory
00002 #include <stdio.h>
00003 #include <assert.h>
00004 #include <math.h>
00005 #include "MC4Vector.H"
00006 #include "HEPParticle.H"
00007 #include "TH1.h"
00008 #include "Setup.H"
00009 #include "TObjArray.h"
00010 #include "TMath.h"
00011 #define PI 3.141592653
00012 
00013 inline bool ifSofty(int Id,int nparams, double *params){
00014   if (nparams<5 && Id==22) return true;   // to remove photons only
00015   for (int i=nparams-1; i>3; i--)
00016     if (Id==params[i]) return true;       // to remove all what is in params from nparams down to 4
00017   return false;
00018 }
00019 
00020 // very similar to  MC_FillUserHistogram from Generate.cxx
00021 inline void fillUserHisto(char *name,double val, double weight=1.0, 
00022                           double min=Setup::bin_min[0][0], 
00023                           double max=Setup::bin_max[0][0]){
00024 
00025     TH1D *h=(TH1D*)(Setup::user_histograms->FindObject(name));
00026     if(!h){
00027       h=new TH1D(name,name,Setup::nbins[0][0],min,max);
00028       if(!h) return;
00029       Setup::user_histograms->Add(h);
00030       //      printf("user histogram created %s\n", name);
00031 }
00032     h->Fill(val,weight);
00033 
00034 }
00035 double angle(double X, double Y){
00036 
00037  
00038     double  an=0.0; 
00039     double R=sqrt(X*X+Y*Y); 
00040     //    if(R<pow(10,-20)) printf(" angle this time X %f\n", 10000*X);
00041     if(R<pow(10,-20)) return an;
00042 
00043     if(TMath::Abs(X)/R<0.8) 
00044     { 
00045         an=acos(X/R);
00046         if(Y<0 && an>0) an=-an;
00047         if(Y>0 && an<0) an=-an; 
00048     }
00049     else 
00050     { 
00051       an=asin(Y/R);
00052       if(X<0 && an>=0.)  an=PI-an; 
00053       else if(X<0.)      an=-PI-an; 
00054        
00055     } 
00056   return an;
00057 }
00058 
00059 int counter=0;
00060 
00061 int ZeeAnalysis(HEPParticle *mother,HEPParticleList *stableDaughters, int nparams, double *params)
00062 {
00063   // THIS IS EXAMPLE of userTreeAnalysis. It acts on list of particle X (under our MC-test) final
00064   // daughters. Of course in real life many options may need to be introduced by the user.
00065     // PARAMETERS:
00066     // params[0] threshold on Energy (or p_T) of particle expressed as a fraction of mother's 
00067     //           Energy (or p_T) in mother's frame. If not specified - default is 0.05 
00068     //           Note that setting it to zero is possible.
00069     // params[1] maximum number of left soft-suspected particles. 
00070     //           0: (default) all listed particles are removed, even if hard
00071     //           1: 2:  one two etc  removable particles are kept (at most)
00072     //          -1: this option is off, all hard particles stay.
00073     // params[2] control tag on discriminating particles, 
00074     //           0: (default) Energy in decaying particle rest frame
00075     //           1: Energy in lab frame
00076     //           2: p_T with respect to beam direction in lab frame.
00077     // params[3] type of action to be applied on removed daughters 
00078     //           0: (default) nothing is done, they are lost
00079     //           1:  algorithm as in studies on PHOTOS is used, see papers P. Golonka and Z. Was.
00080     // params[4] from this paramter on PDG id-s of particles to be removed can be added. 
00081     //           Default is that only photons are removed.
00082     // To get to origin of our particle (mother) one need to go after UserEventAnalysis,
00083     // instructive example will be given later. 
00084     assert(mother!=0);
00085     assert(stableDaughters!=0);
00086 
00087 
00088     double threshold=0.05, leftmax=0.0, selector=0.0, actLost=0.0; 
00089     if (nparams>0 && params==0)  return -1;
00090     if (nparams>0) threshold=params[0]; 
00091     if (nparams>1) leftmax  =params[1]; 
00092     if (nparams>2) selector =params[2]; 
00093     if (nparams>3) actLost  =params[3]; 
00094 
00095     
00096     HEPParticleList *lostDaughters=new HEPParticleList;    
00097     double pt_limit=threshold*mother->GetM();
00098 
00099     HEPParticleListIterator daughters (*stableDaughters);
00100     int nphot=0;
00101     double ephot=pow(10,22);
00102     bool redo=false;
00103     for (HEPParticle *part=daughters.first(); part!=0; (redo)? part=part:part=daughters.next() ) {
00104       if(redo) redo=false;
00105         MC4Vector d4(part->GetE(),part->GetPx(),part->GetPy(),part->GetPz(),part->GetM());
00106         // boost to mother's frame:
00107         d4.Boost(mother->GetPx(),mother->GetPy(),mother->GetPz(),mother->GetE(),mother->GetM());
00108 
00109 
00110 
00111         double p_t=d4.GetX0();  // default
00112         switch ((int) selector)
00113         {
00114         case 1: p_t=part->GetE(); break;
00115         case 2: p_t=d4.Xt(); break;
00116         }
00117 
00118 
00119         if(ifSofty(part->GetPDGId(),nparams,params)) nphot++;
00120         if ( ifSofty(part->GetPDGId(),nparams,params) && p_t < pt_limit) {
00121           //      printf("Odrzucamy! %i\n",count++);
00122           nphot=0;
00123             lostDaughters->push_back(part);
00124             stableDaughters->remove(part);
00125             part=daughters.first();
00126             redo=true;
00127             continue;
00128         }
00129         if( ifSofty(part->GetPDGId(),nparams,params) && ephot>p_t) ephot=p_t;
00130     }
00131     while(nphot>(int) leftmax)
00132     {
00133       double ephot1=pow(10,22);
00134       redo=false;
00135       for (HEPParticle *part=daughters.first(); part!=0; (redo)? part=part:part=daughters.next() ) {
00136         if(redo) redo=false;
00137         MC4Vector d4(part->GetE(),part->GetPx(),part->GetPy(),part->GetPz(),part->GetM());
00138         // boost to mother's frame:
00139         d4.Boost(mother->GetPx(),mother->GetPy(),mother->GetPz(),mother->GetE(),mother->GetM());
00140 
00141 
00142         double p_t=d4.GetX0();  // default
00143         switch ((int) selector)
00144         {
00145         case 1: p_t=part->GetE(); break;
00146         case 2: p_t=d4.Xt(); break;
00147         }
00148 
00149 
00150         if ( ifSofty(part->GetPDGId(),nparams,params) && p_t == ephot) {
00151 
00152           nphot--;
00153             lostDaughters->push_back(part);
00154             stableDaughters->remove(part);
00155             ephot1=pow(10,22);
00156              part=daughters.first();
00157              redo=true;
00158              continue;
00159         }
00160         if( ifSofty(part->GetPDGId(),nparams,params) && ephot1>p_t) ephot1=p_t;
00161       }
00162       ephot=ephot1;
00163 
00164     }
00165 
00166     // ##############################################################
00167     // here functionality of removig some daughters is finished
00168     // now we reconsider what to do with them
00169     // ##############################################################
00170     // delete lostDaughters;
00171     // return 1; 
00172 
00173 
00174     // ##############################################################
00175     // Now: What to do with lost daughters?
00176     // ##############################################################
00177     //    lostDaughters->ls();
00178     int version=(int) actLost;
00179         
00180     switch(version)
00181       {
00182       case 1:  // add lost to charged
00183         {
00184           HEPParticleListIterator lost (*lostDaughters);
00185           for (HEPParticle *lostpart=lost.first(); lostpart!=0; lostpart=lost.next() ) {
00186             HEPParticle* Best=0;
00187             double searchvirt=pow(10.0,22);
00188             MC4Vector VV;
00189             for( HEPParticle *part=daughters.first(); part!=0; part=daughters.next() ){
00190               if(part->GetCharge()==0) continue;
00191               VV=lostpart->GetP4()+part->GetP4();
00192               VV.AdjustM();
00193               if (VV.GetM()<searchvirt) {searchvirt=VV.GetM(); Best=part;}
00194             }
00195             if(Best) {
00196               Best->SetPx(Best->GetPx()+lostpart->GetPx());
00197               Best->SetPy(Best->GetPy()+lostpart->GetPy());
00198               Best->SetPz(Best->GetPz()+lostpart->GetPz());
00199               Best->SetE (Best->GetE ()+lostpart->GetE ());
00200             }
00201           }
00202         break;
00203         }
00204       default: break; // do nothing
00205       }
00206 
00207     delete lostDaughters;
00208 
00209 
00210         /*
00211                 Create e+e- inv. mass histogram
00212         */
00213         double px=0,py=0,pz=0,e=0;
00214 
00215         //loop over all daughters
00216         for(HEPParticle *part=daughters.first();part!=0;part=daughters.next())
00217         {
00218                 if(abs(part->GetPDGId())!=11) continue;
00219                 //Sum e+ e- 4-vectors
00220                 px+=part->GetPx();
00221                 py+=part->GetPy();
00222                 pz+=part->GetPz();
00223                 e+=part->GetE();
00224         }
00225         double mass=sqrt(e*e-px*px-py*py-pz*pz);
00226         char mm[] = "e+e-mass_distribution";
00227         fillUserHisto(mm,mass,1.0,40.0,140.0);
00228 
00229         //For ZeeAnalysis - other histograms are not needed
00230     bool activateUserHist=false;
00231     if(!activateUserHist) return 1;
00232 
00233     // segmet of user defined histograms
00234 
00235     double X=mother->GetPx(), Y=mother->GetPy(), Z=mother->GetPz(); 
00236     // double E=mother->GetE(),  MM=mother->GetM();
00237 
00238 
00239     double pt=sqrt(X*X+Y*Y);
00240     double eta=log((sqrt(pt*pt+Z*Z)+TMath::Abs(Z))/pt);
00241     if(Z<0 && eta>0) eta=-eta;
00242     if(Z>0 && eta<0) eta=-eta;
00243     double phi=angle(X,Y);
00244     char hist1[] = "mother-PT";
00245     char hist2[] = "mother-eta";
00246     char hist3[] = "mother-phi";
00247     fillUserHisto(hist1,pt,1.0,0.0,100.0);
00248     fillUserHisto(hist2,eta,1.0,-8.0,8.0);
00249     fillUserHisto(hist3,phi,1.0,-PI,PI);
00250     return 1;
00251 }
Generated on Sun Oct 20 20:23:56 2013 for C++InterfacetoPHOTOS by  doxygen 1.6.3