ROOT logo
#include "hparticlecutrange.h"
#include <iostream>
#include <iomanip>

#include "TDirectory.h"
#include "TROOT.h"

using namespace std;

ClassImp(HParticleCutRange)

//_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();
//
//
//}
////////////////////////////////////////////////////////////////////////////////

HParticleCutRange::HParticleCutRange(TString name,
				     Int_t num,Double_t l,Double_t u)
{
    SetName(name);
    TDirectory* saveDir = gDirectory;
    gROOT->cd();
    fCut = new TCutG();
    fCut->SetName(name);
    fLowF1  = NULL;
    fUpF1   = NULL;
    fF1Mode = -1;
    flow        = l;
    fup         = u;
    fCutNumber  = num;
    fmaxVersion = 0;
    fbInverseCut= kFALSE;
    setMaxCut(4);
    saveDir->cd();
}

HParticleCutRange::~HParticleCutRange()
{
    ;
}


void HParticleCutRange::setTF1(TF1* flowF, TF1* fupF,TString mode)
{
    // 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


    fF1Mode= -1;
    fLowF1 = NULL;
    fUpF1  = NULL;

    if(mode.CompareTo("low+up") == 0){
	if(flowF == 0 || fupF == 0){
	    Error("setTF1()","Both TF1 needs to be !=0 in mode \"low+up\"");
	    return;
	} else {
	    fF1Mode = 0;
	    fLowF1  = flowF;
	    fUpF1   = fupF;
	}
    } else if(mode.CompareTo("low") == 0){
	if(flowF == 0 || fupF != 0){
	    Error("setTF1()","Lower TF1 needs to be set and upper TF1 == 0  in mode \"low\"");
	    return;

	} else {
            fF1Mode=1;
	    fLowF1 = flowF;
	}
    } else if(mode.CompareTo("up") == 0){
	if(flowF != 0 || fupF == 0){
	    Error("setTF1()","Upper TF1 needs to be set and lower TF1 == 0  in mode \"up\"");
	    return;

	} else {
            fF1Mode=2;
	    fUpF1 = fupF;
	}
    } else if(mode.CompareTo("mean+width") == 0){
	if(flowF == 0 || fupF == 0){
	    Error("setTF1()","Both TF1 needs to be !=0 in mode \"mean+width\"");
	    return;
	} else {
            fF1Mode=3;
	    fLowF1 = flowF;
	    fUpF1 = fupF;
	}
    } else {
	Error("setTF1()","Unknown mode =%s",mode.Data());
    }
}

void   HParticleCutRange::setMaxCut(UInt_t max)
{
    // set the max number of versions of stst counters
    fmaxCut     = max;
    if(fmaxCut == 0) fmaxCut = 1;
    resetCounter();
}

Bool_t HParticleCutRange::eval(Double_t var,UInt_t version){

    // returns kTRUE if the object fullfills the condition.
    // counters for calls and failed conditions are filled.
    // version is used to fill the stat for a given version
    // (0 up to 9 incl )
    if(version<fmaxCut){
	fctCall[version]++;
        if(version>fmaxVersion) fmaxVersion = version;
    }
    if( (!fbInverseCut &&  (var > fup || var < flow)) ||
        ( fbInverseCut && !(var > fup || var < flow))
      )
       {
	if(version<fmaxCut)fctFail[version]++;
	return kFALSE;
    } else return kTRUE;
}

Bool_t HParticleCutRange::evalG(Double_t var,Double_t var2,UInt_t version){

    // returns kTRUE if the object fullfills the condition.
    // counters for calls and failed conditions are filled.
    // version is used to fill the stat for a given version

    if(fCut->GetN() == 0) {
	Warning("evalG(var,var)","Cut %s has no set TCutG !",GetName());
        return kFALSE;
    }
    if(version<fmaxCut){
	fctCall[version]++;
        if(version>fmaxVersion) fmaxVersion = version;
    }
    if( (!fbInverseCut && !fCut->IsInside(var,var2)) ||
        ( fbInverseCut &&  fCut->IsInside(var,var2))
      )
       {
	if(version<fmaxCut)fctFail[version]++;
	return kFALSE;
    } else return kTRUE;
}

Bool_t HParticleCutRange::evalF(Double_t var,Double_t var2,UInt_t version){

    // returns kTRUE if the object fullfills the condition.
    // counters for calls and failed conditions are filled.
    // version is used to fill the stat for a given version

    if(fF1Mode == -1) {
	Warning("evalF(var,var)","TF1 Cut %s has no set !",GetName());
        return kFALSE;
    }

    if(version<fmaxCut){
	fctCall[version]++;
        if(version>fmaxVersion) fmaxVersion = version;
    }

    Double_t lowCut = -1e30;
    Double_t upCut  =  1e30;

    if(fF1Mode == 0) {    // low+up
	lowCut = fLowF1->Eval(var2);
	upCut  = fUpF1 ->Eval(var2);
    } else if (fF1Mode == 1 ) { // low
	lowCut = fLowF1->Eval(var2);
    } else if (fF1Mode == 2) { // up
	upCut = fUpF1->Eval(var2);
    } else {  // mean+width
	Double_t width = fUpF1->Eval(var2);
	Double_t mean  = fLowF1->Eval(var2);
	lowCut = mean - width;
	upCut  = mean + width;
    }


    if( (!fbInverseCut && !( var<upCut && var>lowCut) )  ||
        ( fbInverseCut &&  ( var<upCut && var>lowCut))
      )
       {
	if(version<fmaxCut)fctFail[version]++;
	return kFALSE;
    } else return kTRUE;
}

void HParticleCutRange::print()
{
    Int_t p = cout.precision();
    std::ios_base::fmtflags fl =cout.flags();
    cout<<"CutNumber : "<<setw(3)<<fCutNumber
	<<" name : "<<setw(20)<<GetName()
	<<", lowCut  : "<<setw(10)<<flow
	<<", UpCut  : "<<setw(10)<<fup<<", cut [%] : "<<flush;
    for(UInt_t i = 0; i < fmaxCut ; i++){ cout<<setw(5)<<fixed<<right<<setprecision(1)<<getCutRate(i)<<" "<<flush;}
    cout<<setprecision(p)<<endl;
    cout.flags(fl);
}

void HParticleCutRange::resetCounter()
{
    fctFail.clear();
    fctCall.clear();
    fctFail.assign(fmaxCut,0);
    fctCall.assign(fmaxCut,0);
}

Float_t HParticleCutRange::getCutRate(UInt_t version)
{
    return (version <fmaxCut && fctCall[version] > 0) ? (fctFail[version]/(Float_t)fctCall[version])*100.:0;
}

UInt_t HParticleCutRange::getNCall(UInt_t version)
{
    return (version <fmaxCut ) ? fctCall[version]:0;
}

UInt_t HParticleCutRange::getNFail(UInt_t version)
{
    return (version <fmaxCut ) ? fctFail[version]:0;
}

 hparticlecutrange.cc:1
 hparticlecutrange.cc:2
 hparticlecutrange.cc:3
 hparticlecutrange.cc:4
 hparticlecutrange.cc:5
 hparticlecutrange.cc:6
 hparticlecutrange.cc:7
 hparticlecutrange.cc:8
 hparticlecutrange.cc:9
 hparticlecutrange.cc:10
 hparticlecutrange.cc:11
 hparticlecutrange.cc:12
 hparticlecutrange.cc:13
 hparticlecutrange.cc:14
 hparticlecutrange.cc:15
 hparticlecutrange.cc:16
 hparticlecutrange.cc:17
 hparticlecutrange.cc:18
 hparticlecutrange.cc:19
 hparticlecutrange.cc:20
 hparticlecutrange.cc:21
 hparticlecutrange.cc:22
 hparticlecutrange.cc:23
 hparticlecutrange.cc:24
 hparticlecutrange.cc:25
 hparticlecutrange.cc:26
 hparticlecutrange.cc:27
 hparticlecutrange.cc:28
 hparticlecutrange.cc:29
 hparticlecutrange.cc:30
 hparticlecutrange.cc:31
 hparticlecutrange.cc:32
 hparticlecutrange.cc:33
 hparticlecutrange.cc:34
 hparticlecutrange.cc:35
 hparticlecutrange.cc:36
 hparticlecutrange.cc:37
 hparticlecutrange.cc:38
 hparticlecutrange.cc:39
 hparticlecutrange.cc:40
 hparticlecutrange.cc:41
 hparticlecutrange.cc:42
 hparticlecutrange.cc:43
 hparticlecutrange.cc:44
 hparticlecutrange.cc:45
 hparticlecutrange.cc:46
 hparticlecutrange.cc:47
 hparticlecutrange.cc:48
 hparticlecutrange.cc:49
 hparticlecutrange.cc:50
 hparticlecutrange.cc:51
 hparticlecutrange.cc:52
 hparticlecutrange.cc:53
 hparticlecutrange.cc:54
 hparticlecutrange.cc:55
 hparticlecutrange.cc:56
 hparticlecutrange.cc:57
 hparticlecutrange.cc:58
 hparticlecutrange.cc:59
 hparticlecutrange.cc:60
 hparticlecutrange.cc:61
 hparticlecutrange.cc:62
 hparticlecutrange.cc:63
 hparticlecutrange.cc:64
 hparticlecutrange.cc:65
 hparticlecutrange.cc:66
 hparticlecutrange.cc:67
 hparticlecutrange.cc:68
 hparticlecutrange.cc:69
 hparticlecutrange.cc:70
 hparticlecutrange.cc:71
 hparticlecutrange.cc:72
 hparticlecutrange.cc:73
 hparticlecutrange.cc:74
 hparticlecutrange.cc:75
 hparticlecutrange.cc:76
 hparticlecutrange.cc:77
 hparticlecutrange.cc:78
 hparticlecutrange.cc:79
 hparticlecutrange.cc:80
 hparticlecutrange.cc:81
 hparticlecutrange.cc:82
 hparticlecutrange.cc:83
 hparticlecutrange.cc:84
 hparticlecutrange.cc:85
 hparticlecutrange.cc:86
 hparticlecutrange.cc:87
 hparticlecutrange.cc:88
 hparticlecutrange.cc:89
 hparticlecutrange.cc:90
 hparticlecutrange.cc:91
 hparticlecutrange.cc:92
 hparticlecutrange.cc:93
 hparticlecutrange.cc:94
 hparticlecutrange.cc:95
 hparticlecutrange.cc:96
 hparticlecutrange.cc:97
 hparticlecutrange.cc:98
 hparticlecutrange.cc:99
 hparticlecutrange.cc:100
 hparticlecutrange.cc:101
 hparticlecutrange.cc:102
 hparticlecutrange.cc:103
 hparticlecutrange.cc:104
 hparticlecutrange.cc:105
 hparticlecutrange.cc:106
 hparticlecutrange.cc:107
 hparticlecutrange.cc:108
 hparticlecutrange.cc:109
 hparticlecutrange.cc:110
 hparticlecutrange.cc:111
 hparticlecutrange.cc:112
 hparticlecutrange.cc:113
 hparticlecutrange.cc:114
 hparticlecutrange.cc:115
 hparticlecutrange.cc:116
 hparticlecutrange.cc:117
 hparticlecutrange.cc:118
 hparticlecutrange.cc:119
 hparticlecutrange.cc:120
 hparticlecutrange.cc:121
 hparticlecutrange.cc:122
 hparticlecutrange.cc:123
 hparticlecutrange.cc:124
 hparticlecutrange.cc:125
 hparticlecutrange.cc:126
 hparticlecutrange.cc:127
 hparticlecutrange.cc:128
 hparticlecutrange.cc:129
 hparticlecutrange.cc:130
 hparticlecutrange.cc:131
 hparticlecutrange.cc:132
 hparticlecutrange.cc:133
 hparticlecutrange.cc:134
 hparticlecutrange.cc:135
 hparticlecutrange.cc:136
 hparticlecutrange.cc:137
 hparticlecutrange.cc:138
 hparticlecutrange.cc:139
 hparticlecutrange.cc:140
 hparticlecutrange.cc:141
 hparticlecutrange.cc:142
 hparticlecutrange.cc:143
 hparticlecutrange.cc:144
 hparticlecutrange.cc:145
 hparticlecutrange.cc:146
 hparticlecutrange.cc:147
 hparticlecutrange.cc:148
 hparticlecutrange.cc:149
 hparticlecutrange.cc:150
 hparticlecutrange.cc:151
 hparticlecutrange.cc:152
 hparticlecutrange.cc:153
 hparticlecutrange.cc:154
 hparticlecutrange.cc:155
 hparticlecutrange.cc:156
 hparticlecutrange.cc:157
 hparticlecutrange.cc:158
 hparticlecutrange.cc:159
 hparticlecutrange.cc:160
 hparticlecutrange.cc:161
 hparticlecutrange.cc:162
 hparticlecutrange.cc:163
 hparticlecutrange.cc:164
 hparticlecutrange.cc:165
 hparticlecutrange.cc:166
 hparticlecutrange.cc:167
 hparticlecutrange.cc:168
 hparticlecutrange.cc:169
 hparticlecutrange.cc:170
 hparticlecutrange.cc:171
 hparticlecutrange.cc:172
 hparticlecutrange.cc:173
 hparticlecutrange.cc:174
 hparticlecutrange.cc:175
 hparticlecutrange.cc:176
 hparticlecutrange.cc:177
 hparticlecutrange.cc:178
 hparticlecutrange.cc:179
 hparticlecutrange.cc:180
 hparticlecutrange.cc:181
 hparticlecutrange.cc:182
 hparticlecutrange.cc:183
 hparticlecutrange.cc:184
 hparticlecutrange.cc:185
 hparticlecutrange.cc:186
 hparticlecutrange.cc:187
 hparticlecutrange.cc:188
 hparticlecutrange.cc:189
 hparticlecutrange.cc:190
 hparticlecutrange.cc:191
 hparticlecutrange.cc:192
 hparticlecutrange.cc:193
 hparticlecutrange.cc:194
 hparticlecutrange.cc:195
 hparticlecutrange.cc:196
 hparticlecutrange.cc:197
 hparticlecutrange.cc:198
 hparticlecutrange.cc:199
 hparticlecutrange.cc:200
 hparticlecutrange.cc:201
 hparticlecutrange.cc:202
 hparticlecutrange.cc:203
 hparticlecutrange.cc:204
 hparticlecutrange.cc:205
 hparticlecutrange.cc:206
 hparticlecutrange.cc:207
 hparticlecutrange.cc:208
 hparticlecutrange.cc:209
 hparticlecutrange.cc:210
 hparticlecutrange.cc:211
 hparticlecutrange.cc:212
 hparticlecutrange.cc:213
 hparticlecutrange.cc:214
 hparticlecutrange.cc:215
 hparticlecutrange.cc:216
 hparticlecutrange.cc:217
 hparticlecutrange.cc:218
 hparticlecutrange.cc:219
 hparticlecutrange.cc:220
 hparticlecutrange.cc:221
 hparticlecutrange.cc:222
 hparticlecutrange.cc:223
 hparticlecutrange.cc:224
 hparticlecutrange.cc:225
 hparticlecutrange.cc:226
 hparticlecutrange.cc:227
 hparticlecutrange.cc:228
 hparticlecutrange.cc:229
 hparticlecutrange.cc:230
 hparticlecutrange.cc:231
 hparticlecutrange.cc:232
 hparticlecutrange.cc:233
 hparticlecutrange.cc:234
 hparticlecutrange.cc:235
 hparticlecutrange.cc:236
 hparticlecutrange.cc:237
 hparticlecutrange.cc:238
 hparticlecutrange.cc:239
 hparticlecutrange.cc:240
 hparticlecutrange.cc:241
 hparticlecutrange.cc:242
 hparticlecutrange.cc:243
 hparticlecutrange.cc:244
 hparticlecutrange.cc:245
 hparticlecutrange.cc:246
 hparticlecutrange.cc:247
 hparticlecutrange.cc:248
 hparticlecutrange.cc:249
 hparticlecutrange.cc:250
 hparticlecutrange.cc:251
 hparticlecutrange.cc:252
 hparticlecutrange.cc:253
 hparticlecutrange.cc:254
 hparticlecutrange.cc:255
 hparticlecutrange.cc:256
 hparticlecutrange.cc:257
 hparticlecutrange.cc:258
 hparticlecutrange.cc:259
 hparticlecutrange.cc:260
 hparticlecutrange.cc:261
 hparticlecutrange.cc:262
 hparticlecutrange.cc:263
 hparticlecutrange.cc:264
 hparticlecutrange.cc:265
 hparticlecutrange.cc:266
 hparticlecutrange.cc:267
 hparticlecutrange.cc:268
 hparticlecutrange.cc:269
 hparticlecutrange.cc:270
 hparticlecutrange.cc:271
 hparticlecutrange.cc:272
 hparticlecutrange.cc:273
 hparticlecutrange.cc:274
 hparticlecutrange.cc:275
 hparticlecutrange.cc:276
 hparticlecutrange.cc:277
 hparticlecutrange.cc:278
 hparticlecutrange.cc:279
 hparticlecutrange.cc:280
 hparticlecutrange.cc:281
 hparticlecutrange.cc:282
 hparticlecutrange.cc:283
 hparticlecutrange.cc:284
 hparticlecutrange.cc:285
 hparticlecutrange.cc:286
 hparticlecutrange.cc:287
 hparticlecutrange.cc:288
 hparticlecutrange.cc:289
 hparticlecutrange.cc:290
 hparticlecutrange.cc:291
 hparticlecutrange.cc:292
 hparticlecutrange.cc:293
 hparticlecutrange.cc:294
 hparticlecutrange.cc:295
 hparticlecutrange.cc:296
 hparticlecutrange.cc:297
 hparticlecutrange.cc:298
 hparticlecutrange.cc:299
 hparticlecutrange.cc:300
 hparticlecutrange.cc:301
 hparticlecutrange.cc:302
 hparticlecutrange.cc:303
 hparticlecutrange.cc:304
 hparticlecutrange.cc:305
 hparticlecutrange.cc:306
 hparticlecutrange.cc:307
 hparticlecutrange.cc:308
 hparticlecutrange.cc:309
 hparticlecutrange.cc:310
 hparticlecutrange.cc:311
 hparticlecutrange.cc:312
 hparticlecutrange.cc:313
 hparticlecutrange.cc:314
 hparticlecutrange.cc:315
 hparticlecutrange.cc:316
 hparticlecutrange.cc:317
 hparticlecutrange.cc:318
 hparticlecutrange.cc:319
 hparticlecutrange.cc:320
 hparticlecutrange.cc:321
 hparticlecutrange.cc:322
 hparticlecutrange.cc:323
 hparticlecutrange.cc:324
 hparticlecutrange.cc:325
 hparticlecutrange.cc:326
 hparticlecutrange.cc:327
 hparticlecutrange.cc:328
 hparticlecutrange.cc:329
 hparticlecutrange.cc:330
 hparticlecutrange.cc:331
 hparticlecutrange.cc:332
 hparticlecutrange.cc:333
 hparticlecutrange.cc:334
 hparticlecutrange.cc:335
 hparticlecutrange.cc:336
 hparticlecutrange.cc:337
 hparticlecutrange.cc:338
 hparticlecutrange.cc:339
 hparticlecutrange.cc:340
 hparticlecutrange.cc:341
 hparticlecutrange.cc:342
 hparticlecutrange.cc:343
 hparticlecutrange.cc:344
 hparticlecutrange.cc:345
 hparticlecutrange.cc:346
 hparticlecutrange.cc:347
 hparticlecutrange.cc:348
 hparticlecutrange.cc:349
 hparticlecutrange.cc:350
 hparticlecutrange.cc:351
 hparticlecutrange.cc:352
 hparticlecutrange.cc:353
 hparticlecutrange.cc:354
 hparticlecutrange.cc:355
 hparticlecutrange.cc:356
 hparticlecutrange.cc:357
 hparticlecutrange.cc:358
 hparticlecutrange.cc:359
 hparticlecutrange.cc:360
 hparticlecutrange.cc:361
 hparticlecutrange.cc:362
 hparticlecutrange.cc:363
 hparticlecutrange.cc:364
 hparticlecutrange.cc:365
 hparticlecutrange.cc:366
 hparticlecutrange.cc:367
 hparticlecutrange.cc:368
 hparticlecutrange.cc:369
 hparticlecutrange.cc:370
 hparticlecutrange.cc:371
 hparticlecutrange.cc:372
 hparticlecutrange.cc:373
 hparticlecutrange.cc:374
 hparticlecutrange.cc:375
 hparticlecutrange.cc:376
 hparticlecutrange.cc:377
 hparticlecutrange.cc:378
 hparticlecutrange.cc:379
 hparticlecutrange.cc:380
 hparticlecutrange.cc:381
 hparticlecutrange.cc:382
 hparticlecutrange.cc:383
 hparticlecutrange.cc:384
 hparticlecutrange.cc:385
 hparticlecutrange.cc:386
 hparticlecutrange.cc:387
 hparticlecutrange.cc:388
 hparticlecutrange.cc:389
 hparticlecutrange.cc:390
 hparticlecutrange.cc:391
 hparticlecutrange.cc:392
 hparticlecutrange.cc:393
 hparticlecutrange.cc:394
 hparticlecutrange.cc:395
 hparticlecutrange.cc:396
 hparticlecutrange.cc:397
 hparticlecutrange.cc:398
 hparticlecutrange.cc:399
 hparticlecutrange.cc:400
 hparticlecutrange.cc:401
 hparticlecutrange.cc:402