ROOT logo
HYDRA - THE HADES ANALYSIS PACKAGE » (UNKNOWN) » HParticleCutRange

class HParticleCutRange: public TNamed

_HADES_CLASS_DESCRIPTION



 HParticleCutRange

 Helper class for cuts.
 A cut has a name (should be unique), a cut number (better unique)
 and a condition as argument of creation. The cut object owns counters
 for the number of calls and the number of failed evaluations. The
 The counters cand be used to monitor the efficiency of the cut
 If the TCutG object belonging to the cut is set  evalG(var,var,version)
 will evaluate the graph. cut. As option TF1 objects cann be used to
 cut (see the exapmple macro below)
 The cuts can be inversed by setInverse(Bool_t). In fact
 cuts will than select outside the band low-up range and
 for the graph cut isInside(var,var)==kFALSE will be used for selection.


 USAGE:
  // outside eventloop:
  HParticleCutRange cut1("BetaCut",1,0.9,1.1);
  // will create a cut object selecting beta range and

  //  inside eventloop:
  if(cut1->eval(beta,0)) { // true if condition is fullfilled, stat counted for version 0

  }
  if(cut1->eval(beta,1)) { // true if condition is fullfilled, stat counted for version 1

  }


  // after eventloop
  cut1.print(); // show cut name,number,cut statistics and range

void testCutRange()
{
    // example macro for HParticleCutRange demonstrating the
    // different features. The cut will be applied on the
    // mass of an particle in this example.

    // Cuts can be applied in 3 ways:

    // CUTMODE:
    // 1. 1-dim cuts (upper+lower limit)
    // 2. 2-dim cut using TCutG
    // 3. 2-dim cut using TF1

    // In TCutG or TF1 mode the example cuts in mass depends on the momentum.
    // The TCutG is owned by the cut and the user just has to set the
    // points. TF1 have to be contructed outside the clas due to the
    // available contructors of TF1. The user ahs to set the pointers
    // to the TF1. Remark: if the TF1 is contructed from a user function
    // the cut will not work when read back from root file.

    // TF1MODE:
    // There are 4 options :
    // 1. low+up  (1 seperate TF1 for lower and upper cut)
    // 2. low     (1 TF1 for lower cut)
    // 3. up      (1 TF1 for upper cut)
    // 3. mean+width (1 TF1 for mean of the cut, 1 TF1 for the width of the cut)


    Bool_t   invert   = kFALSE;  // invert cut : default: kFALSE
    Int_t    cutmode  = 2; //0=range,1=cutG,2=TF1
    TString  tf1mode  = "low+up"; // TF1 mode: low+up, low,up,mean+width

    Double_t slope    = .10;

    Double_t mass_pion   = HPhysicsConstants::mass(8);
    Double_t mass_proton = HPhysicsConstants::mass(14);


    HParticleCutRange cut1("massCut",1,mass_proton-50,mass_proton+50);
    cut1.setInverse(invert);


    //--------------------------------------------------
    // TF1 mode
    // Cuts in mass depend on mommentum
    // low+up , low , up mode
    TF1* flow = new TF1("flow","pol1",0,1000);
    flow->SetParameter(0,mass_pion);
    flow->SetParameter(1,-slope);

    TF1* fup = new TF1("fup","pol1"  ,0,1000);
    fup->SetParameter(0,mass_pion);
    fup->SetParameter(1,slope);
    fup->SetLineColor(4);

    // mean+width mode
    TF1* fmean = new TF1("fmean","pol1",0,1000);
    fmean->SetParameter(0,mass_pion);
    fmean->SetParameter(1,0);

    TF1* fwidth = new TF1("fwidth","pol1",0,1000);
    fwidth->SetParameter(0,0);
    fwidth->SetParameter(1,slope);
    fwidth->SetLineColor(4);

    if(tf1mode=="low+up")    cut1.setTF1(flow,fup,tf1mode);
    if(tf1mode=="low")       cut1.setTF1(flow,  0,tf1mode);
    if(tf1mode=="up")        cut1.setTF1(   0,fup,tf1mode);
    if(tf1mode=="mean+width")cut1.setTF1(flow,fup,tf1mode);
    //--------------------------------------------------


    //--------------------------------------------------
    // TCutG mode
    cut1.getTCutG()->SetPoint(0,   0,mass_proton);
    cut1.getTCutG()->SetPoint(1,1000,mass_proton - slope*1000);
    cut1.getTCutG()->SetPoint(2,1000,mass_proton + slope*1000);
    cut1.getTCutG()->SetPoint(3,   0,mass_proton);
    //--------------------------------------------------



    TH2F* hmomvsmass = new TH2F("hmassvsmom","mom vs mass",200,0,2000,200,0,1200);
    hmomvsmass->SetXTitle("mom [MeV/c]");
    hmomvsmass->SetYTitle("mass [MeV/c^{2}]");

    //--------------------------------------------------
    // create some dummy spectrum
    Double_t masses [] = {mass_pion,mass_proton};

    for(Int_t i = 0 ; i <10000; i ++){
        Double_t mom     =  gRandom->Rndm()*2000;
	Double_t massres =  mom*0.1;
        Int_t    index   = (Int_t) (gRandom->Rndm()*2);
        Double_t mass    = gRandom->Gaus(masses[index],massres);

	if(cutmode == 0 && cut1.eval(mass)) {
	    hmomvsmass->Fill(mom,mass);
	}

	if(cutmode == 1 && cut1.evalG(mom,mass)) {
	    hmomvsmass->Fill(mom,mass);
	}

	if(cutmode == 2 && cut1.evalF(mass,mom)) {
	    hmomvsmass->Fill(mom,mass);
	}
    }
    //--------------------------------------------------


    cut1.print();

    TCanvas* result = new TCanvas("result","result",1000,800);

    hmomvsmass->Draw("colz");


    if(cutmode == 1){
        cut1.getTCutG()->Draw("same");
    }

    if(cutmode == 2){
	if(cut1.getTF1Low())cut1.getTF1Low()->DrawCopy("same");
	if(cut1.getTF1Up() )cut1.getTF1Up() ->DrawCopy("same");
    }

    TFile* cutsfile = new TFile("myCuts.root","RECREATE");
    cut1.Write();
    cutsfile->Save();
    cutsfile->Close();


}

Function Members (Methods)

public:
HParticleCutRange(const HParticleCutRange&)
HParticleCutRange(TString name = "cut", Int_t num = -1, Double_t l = -1e30, Double_t u = 1e30)
virtual~HParticleCutRange()
voidTObject::AbstractMethod(const char* method) const
virtual voidTObject::AppendPad(Option_t* option = "")
virtual voidTObject::Browse(TBrowser* b)
static TClass*Class()
virtual const char*TObject::ClassName() const
virtual voidTNamed::Clear(Option_t* option = "")
virtual TObject*TNamed::Clone(const char* newname = "") const
virtual Int_tTNamed::Compare(const TObject* obj) const
virtual voidTNamed::Copy(TObject& named) const
virtual voidTObject::Delete(Option_t* option = "")MENU
virtual Int_tTObject::DistancetoPrimitive(Int_t px, Int_t py)
virtual voidTObject::Draw(Option_t* option = "")
virtual voidTObject::DrawClass() constMENU
virtual TObject*TObject::DrawClone(Option_t* option = "") constMENU
virtual voidTObject::Dump() constMENU
virtual voidTObject::Error(const char* method, const char* msgfmt) const
Bool_teval(Double_t var, UInt_t version = 0)
Bool_tevalF(Double_t var, Double_t var2, UInt_t version = 0)
Bool_tevalG(Double_t var, Double_t var2, UInt_t version = 0)
virtual voidTObject::Execute(const char* method, const char* params, Int_t* error = 0)
virtual voidTObject::Execute(TMethod* method, TObjArray* params, Int_t* error = 0)
virtual voidTObject::ExecuteEvent(Int_t event, Int_t px, Int_t py)
virtual voidTObject::Fatal(const char* method, const char* msgfmt) const
virtual voidTNamed::FillBuffer(char*& buffer)
virtual TObject*TObject::FindObject(const char* name) const
virtual TObject*TObject::FindObject(const TObject* obj) const
Int_tgetCutNumber()
Float_tgetCutRate(UInt_t version = 0)
virtual Option_t*TObject::GetDrawOption() const
static Long_tTObject::GetDtorOnly()
virtual const char*TObject::GetIconName() const
Double_tgetLowCut()
UInt_tgetMaxVersion()
virtual const char*TNamed::GetName() const
UInt_tgetNCall(UInt_t version = 0)
UInt_tgetNFail(UInt_t version = 0)
virtual char*TObject::GetObjectInfo(Int_t px, Int_t py) const
static Bool_tTObject::GetObjectStat()
virtual Option_t*TObject::GetOption() const
TCutG*getTCutG()
TF1*getTF1Low()
TF1*getTF1Up()
virtual const char*TNamed::GetTitle() const
virtual UInt_tTObject::GetUniqueID() const
Double_tgetUpCut()
virtual Bool_tTObject::HandleTimer(TTimer* timer)
virtual ULong_tTNamed::Hash() const
virtual voidTObject::Info(const char* method, const char* msgfmt) const
virtual Bool_tTObject::InheritsFrom(const char* classname) const
virtual Bool_tTObject::InheritsFrom(const TClass* cl) const
virtual voidTObject::Inspect() constMENU
voidTObject::InvertBit(UInt_t f)
virtual TClass*IsA() const
virtual Bool_tTObject::IsEqual(const TObject* obj) const
virtual Bool_tTObject::IsFolder() const
Bool_tTObject::IsOnHeap() const
virtual Bool_tTNamed::IsSortable() const
Bool_tTObject::IsZombie() const
virtual voidTNamed::ls(Option_t* option = "") const
voidTObject::MayNotUse(const char* method) const
virtual Bool_tTObject::Notify()
voidTObject::Obsolete(const char* method, const char* asOfVers, const char* removedFromVers) const
static voidTObject::operator delete(void* ptr)
static voidTObject::operator delete(void* ptr, void* vp)
static voidTObject::operator delete[](void* ptr)
static voidTObject::operator delete[](void* ptr, void* vp)
void*TObject::operator new(size_t sz)
void*TObject::operator new(size_t sz, void* vp)
void*TObject::operator new[](size_t sz)
void*TObject::operator new[](size_t sz, void* vp)
HParticleCutRange&operator=(const HParticleCutRange&)
virtual voidTObject::Paint(Option_t* option = "")
virtual voidTObject::Pop()
voidprint()
virtual voidTNamed::Print(Option_t* option = "") const
virtual Int_tTObject::Read(const char* name)
virtual voidTObject::RecursiveRemove(TObject* obj)
voidTObject::ResetBit(UInt_t f)
voidresetCounter()
virtual voidTObject::SaveAs(const char* filename = "", Option_t* option = "") constMENU
virtual voidTObject::SavePrimitive(ostream& out, Option_t* option = "")
voidset(Double_t l = -1e30, Double_t u = 1e30)
voidTObject::SetBit(UInt_t f)
voidTObject::SetBit(UInt_t f, Bool_t set)
virtual voidTObject::SetDrawOption(Option_t* option = "")MENU
static voidTObject::SetDtorOnly(void* obj)
voidsetInverse(Bool_t inv)
voidsetMaxCut(UInt_t max = 4)
virtual voidTNamed::SetName(const char* name)MENU
virtual voidTNamed::SetNameTitle(const char* name, const char* title)
static voidTObject::SetObjectStat(Bool_t stat)
voidsetTF1(TF1* flow, TF1* fup = 0, TString mode = "low+up")
virtual voidTNamed::SetTitle(const char* title = "")MENU
virtual voidTObject::SetUniqueID(UInt_t uid)
virtual voidShowMembers(TMemberInspector&)
virtual Int_tTNamed::Sizeof() const
virtual voidStreamer(TBuffer&)
voidStreamerNVirtual(TBuffer& ClassDef_StreamerNVirtual_b)
virtual voidTObject::SysError(const char* method, const char* msgfmt) const
Bool_tTObject::TestBit(UInt_t f) const
Int_tTObject::TestBits(UInt_t f) const
virtual voidTObject::UseCurrentStyle()
virtual voidTObject::Warning(const char* method, const char* msgfmt) const
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0)
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0) const
protected:
virtual voidTObject::DoError(int level, const char* location, const char* fmt, va_list va) const
voidTObject::MakeZombie()

Data Members

protected:
TStringTNamed::fNameobject identifier
TStringTNamed::fTitleobject title
private:
TCutG*fCutgraph cut object
Int_tfCutNumbernumber of the cut
Int_tfF1Mode0 = "low+up" , 1="low" , 2="up" , 3="mean+width"
TF1*fLowF1TFormular for low cut
TF1*fUpF1TFormular for upper cut
Bool_tfbInverseCutdefault kFALSE , kTRUE will invert selection
vector<ULong64_t>fctCallcount all call to the evaluation
vector<ULong64_t>fctFailcount evaluation == kFALSE
Double_tflowlower cut
UInt_tfmaxCuthow many version can be used at max
UInt_tfmaxVersionhow many version have been used
Double_tfupupper cut

Class Charts

Inheritance Inherited Members Includes Libraries
Class Charts

Function documentation

HParticleCutRange(TString name = "cut", Int_t num = -1, Double_t l = -1e30, Double_t u = 1e30)
~HParticleCutRange()
void setTF1(TF1* flow, TF1* fup = 0, TString mode = "low+up")
 Set the TF1 objects for low and up cut. The TF1
 will be used by evalF(var,var2) for cutting where
 the cuts in var will be a function of var2.
 There are the following modes implemented:

 mode = "low+up" : lower and upper cur are defined
                   by independent TF1 (band pass)
      = "low"    : low TF1 is used to defined a cut
                   for var < lowcut(var2) ( high pass)
      = "up"     : up TF1 is used to defined a cut
                   for var > upcut(var2) ( low pass)
      = "mean+width" : is used to define a cut
                       by mean(var2)(low TF1) + width(var2)(up TF1)
                       for var > mean + width && var < mean - width  ( band pass)
 the cuts can be inverted as the other cuts
void setMaxCut(UInt_t max = 4)
 set the max number of versions of stst counters
Bool_t eval(Double_t var, UInt_t version = 0)
Bool_t evalG(Double_t var, Double_t var2, UInt_t version = 0)
Bool_t evalF(Double_t var, Double_t var2, UInt_t version = 0)
void print()
void resetCounter()
Float_t getCutRate(UInt_t version = 0)
UInt_t getNCall(UInt_t version = 0)
UInt_t getNFail(UInt_t version = 0)
HParticleCutRange(TString name = "cut", Int_t num = -1, Double_t l = -1e30, Double_t u = 1e30)
void setInverse(Bool_t inv)
{ fbInverseCut = inv;}
void set(Double_t l = -1e30, Double_t u = 1e30)
{ flow = l; fup = u; }
TF1* getTF1Low()
{ return fLowF1; }
TF1* getTF1Up()
{ return fUpF1 ; }
TCutG* getTCutG()
{ return fCut; }
UInt_t getMaxVersion()
{ return fmaxVersion; }
Int_t getCutNumber()
{ return fCutNumber; }
Double_t getLowCut()
{ return flow;}
Double_t getUpCut()
{ return fup;}