Trimmer Class Reference

#include <Trimmer.h>

List of all members.

Public Member Functions

 Trimmer ()
void AddFiles (string files)
void SetOutputFile (string name)
void Config (Registry &r)
void SetCuts (string type, int level=0)
bool EvaluateCuts (NueRecord *nr)
void RunTrimmer ()
void SetDeltaMSquare (double dm2)
void SetUe3Square (double dUe32)
void SetTheta23 (double t23)
void SetBaseline (double bl)
void RecalculatePOT ()
void SeparatebyRunPeriod ()
void SetRHC ()

Private Attributes

string outFile
bool outSet
vector< string > files
NueAnalysisCuts fCuts
double fDeltaMSquare
double fUe3Square
double fTheta23
double fBaseline
bool fReweight
string cutSet
int cutLevel
bool fOverWritePOT
bool separatebyRunPeriod
bool isRHC


Detailed Description

Definition at line 16 of file Trimmer.h.


Constructor & Destructor Documentation

Trimmer::Trimmer (  ) 

Definition at line 6 of file Trimmer.cxx.

References fBaseline, fDeltaMSquare, fOverWritePOT, fReweight, fTheta23, fUe3Square, isRHC, outSet, and separatebyRunPeriod.

00006                 { 
00007   outSet = false;
00008   fReweight = false;
00009   fBaseline = 735;
00010   fDeltaMSquare = 0.0027;
00011   fUe3Square = 0.025;
00012   fTheta23 = TMath::Pi()/4;
00013 
00014   fOverWritePOT = false;
00015   separatebyRunPeriod = false;
00016   isRHC = false;
00017 }


Member Function Documentation

void Trimmer::AddFiles ( string  files  ) 

Definition at line 19 of file Trimmer.cxx.

References files.

00020 {
00021    files.push_back(infiles);
00022 }

void Trimmer::Config ( Registry r  ) 

Definition at line 30 of file Trimmer.cxx.

References NueAnalysisCuts::Config(), and fCuts.

00031 {
00032   fCuts.Config(r);
00033 }

bool Trimmer::EvaluateCuts ( NueRecord nr  ) 

Definition at line 436 of file Trimmer.cxx.

References cutLevel, cutSet, NueStandard::GetPIDValue(), NueStandard::IsGoodRun(), NueStandard::IsInFid(), Selection::kANN2PE, Selection::kBasic, Selection::kCC, Selection::kDataQual, Selection::kFid, Selection::kMCNN, Selection::kPre, NueStandard::PassesMRCCFiducial(), NueStandard::PassesMRCCPreSelection(), NueStandard::PassesMREFiducial(), NueStandard::PassesNonHEPreSelection(), NueStandard::PassesPreSelectionTrackCuts(), NueStandard::PassesSelection(), NueStandard::PassesSysPreSelection(), and NueStandard::PassesSysPreSelectionNoHE().

Referenced by RunTrimmer().

00437 {
00438         //just for merging files!
00439    if(cutSet == "Merge")
00440         return true;
00441 
00442    if(cutSet == "Fiducial") 
00443      return NueStandard::IsInFid(nr);
00444 
00445    if(cutSet == "Standard"){
00446       bool ret;
00447       switch(cutLevel){
00448         case 0: return NueStandard::PassesSelection(nr, Selection::kDataQual);
00449         case 1: return NueStandard::PassesSelection(nr, Selection::kFid);
00450         case 2: return NueStandard::PassesSelection(nr, Selection::kBasic);
00451         case 3: 
00452              ret = NueStandard::PassesSelection(nr, Selection::kBasic);
00453              ret = ret && NueStandard::PassesPreSelectionTrackCuts(nr);
00454              return ret;
00455         case 4: 
00456              ret = NueStandard::PassesSelection(nr, Selection::kBasic);
00457              ret = ret &&  NueStandard::PassesNonHEPreSelection(nr);
00458              return ret;
00459         case 5: return NueStandard::PassesSelection(nr, Selection::kPre);
00460       }
00461    }
00462    if(cutSet == "MRE" || cutSet == "MRCC"){
00463       bool dq = NueStandard::PassesSelection(nr, Selection::kDataQual);
00464       bool mrefid = NueStandard::PassesMREFiducial(nr);
00465       bool ret;
00466       switch(cutLevel){
00467         case 0: return dq && NueStandard::IsInFid(nr);
00468         case 1: return mrefid && dq && NueStandard::IsInFid(nr);
00469         case 2: 
00470           return NueStandard::PassesSelection(nr, Selection::kFid);//includes mrcc presel & mrcc fiducial
00471         case 3: return NueStandard::PassesSelection(nr,Selection::kBasic);
00472         case 4:
00473              ret = NueStandard::PassesSelection(nr,Selection::kBasic);
00474              ret = ret && NueStandard::PassesPreSelectionTrackCuts(nr);
00475              return ret;
00476         case 5:
00477              ret = NueStandard::PassesSelection(nr, Selection::kBasic);
00478              ret = ret &&  NueStandard::PassesNonHEPreSelection(nr);
00479              return ret;
00480         case 6: 
00481          return NueStandard::PassesSelection(nr, Selection::kPre);
00482       }
00483    }
00484    if(cutSet == "Presel")
00485       return NueStandard::PassesSelection(nr, Selection::kPre);
00486 
00487    if(cutSet == "PreselNoHE")
00488       return  NueStandard::PassesNonHEPreSelection(nr);
00489 
00490    if(cutSet == "Systematic"){
00491       // O - standard systematic envelope
00492       // 1 - standard systematic w/out HE cut
00493       if(!NueStandard::IsInFid(nr)) return false;
00494 
00495       switch(cutLevel){
00496         case 0: return NueStandard::PassesSysPreSelection(nr);
00497         case 1: return NueStandard::PassesSysPreSelectionNoHE(nr);
00498       }
00499    }
00500 
00501 
00502    if(cutSet == "CC")
00503       return NueStandard::PassesSelection(nr, Selection::kCC);
00504    
00505    if(cutSet == "GoodRun")
00506      return NueStandard::IsGoodRun(nr);
00507    
00508    if(cutSet == "AntiANN11")
00509    {
00510      if(NueStandard::PassesSelection(nr, Selection::kPre) && NueStandard::GetPIDValue(nr,Selection::kANN2PE)<0.55)
00511      {
00512        return true;
00513      }
00514      else
00515      {
00516        return false;
00517      }
00518    }
00519    
00520    if(cutSet == "AntiMCNN")
00521    {
00522      if(NueStandard::PassesSelection(nr, Selection::kPre) && NueStandard::GetPIDValue(nr,Selection::kMCNN)<0.55)
00523      {
00524        return true;
00525      }
00526      else
00527      {
00528        return false;
00529      }
00530    }
00531    if(cutSet=="MRCuts")
00532    {
00533      if(NueStandard::PassesSelection(nr, Selection::kDataQual) && NueStandard::PassesMRCCFiducial(nr) && NueStandard::PassesMRCCPreSelection(nr))
00534      {
00535        return true;
00536      }
00537      return false;
00538    }
00539    
00540   cout<<"Invalid Cut Level for "<<cutSet<<": "<<cutLevel<<endl;
00541   return false;
00542 }

void Trimmer::RecalculatePOT (  )  [inline]

Definition at line 32 of file Trimmer.h.

References fOverWritePOT.

00032 {fOverWritePOT = true;};

void Trimmer::RunTrimmer (  ) 

Definition at line 40 of file Trimmer.cxx.

References base, NuePOT::beamtype, beamtype, count, det, EvaluateCuts(), fDeltaMSquare, files, fOverWritePOT, fReweight, fTheta23, fUe3Square, NueHeader::GetBeamType(), VldContext::GetDetector(), RecRecordImp< T >::GetHeader(), NueStandard::GetIntensityBeamWeight(), MCInfo::GetMCPoT(), NueStandard::GetNDDataWeights(), NueHeader::GetRelease(), RecHeader::GetVldContext(), DataUtil::IsGoodData(), isMC, isRHC, NueStandard::IsRun1(), NueStandard::IsRun2(), NueStandard::IsRun3(), Detector::kFar, SimFlag::kMC, Detector::kNear, NueStandard::ModifyANNPID(), NuePOT::nruns, NuePOT::nsnarls, NueConvention::NueEnergyCorrection(), NueConvention::Oscillate(), outFile, outSet, NueStandard::PassesPOTStandards(), NuePOT::pot, NuePOT::pot_nocut, NueConvention::RHCNueEnergyCorrection(), and separatebyRunPeriod.

00041 {
00042   
00043   NueRecord *nr = new NueRecord();
00044   NuePOT *np = new NuePOT();
00045 
00046   TChain *chain = new TChain("ana_nue");
00047   chain->SetBranchAddress("NueRecord",&nr);
00048                                                                         
00049   TChain *pchain = new TChain("pottree");
00050   pchain->SetBranchAddress("NuePOT", &np);
00051  
00052   for(unsigned int i = 0; i < files.size(); i++)
00053   {
00054       chain->Add(files[i].c_str());
00055       pchain->Add(files[i].c_str());
00056   }
00057   
00058   if(!outSet){
00059      string file;
00060                                                                         
00061      if(pchain->GetEntries() > 0){
00062         pchain->GetEntry(0);
00063         file = pchain->GetFile()->GetName();
00064      }
00065      string minifile = file.substr(file.find_last_of("/")+1, file.find_last_of(".root")-file.find_last_of("/") - 5);
00066      minifile += "-Trim.root";
00067   
00068      outFile = minifile;
00069   }
00070   if(outFile == "-Trim.root"){
00071      cout<<"No input file found"<<endl; 
00072      return;
00073   }
00074   
00075   string out1,out2,out3;
00076   if(separatebyRunPeriod)
00077   {
00078     string base = outFile.substr(0,outFile.size()-5);//strip off '.root' and add label for run period
00079     out1 = base + "_RunPeriod1.root";
00080     out2 = base + "_RunPeriod2.root";
00081     out3 = base + "_RunPeriod3.root";
00082     cout<<"Setting output to:"<<endl;
00083     cout<<out1<<endl;
00084     cout<<out2<<endl;
00085     cout<<out3<<endl;
00086   }
00087   else
00088   {
00089     cout<<"Setting output to "<<outFile<<endl;
00090     out1 = outFile;
00091   }
00092 
00093   // And now the actual looping
00094   
00095   TFile *f1 = new TFile(out1.c_str(),"RECREATE");
00096   TFile *f2=0,*f3=0;
00097   if(separatebyRunPeriod)
00098   {
00099     f2 = new TFile(out2.c_str(),"RECREATE");
00100     f3 = new TFile(out3.c_str(),"RECREATE");
00101   }
00102   
00103   f1->cd();
00104   TTree *tree1 = new TTree("ana_nue","ana_nue");
00105   TTree::SetBranchStyle(1);
00106   TBranch* br1 = tree1->Branch("NueRecord","NueRecord", &nr );
00107   br1->SetAutoDelete(kFALSE);
00108   
00109   TTree *tree2=0;
00110   TTree *tree3=0;
00111   TBranch* br2=0;
00112   TBranch* br3=0;
00113   if(separatebyRunPeriod)
00114   {
00115     f2->cd();
00116     tree2 = new TTree("ana_nue","ana_nue");
00117     br2 = tree2->Branch("NueRecord","NueRecord", &nr );
00118     br2->SetAutoDelete(kFALSE);
00119     
00120     f3->cd();
00121     tree3 = new TTree("ana_nue","ana_nue");
00122     br3 = tree3->Branch("NueRecord","NueRecord", &nr );
00123     br3->SetAutoDelete(kFALSE);
00124   }
00125   
00126   bool isRunPeriod1=false;
00127   bool isRunPeriod2=false;
00128   bool isRunPeriod3=false;
00129   
00130   Int_t n = chain->GetEntries();
00131   int count = 0;
00132 
00133   bool isMC = false;
00134   
00135   //for recalculating the pot tree for each run period
00136   double TotPOT1 = 0.0,TotPOT2 = 0.0,TotPOT3 = 0.0;
00137   double TotPOT1_nocut = 0.0,TotPOT2_nocut = 0.0,TotPOT3_nocut = 0.0;
00138   Int_t nruns=0;
00139   Int_t nsnarls1=0,nsnarls2=0,nsnarls3=0;
00140   
00141   int lastrun=-1;
00142 
00143   double potweight;
00144   double nddataweight;
00145   
00146   chain->GetEntry(0);
00147   ReleaseType::Release_t rel = nr->GetHeader().GetRelease();
00148   BeamType::BeamType_t beamtype = nr->GetHeader().GetBeamType();
00149   Detector::Detector_t det = nr->GetHeader().GetVldContext().GetDetector();
00150   
00151   for(int i=0;i<n;i++){
00152     if(i%10000==0) cout << 100*float(i)/float(n)
00153                        << "% done" << endl;
00154     chain->GetEvent(i);
00155     
00156     NueConvention::NueEnergyCorrection(nr);
00157     if(isRHC) NueConvention::RHCNueEnergyCorrection(nr);
00158     
00159     nr->anainfo.isRHC = isRHC;
00160     
00161     NueStandard::ModifyANNPID(nr);
00162     
00163     nr->eventq.passFarDetQuality = DataUtil::IsGoodData(nr->GetHeader().GetVldContext());
00164     nr->eventq.passNearDetQuality = nr->eventq.passFarDetQuality;
00165     
00166     if(i == 0){
00167       SimFlag::SimFlag_t sim = nr->GetHeader().GetVldContext().GetSimFlag();
00168       isMC = (sim == SimFlag::kMC);
00169     }
00170     
00171     if(separatebyRunPeriod && !isMC)
00172     {
00173       cout<<"'separatebyRunPeriod' option is not meant to be used for data - pottree would not be calculated correctly for data.  Qutting..."<<endl;
00174       return;
00175     }
00176     
00177     //sanity check
00178     if(NueStandard::PassesPOTStandards(nr))
00179     {
00180       if(nr->GetHeader().GetEventNo() == 0 || nr->GetHeader().GetEvents() == 0)
00181       {
00182         if(nr->bmon.bI < 0 && !isMC) cout<<"Unexpected behavior"<<endl;
00183       }
00184     }
00185    
00186     isRunPeriod1=false;
00187     isRunPeriod2=false;
00188     isRunPeriod3=false;
00189     if(separatebyRunPeriod)
00190     {
00191       if(NueStandard::IsRun1(nr))
00192       {
00193         isRunPeriod1=true;
00194       }
00195       if(NueStandard::IsRun2(nr))
00196       {
00197         isRunPeriod2=true;
00198       }
00199       if(NueStandard::IsRun3(nr))
00200       {
00201         isRunPeriod3=true;
00202       }
00203     }
00204     
00205     if(separatebyRunPeriod)//recalculating pottree; this part works only for MC
00206     {
00207       if(nr->GetHeader().GetEventNo() == 0 || nr->GetHeader().GetEvents() == 0)//for each snarl
00208       {
00209         if(isRunPeriod1)
00210         {
00211           if(det==Detector::kNear)//this only works for ND MC
00212           {
00213             TotPOT1_nocut += MCInfo::GetMCPoT(det, beamtype, rel);
00214             TotPOT1 += MCInfo::GetMCPoT(det, beamtype, rel);
00215           }
00216           
00217           nsnarls1++;//count the snarls per run period
00218         }
00219         if(isRunPeriod2)
00220         {
00221           if(det==Detector::kNear)//this only works for ND MC
00222           {
00223             TotPOT2_nocut += MCInfo::GetMCPoT(det, beamtype, rel);
00224             TotPOT2 += MCInfo::GetMCPoT(det, beamtype, rel);
00225           }
00226           
00227           nsnarls2++;
00228         }
00229         if(isRunPeriod3)
00230         {
00231           if(det==Detector::kNear)//this only works for ND MC
00232           {
00233             potweight = NueStandard::GetIntensityBeamWeight(nr);
00234             TotPOT3_nocut += (potweight*MCInfo::GetMCPoT(det, beamtype, rel));
00235             TotPOT3 += (potweight*MCInfo::GetMCPoT(det, beamtype, rel));
00236           }
00237           
00238           nsnarls3++;
00239         }
00240       }
00241       
00242       if(nr->GetHeader().GetRun()!=lastrun)//counting runs - BUT not counting separately for each run period!!
00243       {
00244         lastrun = nr->GetHeader().GetRun();
00245         nruns++;
00246         
00247         if(det==Detector::kFar)//in FD MC, each run is a single file, and GetMCPoT returns pot/file
00248         {
00249           if(isRunPeriod1)
00250           {
00251             TotPOT1_nocut += MCInfo::GetMCPoT(det, beamtype, rel);
00252             TotPOT1 += MCInfo::GetMCPoT(det, beamtype, rel);
00253           }
00254           if(isRunPeriod2)
00255           {
00256             TotPOT2_nocut += MCInfo::GetMCPoT(det, beamtype, rel);
00257             TotPOT2 += MCInfo::GetMCPoT(det, beamtype, rel);
00258           }
00259           if(isRunPeriod3)
00260           {
00261             TotPOT3_nocut += MCInfo::GetMCPoT(det, beamtype, rel);
00262             TotPOT3 += MCInfo::GetMCPoT(det, beamtype, rel);
00263           }
00264         }
00265       }
00266     }
00267     
00268     if(fOverWritePOT && !isMC)//recalculating POT for data
00269     {
00270       if(NueStandard::PassesPOTStandards(nr))
00271       {
00272         if(nr->GetHeader().GetEventNo() == 0 || nr->GetHeader().GetEvents() == 0)
00273         {
00274           if(det==Detector::kNear)
00275           {
00276             nddataweight = NueStandard::GetNDDataWeights(nr);
00277             TotPOT1 += (nr->bmon.bI*nddataweight);
00278           }
00279           else//FD
00280           {
00281             TotPOT1 += nr->bmon.bI;
00282           }
00283         }
00284       }
00285     }
00286     
00287     if(EvaluateCuts(nr)){
00288        if(fReweight && nr->mctrue.nuFlavor > -10)
00289        {
00290           int nuFlavor = nr->mctrue.nuFlavor;
00291           int  nonOsc = nr->mctrue.nonOscNuFlavor;
00292           float energy = nr->mctrue.nuEnergy;
00293                                                                                 
00294           Float_t newWeight = NueConvention::Oscillate(nuFlavor, nonOsc, energy,735, fDeltaMSquare, fTheta23, fUe3Square);
00295                                                                                 
00296           nr->mctrue.Ue3Squared = fUe3Square;
00297           nr->mctrue.DeltamSquared23 = fDeltaMSquare;
00298           nr->mctrue.Theta23  = fTheta23;
00299           nr->mctrue.fOscProb = newWeight;
00300        }
00301        
00302        if(separatebyRunPeriod)
00303        {
00304          if(isRunPeriod1)
00305          {
00306            tree1->Fill();
00307          }
00308          if(isRunPeriod2)
00309          {
00310            tree2->Fill();
00311          }
00312          if(isRunPeriod3)
00313          {
00314            tree3->Fill();
00315          }
00316        }
00317        else
00318        {
00319          tree1->Fill();
00320        }
00321        count++;
00322     }
00323   }
00324   cout<<count<<" of "<<n<<" entries were passed"<<endl;
00325   
00326   f1->cd();
00327   NuePOT *total1 = new NuePOT();
00328   TTree *ptree1 = new TTree("pottree","pottree");
00329   TBranch* brp1 = ptree1->Branch("NuePOT","NuePOT", &total1 );
00330   brp1->SetAutoDelete(kFALSE);
00331   
00332   NuePOT *total2=0,*total3=0;
00333   TTree *ptree2=0,*ptree3=0;
00334   TBranch *brp2=0,*brp3=0;
00335   if(separatebyRunPeriod)
00336   {
00337     f2->cd();
00338     total2 = new NuePOT();
00339     ptree2 = new TTree("pottree","pottree");
00340     brp2 = ptree2->Branch("NuePOT","NuePOT", &total2 );
00341     brp2->SetAutoDelete(kFALSE);
00342     
00343     f3->cd();
00344     total3 = new NuePOT();
00345     ptree3 = new TTree("pottree","pottree");
00346     brp3 = ptree3->Branch("NuePOT","NuePOT", &total3 );
00347     brp3->SetAutoDelete(kFALSE);
00348   }
00349   
00350   if(pchain->GetEntries() > 0){
00351     pchain->GetEntry(0);
00352     total1->beamtype = np->beamtype;
00353     if(separatebyRunPeriod)
00354     {
00355       total2->beamtype = np->beamtype;
00356       total3->beamtype = np->beamtype;
00357     }
00358   }
00359   
00360   for(int i = 0; i < pchain->GetEntries(); i++)
00361   {
00362      pchain->GetEntry(i);
00363      
00364      if(total1->beamtype != np->beamtype)
00365      {
00366        cerr<<"You are merging files of different beamtype - BAD"<<endl;
00367      }
00368      
00369      total1->nruns +=  np->nruns;
00370      total1->nsnarls +=  np->nsnarls;
00371      total1->pot += np->pot;
00372      total1->pot_nocut += np->pot_nocut;
00373   }
00374 
00375   if(fOverWritePOT && !isMC) total1->pot = TotPOT1;
00376   
00377   if(separatebyRunPeriod)//recalculating pottree
00378   {
00379     total1->pot = TotPOT1;
00380     total2->pot = TotPOT2;
00381     total3->pot = TotPOT3;
00382     
00383     total1->pot_nocut = TotPOT1_nocut;
00384     total2->pot_nocut = TotPOT2_nocut;
00385     total3->pot_nocut = TotPOT3_nocut;
00386     
00387     total1->nruns = nruns;
00388     total2->nruns = nruns;
00389     total3->nruns = nruns;
00390     
00391     total1->nsnarls = nsnarls1;
00392     total2->nsnarls = nsnarls2;
00393     total3->nsnarls = nsnarls3;
00394   }
00395   
00396   //fill the pottree
00397   ptree1->Fill();
00398   if(separatebyRunPeriod)
00399   {
00400     ptree2->Fill();
00401     ptree3->Fill();
00402   }
00403   
00404   //save the file(s)
00405   f1->cd();
00406   tree1->Write("ana_nue",TObject::kWriteDelete);
00407   ptree1->Write("pottree",TObject::kWriteDelete);
00408   if(separatebyRunPeriod)
00409   {
00410     f2->cd();
00411     tree2->Write("ana_nue",TObject::kWriteDelete);
00412     ptree2->Write("pottree",TObject::kWriteDelete);
00413     
00414     f3->cd();
00415     tree3->Write("ana_nue",TObject::kWriteDelete);
00416     ptree3->Write("pottree",TObject::kWriteDelete);
00417   }
00418   
00419   f1->Close();
00420   if(separatebyRunPeriod)
00421   {
00422     f2->Close();
00423     f3->Close();
00424   }
00425   
00426   
00427   cout<<"Trimming completed with "<<np->pot<<"x10^12 POT"<<endl;
00428 }

void Trimmer::SeparatebyRunPeriod (  )  [inline]

Definition at line 33 of file Trimmer.h.

References separatebyRunPeriod.

00033 {separatebyRunPeriod = true;};

void Trimmer::SetBaseline ( double  bl  ) 

Definition at line 38 of file Trimmer.cxx.

References fBaseline, and fReweight.

00038 { fBaseline = bl; fReweight = true;}

void Trimmer::SetCuts ( string  type,
int  level = 0 
)

Definition at line 430 of file Trimmer.cxx.

References cutLevel, and cutSet.

00431 {
00432    cutSet = type;
00433    cutLevel = level;   
00434 }

void Trimmer::SetDeltaMSquare ( double  dm2  ) 

Definition at line 35 of file Trimmer.cxx.

References fDeltaMSquare, and fReweight.

00035 { fDeltaMSquare = dm2; fReweight = true;}

void Trimmer::SetOutputFile ( string  name  ) 

Definition at line 24 of file Trimmer.cxx.

References outFile, and outSet.

00025 {
00026    outSet = true;
00027    outFile = file;
00028 }

void Trimmer::SetRHC (  )  [inline]

Definition at line 34 of file Trimmer.h.

References isRHC.

00034 {isRHC = true;};

void Trimmer::SetTheta23 ( double  t23  ) 

Definition at line 37 of file Trimmer.cxx.

References fReweight, and fTheta23.

00037 { fTheta23 = t23; fReweight = true;}

void Trimmer::SetUe3Square ( double  dUe32  ) 

Definition at line 36 of file Trimmer.cxx.

References fReweight, and fUe3Square.

00036 { fUe3Square = dUe32; fReweight = true;}


Member Data Documentation

int Trimmer::cutLevel [private]

Definition at line 49 of file Trimmer.h.

Referenced by EvaluateCuts(), and SetCuts().

string Trimmer::cutSet [private]

Definition at line 48 of file Trimmer.h.

Referenced by EvaluateCuts(), and SetCuts().

double Trimmer::fBaseline [private]

Definition at line 45 of file Trimmer.h.

Referenced by SetBaseline(), and Trimmer().

NueAnalysisCuts Trimmer::fCuts [private]

Definition at line 40 of file Trimmer.h.

Referenced by Config().

double Trimmer::fDeltaMSquare [private]

Definition at line 42 of file Trimmer.h.

Referenced by RunTrimmer(), SetDeltaMSquare(), and Trimmer().

vector<string> Trimmer::files [private]

Definition at line 39 of file Trimmer.h.

Referenced by AddFiles(), and RunTrimmer().

bool Trimmer::fOverWritePOT [private]

Definition at line 51 of file Trimmer.h.

Referenced by RecalculatePOT(), RunTrimmer(), and Trimmer().

bool Trimmer::fReweight [private]

Definition at line 46 of file Trimmer.h.

Referenced by RunTrimmer(), SetBaseline(), SetDeltaMSquare(), SetTheta23(), SetUe3Square(), and Trimmer().

double Trimmer::fTheta23 [private]

Definition at line 44 of file Trimmer.h.

Referenced by RunTrimmer(), SetTheta23(), and Trimmer().

double Trimmer::fUe3Square [private]

Definition at line 43 of file Trimmer.h.

Referenced by RunTrimmer(), SetUe3Square(), and Trimmer().

bool Trimmer::isRHC [private]

Definition at line 53 of file Trimmer.h.

Referenced by RunTrimmer(), SetRHC(), and Trimmer().

string Trimmer::outFile [private]

Definition at line 34 of file Trimmer.h.

Referenced by RunTrimmer(), and SetOutputFile().

bool Trimmer::outSet [private]

Definition at line 38 of file Trimmer.h.

Referenced by RunTrimmer(), SetOutputFile(), and Trimmer().

bool Trimmer::separatebyRunPeriod [private]

Definition at line 52 of file Trimmer.h.

Referenced by RunTrimmer(), SeparatebyRunPeriod(), and Trimmer().


The documentation for this class was generated from the following files:
Generated on Thu Apr 10 23:03:51 2014 for loon by  doxygen 1.4.7