FEDRA emulsion software from the OPERA Collaboration
EdbShowerAlgESimple Class Reference

#include <EdbShowerAlg.h>

Inheritance diagram for EdbShowerAlgESimple:
Collaboration diagram for EdbShowerAlgESimple:

Public Member Functions

 ClassDef (EdbShowerAlgESimple, 1)
 
void CreateANN ()
 
void DoRun ()
 
void DoRun (EdbTrackP *shower)
 
void DoRun (TObjArray *trackarray)
 
void DumpNeuralNetworkWeight (TString weight, Int_t ANNType=0)
 
 EdbShowerAlgESimple ()
 
 EdbShowerAlgESimple (EdbTrackP *track)
 in libShower More...
 
 EdbShowerAlgESimple (TObjArray *RecoShowerArray)
 
Int_t FindClosestEfficiencyParametrization (Double_t TestAngle, Double_t ReferenceEff)
 
TF1 * GetEffFunc_all ()
 
TF1 * GetEffFunc_edefault ()
 
TF1 * GetEffFunc_elletroni ()
 
TF1 * GetEffFunc_LowEff ()
 
TF1 * GetEffFunc_MiddleFix ()
 
TF1 * GetEffFunc_neuchmicro ()
 
TF1 * GetEffFunc_UserEfficiency ()
 
TSpline3 * GetEfficiencyParametrisation ()
 
Float_t GetEnergy ()
 
Float_t GetEnergy (EdbTrackP *track)
 
TArrayF * GetEnergyArray ()
 
TArrayF * GetEnergyArrayUnCorrected ()
 
TMultiLayerPerceptron * GetNeuralNetwork (Int_t ANNType=0)
 
void GetNplIndexNr (Int_t sh_npl, Int_t &check_Npl_index, Int_t ePlateNumberType)
 
void GetPara (EdbTrackP *track)
 
Int_t GetRecoShowerArrayN () const
 
void GetSpecifications ()
 
Int_t GetSpecType (Int_t SpecificationType)
 
void Help ()
 
void InitStrings ()
 
void LoadSpecificationWeightFile ()
 
void Print ()
 
void PrintEfficiencyParametrisation ()
 
void PrintSpecifications ()
 
void ReadCorrectionFactors (TString weigthstring, Float_t &p0, Float_t &p1)
 
void ReadTables ()
 
void ReadTables_Energy ()
 
void SetCalibrationOffset (Float_t CalibrationOffset)
 
void SetCalibrationSlope (Float_t CalibrationSlope)
 
void SetEfficiencyParametrisation (TSpline3 *EfficiencyParametrisation)
 
void SetEfficiencyParametrisationAngles ()
 
void SetEfficiencyParametrisationValues (Double_t *Angles, Double_t *EffValuesAtAngles)
 
void SetForceSpecificationReload ()
 
void SetPlateNumber (Int_t PlateNumber)
 
void SetPlateNumberType (Int_t PlateNumberType)
 
void SetRecoShowerArray (TObjArray *RecoShowerArray)
 
void SetRecoShowerArrayN (Int_t RecoShowerArrayN)
 
void SetSpecifications (Int_t sp0, Int_t sp1, Int_t sp2, Int_t sp3, Int_t sp4, Int_t sp5)
 
void SetSpecificationType (Int_t SpecificationType, Int_t SpecificationTypeVal)
 
void SetWeightFileString (TString weightstring)
 
void TrainNeuralNetwork (TString weight, Int_t ANNType=0)
 
void UnSetForceSpecificationReload ()
 
void Update ()
 
void WriteNewRootFile ()
 
void WriteNewRootFile (TString sourcefilename)
 
void WriteNewRootFile (TString sourcefilename, TString treename)
 
virtual ~EdbShowerAlgESimple ()
 virtual constructor due to inherited class More...
 

Protected Member Functions

void Init ()
 
void Set0 ()
 

Protected Attributes

TString ANN_Layout
 
TMultiLayerPerceptron * ANN_MLP
 
TMultiLayerPerceptron * ANN_MLP_ARRAY [15]
 
Int_t ANN_n_InputNeurons
 
Int_t ANN_n_InputNeurons_ARRAY [15]
 
Int_t ANN_nPlates_ARRAY [15]
 
TString ANN_WeightFile_ARRAY [15]
 
TTree * ANNTree
 
Float_t eANN_MLP_CORR_0 [15]
 
Float_t eANN_MLP_CORR_1 [15]
 
Float_t eCalibrationOffset
 
Float_t eCalibrationSlope
 
TSpline3 * eEfficiencyParametrisation
 
Float_t eEnergy
 
TArrayF * eEnergyArray
 
Int_t eEnergyArrayCount
 
TArrayF * eEnergyArraySigmaCorrected
 
TArrayF * eEnergyArrayUnCorrected
 
Float_t eEnergyCorr
 
Float_t eEnergySigmaCorr
 
Float_t eEnergyUnCorr
 
TF1 * EffFunc_all
 
TF1 * EffFunc_edefault
 
TF1 * EffFunc_elletroni
 
TF1 * EffFunc_LowEff
 
TF1 * EffFunc_MiddleFix
 
TF1 * EffFunc_neuchmicro
 
TF1 * EffFunc_UserEfficiency
 
Bool_t eForceSpecificationReload
 
TH1D * eHisto_deltaR
 
TH1D * eHisto_deltaR_mean
 
TH1D * eHisto_deltaR_rms
 
TH1D * eHisto_deltaT
 
TH1D * eHisto_deltaT_mean
 
TH1D * eHisto_deltaT_rms
 
TH1D * eHisto_longprofile
 
TH1D * eHisto_longprofile_av
 
TH1D * eHisto_nbtk
 
TH1D * eHisto_nbtk_av
 
TH1D * eHisto_transprofile
 
TH1D * eHisto_transprofile_av
 
Float_t eParaBT_deltaR_mean
 
Float_t eParaBT_deltaR_rms
 
Float_t eParaBT_deltaT_mean
 
Float_t eParaBT_deltaT_rms
 
Int_t eParalongprofile [57]
 
Int_t eParaName
 
Int_t eParanseg
 
Float_t eParaShowerAxisAngle
 
Int_t ePlateBinning [15]
 
Int_t ePlateNumber
 
Int_t ePlateNumberType
 
TObjArray * eRecoShowerArray
 
Int_t eRecoShowerArrayN
 
Bool_t eSpecificationIsChanged
 
Int_t eSpecificationType [7]
 
TString eSpecificationTypeString [7]
 
TString eSpecificationTypeStringArray [7][7]
 
TObjArray * eSplineArray_Energy_Stat_Electron
 
TObjArray * eSplineArray_Energy_Stat_Gamma
 
TObjArray * eSplineArray_Energy_Sys_Electron
 
TObjArray * eSplineArray_Energy_Sys_Gamma
 
TSpline3 * eSplineCurrent
 
TString eWeightFileString
 
Double_t inANN [70]
 
Double_t outANN
 

Detailed Description



___ Declaration of EdbShowerAlgESimple Class:


Constructor & Destructor Documentation

◆ EdbShowerAlgESimple() [1/3]

EdbShowerAlgESimple::EdbShowerAlgESimple ( )

Default Constructor

2451{
2453 if (gEDBDEBUGLEVEL >1) cout << "EdbShowerAlgESimple::EdbShowerAlgESimple() Default Constructor"<<endl;
2454 Set0();
2455 Init();
2456}
void Init()
Definition: EdbShowerAlg.cxx:2573
void Set0()
Definition: EdbShowerAlg.cxx:2529
gEDBDEBUGLEVEL
Definition: energy.C:7

◆ EdbShowerAlgESimple() [2/3]

EdbShowerAlgESimple::EdbShowerAlgESimple ( EdbTrackP track)

in libShower

CONSTRUCTOR USED FOR libShowRec with EdbShowerP class available.

Default Constructor...

2478{
2480 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::EdbShowerAlgESimple(EdbTrackP* track) Default Constructor"<<endl;
2481 Set0();
2482 Init();
2483 // ... with Shower:
2484 eRecoShowerArray=new TObjArray();
2485 eRecoShowerArray->Add(track);
2487 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::EdbShowerAlgESimple(EdbTrackP* track) Default Constructor...done."<<endl;
2488}
void SetRecoShowerArrayN(Int_t RecoShowerArrayN)
Definition: EdbShowerAlg.h:452
TObjArray * eRecoShowerArray
Definition: EdbShowerAlg.h:332
Definition: bitview.h:14

◆ EdbShowerAlgESimple() [3/3]

EdbShowerAlgESimple::EdbShowerAlgESimple ( TObjArray *  RecoShowerArray)

Default Constructor

2494{
2496 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::EdbShowerAlgESimple(TObjArray* RecoShowerArray) Default Constructor"<<endl;
2497 Set0();
2498 Init();
2499 // ... with ShowerArray:
2501 eRecoShowerArrayN=eRecoShowerArray->GetEntries();
2502 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::EdbShowerAlgESimple(TObjArray* RecoShowerArray) Default Constructor...done."<<endl;
2503}
TObjArray * RecoShowerArray
Definition: Shower_E_FromShowerRoot.C:12
Int_t eRecoShowerArrayN
Definition: EdbShowerAlg.h:333

◆ ~EdbShowerAlgESimple()

EdbShowerAlgESimple::~EdbShowerAlgESimple ( )
virtual

virtual constructor due to inherited class

Default Destructor

2508{
2510 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::~EdbShowerAlgESimple()"<<endl;
2511 // Delete Histograms (on heap):
2512 delete eHisto_nbtk_av;
2513 delete eHisto_longprofile_av;
2515 delete eHisto_deltaR_mean;
2516 delete eHisto_deltaT_mean;
2517 delete eHisto_deltaR_rms;
2518 delete eHisto_deltaT_rms;
2519 delete eHisto_nbtk;
2520 delete eHisto_longprofile;
2521 delete eHisto_transprofile;
2522 delete eHisto_deltaR;
2523 delete eHisto_deltaT;
2524 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::~EdbShowerAlgESimple()...done."<<endl;
2525}
TH1D * eHisto_longprofile
Definition: EdbShowerAlg.h:395
TH1D * eHisto_deltaT
Definition: EdbShowerAlg.h:398
TH1D * eHisto_nbtk
Definition: EdbShowerAlg.h:394
TH1D * eHisto_deltaR
Definition: EdbShowerAlg.h:397
TH1D * eHisto_transprofile_av
Definition: EdbShowerAlg.h:389
TH1D * eHisto_transprofile
Definition: EdbShowerAlg.h:396
TH1D * eHisto_nbtk_av
Definition: EdbShowerAlg.h:387
TH1D * eHisto_deltaT_rms
Definition: EdbShowerAlg.h:393
TH1D * eHisto_longprofile_av
Definition: EdbShowerAlg.h:388
TH1D * eHisto_deltaT_mean
Definition: EdbShowerAlg.h:391
TH1D * eHisto_deltaR_rms
Definition: EdbShowerAlg.h:392
TH1D * eHisto_deltaR_mean
Definition: EdbShowerAlg.h:390

Member Function Documentation

◆ ClassDef()

EdbShowerAlgESimple::ClassDef ( EdbShowerAlgESimple  ,
 
)

◆ CreateANN()

void EdbShowerAlgESimple::CreateANN ( )

Create ANN Tree and MLP:

2692{
2693 if (gEDBDEBUGLEVEL >1) cout << "EdbShowerAlgESimple::CreateANN()"<<endl;
2694
2696 eParaName=2; // Standard Labeling taken over from EdbShowerP
2697 ANNTree = new TTree("ANNTree", "ANNTree");
2698 ANNTree->SetDirectory(0);
2699 ANNTree->Branch("inANN", inANN, "inANN[70]/D");
2700
2701 // Plate Binning: 10,12,14,16,18,20,23,26,29,32,35,40,45,45
2702 // see ePlateBinning[]
2703
2704 for (int k=0; k<15; k++) {
2705 //cout << "creatin ANN " << k << endl;
2706 ANN_Layout="";
2709 for (Int_t i=1; i<ANN_n_InputNeurons; ++i) ANN_Layout += "@inANN["+TString(Form("%d",i))+"],";
2710 ANN_Layout += "@inANN["+TString(Form("%d",ANN_n_InputNeurons))+"]:"+TString(Form("%d",ANN_n_InputNeurons+1))+":"+TString(Form("%d",ANN_n_InputNeurons));
2711 ANN_Layout+=":inANN[0]/1000";
2712 //---------------------------
2713 ANN_MLP_ARRAY[k] = new TMultiLayerPerceptron(ANN_Layout,ANNTree,"","");
2715 ANN_Layout=ANN_MLP_ARRAY[k]->GetStructure();
2716 //---------------------------
2717 if (gEDBDEBUGLEVEL>2) {
2718 ANN_MLP->Print();
2719 cout << ANN_Layout.Data() << endl;
2720 }
2721 }
2722 //---------------------------
2723
2724 if (gEDBDEBUGLEVEL >1) cout << "EdbShowerAlgESimple::CreateANN()...done."<<endl;
2725 return;
2726}
TString ANN_Layout
Definition: EdbShowerAlg.h:369
Double_t inANN[70]
Definition: EdbShowerAlg.h:364
TTree * ANNTree
Definition: EdbShowerAlg.h:363
Int_t ePlateBinning[15]
Definition: EdbShowerAlg.h:358
Int_t ANN_n_InputNeurons
Definition: EdbShowerAlg.h:368
Int_t ANN_n_InputNeurons_ARRAY[15]
Definition: EdbShowerAlg.h:367
TMultiLayerPerceptron * ANN_MLP_ARRAY[15]
Definition: EdbShowerAlg.h:362
Int_t eParaName
Definition: EdbShowerAlg.h:334
TMultiLayerPerceptron * ANN_MLP
Definition: EdbShowerAlg.h:361

◆ DoRun() [1/3]

void EdbShowerAlgESimple::DoRun ( )
2731{
2732 if (gEDBDEBUGLEVEL >1) cout << "EdbShowerAlgESimple::DoRun()" << endl;
2733
2734 for (int i=0; i<eRecoShowerArrayN; i++) {
2735 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() Doing i= " << i << endl;
2736 continue;
2738
2739 //EdbShowerP* shower=(EdbShowerP*) eRecoShowerArray->At(i);
2740 EdbTrackP* shower=(EdbTrackP*) eRecoShowerArray->At(i);
2741
2742 DoRun(shower);
2743
2744 } // (int i=0; i<eRecoShowerArrayN; i++)
2745
2746 if (gEDBDEBUGLEVEL >1) cout << "EdbShowerAlgESimple::DoRun()...done." << endl;
2747 return;
2748}
Int_t eEnergyArrayCount
Definition: EdbShowerAlg.h:423
void DoRun()
Definition: EdbShowerAlg.cxx:2730
Definition: EdbPattern.h:113

◆ DoRun() [2/3]

void EdbShowerAlgESimple::DoRun ( EdbTrackP shower)

Basic run function for energy estimation.
Assumes, that Neural Networks are properly created, weightfiles correctly loaded,
and scanning efficiencies correctly evaluated.

  • Takes a shower, gets its energy parametrisation (i.e. the shower profile)
  • Takes the right ANN, run it with inputvariables from the shower.
  • According to the energy output value, the systematic uncertainties are read from
    a table and then those values are written in the shower storage, or in the
    "treebranch tree" root file.
2781{
2790
2791 if (gEDBDEBUGLEVEL>2) cout << "EdbShowerAlgESimple::DoRun(EdbTrackP* shower)" << endl;
2792
2793 // Get internal variables from the parametrisation filled:
2794 GetPara(shower);
2795
2796
2798
2799 // Check If Efficiency is the one we have or if we have to change/reload
2800 // the ANN weightfiles:
2801 if (GetSpecType(2)!=EffNr&&eForceSpecificationReload==kFALSE) {
2802 cout << "EdbShowerAlgESimple::DoRun() INFO! Calulated Efficiency is more compatible with another one: Change Specifiaction! Call SetSpecificationType(2,"<< EffNr <<")." << endl;
2803 SetSpecificationType(2,EffNr);
2804 }
2805 if (GetSpecType(2)!=EffNr&&eForceSpecificationReload==kTRUE) {
2806 cout << "EdbShowerAlgESimple::DoRun() INFO! Calulated Efficiency would be more compatible with another one: But eForceSpecificationReload is set to kFALSE. So we keep current Efficiency used." << endl;
2807
2808 }
2809
2810
2811
2812
2813 // This is to select the suited ANN to the shower, i.e. the one that matches closest
2814 // the number of plates:
2815 int check_Npl_index =0;
2816 GetNplIndexNr(shower->Npl(),check_Npl_index,ePlateNumberType);
2817
2818
2820 ANN_MLP=ANN_MLP_ARRAY[check_Npl_index];
2821 eWeightFileString=ANN_WeightFile_ARRAY[check_Npl_index];
2822 eCalibrationOffset = eANN_MLP_CORR_0[check_Npl_index];
2823 eCalibrationSlope = eANN_MLP_CORR_1[check_Npl_index];
2824
2825 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() ANN_n_InputNeurons_ARRAY[check_Npl_index]="<< ANN_n_InputNeurons_ARRAY[check_Npl_index] << endl;
2826 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() ANN_MLP_ARRAY[check_Npl_index]="<< ANN_MLP_ARRAY[check_Npl_index] << endl;
2827 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() Using the following layout: " << endl;
2828 if (gEDBDEBUGLEVEL >2) cout << ANN_MLP->GetStructure().Data() << endl;
2829 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() And the following weightfile: " << endl;
2830 if (gEDBDEBUGLEVEL >2) cout << eWeightFileString.Data() << endl;
2831
2832 // Reset InputVariables:
2833 for (int k=0; k<70; k++) {
2834 inANN[k]=0;
2835 }
2836
2837 // inANN[0] is ALWAYS Reseverd for the quantity value to be estimated
2838 // (E,Id,...)
2839 // Test with private variables:
2841 inANN[2]=eParanseg;
2846 for (int ii=0; ii<57; ii++) {
2847 inANN[7+ii]= eParalongprofile[ii];
2848 }
2849
2850 // Fill Tree:
2851 ANNTree->Fill();
2852
2853 Double_t params[70];
2854 Double_t val;
2855
2856 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() Print Inputvariables: " << endl;
2857 for (Int_t j=1; j<=ANN_n_InputNeurons; ++j) {
2858 params[j-1]=inANN[j];
2859 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() : j params[j-1]=inANN[j] " << j<< " " << params[j-1] << endl;
2860 }
2861 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() Print Inputvariables...done." << endl;
2862 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() Evaluate Neural Network Output:" << endl;
2863
2864 // ---------------------------------
2865 // Evaluation of the ANN:
2866 val=(ANN_MLP->Evaluate(0,params));
2867 // ---------------------------------
2868
2869 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() Before Correction: ANN_MLP->Evaluate(0,params)= " << ANN_MLP->Evaluate(0,params) << " (Inputvar=" << inANN[0] << ")." << endl;
2870
2871 if ( ::isnan(val) ) {
2872 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() ERROR! ANN_MLP->Evaluate(0,params) is NAN! Setting value to -1. " << endl;
2873 val=-1;
2874 }
2875 eEnergyUnCorr=val;
2876
2877 // Linear Correction According to the Parameters matching the right weightfile.
2879
2880 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() After Correction: ANN_MLP->Evaluate(0,params)= " << val << " (Inputvar=" << inANN[0] << ")." << endl;
2881
2882 // ---------------------------------
2883 // Evaluation of the statistical andsystematical errors by lookup tables and splines:
2884 // For further info and more details, I will report the methology in the thesis, and at
2885 // some next collaboration meetings.. (As of october 2010,fwm.)
2886 TSpline3* spline;
2887 Float_t val_sigma_stat, val_sigma_sys, val_sigma_tot;
2888 cout << "TSpline3* spline=eSplineArray_Energy_Stat_Gamma[check_Npl_index];" << endl;
2889 spline=(TSpline3*)eSplineArray_Energy_Stat_Gamma->At(check_Npl_index);
2890 cout << "Float_t val_sigma_stat=spline->Eval(val); " << spline->Eval(val) << endl;
2891 val_sigma_stat=spline->Eval(val);
2892
2893 cout << "TSpline3* spline=eSplineArray_Energy_Sys_Gamma[check_Npl_index];" << endl;
2894 spline=(TSpline3*)eSplineArray_Energy_Sys_Gamma->At(check_Npl_index);
2895 cout << "Float_t val_sigma_sys=spline->Eval(val); " << spline->Eval(val) << endl;
2896 val_sigma_sys=spline->Eval(val);
2897
2898 // Addition of errors: stat(+)sys. (+)=quadratic addition.
2899 val_sigma_tot=TMath::Sqrt(val_sigma_stat*val_sigma_stat+val_sigma_sys*val_sigma_sys);
2900 cout << "Float_t val_sigma_tot=... " << val_sigma_tot << endl;
2901
2902 cout << "Doing only statistics and one source of systematicsat the moment! " << endl;
2903 cout << "Notice also that we dont have storage variable in EdbTrackP for the error of P() ... " << endl;
2904 // ---------------------------------
2905
2906 // Quick Estimation of the sigma (from ANN_MEGA_ENERGY)
2907 // This is later to be read from a lookup table....
2908 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() ERROR! Sigma ONLY FROM STATISTICAL UNCERTATINTY NOWNOT correctly set up to now."<<endl;
2909 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() ERROR! Sigma will be fully implemented when lookup tables for all"<<endl;
2910 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() ERROR! systematic uncertainties are availible!"<<endl;
2911 if (gEDBDEBUGLEVEL >1) cout << "EdbShowerAlgESimple::DoRun() Estimated Energy = " << val << " +- " << val_sigma_tot << "..." << endl;
2912 // ---------------------------------
2913
2914
2915 // Finally, set values...
2916 shower->SetP(val);
2917 eEnergy=val;
2918 eEnergyCorr=val;
2919 eEnergySigmaCorr=val_sigma_stat;
2920
2924
2925
2926 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() ...Done." << endl;
2927 return;
2928}
void SetP(float p)
Definition: EdbSegP.h:133
Float_t eANN_MLP_CORR_0[15]
Definition: EdbShowerAlg.h:373
void GetPara(EdbTrackP *track)
Definition: EdbShowerAlg.cxx:3301
void GetNplIndexNr(Int_t sh_npl, Int_t &check_Npl_index, Int_t ePlateNumberType)
Definition: EdbShowerAlg.cxx:2993
Bool_t eForceSpecificationReload
Definition: EdbShowerAlg.h:354
Float_t eParaShowerAxisAngle
Definition: EdbShowerAlg.h:379
Float_t eParaBT_deltaR_mean
Definition: EdbShowerAlg.h:381
TArrayF * eEnergyArray
Definition: EdbShowerAlg.h:415
Int_t ePlateNumberType
Definition: EdbShowerAlg.h:340
TObjArray * eSplineArray_Energy_Stat_Gamma
Definition: EdbShowerAlg.h:427
Float_t eCalibrationOffset
Definition: EdbShowerAlg.h:337
Float_t eParaBT_deltaT_mean
Definition: EdbShowerAlg.h:383
Float_t eEnergySigmaCorr
Definition: EdbShowerAlg.h:421
TObjArray * eSplineArray_Energy_Sys_Gamma
Definition: EdbShowerAlg.h:429
Float_t eEnergy
Definition: EdbShowerAlg.h:418
Float_t eANN_MLP_CORR_1[15]
Definition: EdbShowerAlg.h:374
TSpline3 * eEfficiencyParametrisation
Definition: EdbShowerAlg.h:412
Int_t FindClosestEfficiencyParametrization(Double_t TestAngle, Double_t ReferenceEff)
Definition: EdbShowerAlg.cxx:2933
Float_t eCalibrationSlope
Definition: EdbShowerAlg.h:338
Int_t eParanseg
Definition: EdbShowerAlg.h:380
Float_t eParaBT_deltaT_rms
Definition: EdbShowerAlg.h:384
void SetSpecificationType(Int_t SpecificationType, Int_t SpecificationTypeVal)
Definition: EdbShowerAlg.cxx:3160
TArrayF * eEnergyArraySigmaCorrected
Definition: EdbShowerAlg.h:417
Float_t eEnergyCorr
Definition: EdbShowerAlg.h:419
Float_t eParaBT_deltaR_rms
Definition: EdbShowerAlg.h:382
Int_t eParalongprofile[57]
Definition: EdbShowerAlg.h:385
TArrayF * eEnergyArrayUnCorrected
Definition: EdbShowerAlg.h:416
Float_t eEnergyUnCorr
Definition: EdbShowerAlg.h:420
TString ANN_WeightFile_ARRAY[15]
Definition: EdbShowerAlg.h:370
TString eWeightFileString
Definition: EdbShowerAlg.h:346
Int_t GetSpecType(Int_t SpecificationType)
Definition: EdbShowerAlg.h:527
Int_t Npl() const
Definition: EdbPattern.h:171
Double_t params[3]
Definition: testBGReduction_By_ANN.C:84

◆ DoRun() [3/3]

void EdbShowerAlgESimple::DoRun ( TObjArray *  trackarray)

Runs over the object array of showers.

2754{
2756
2757 if (gEDBDEBUGLEVEL >1) cout << "EdbShowerAlgESimple::DoRun(TObjArray* trackarray)" << endl;
2758
2759 eRecoShowerArrayN=trackarray->GetEntries();
2760 eRecoShowerArray=trackarray;
2761
2762 for (int i=0; i<eRecoShowerArrayN; i++) {
2763 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::DoRun() Doing i= " << i << endl;
2764
2766
2767 //EdbShowerP* shower=(EdbShowerP*) eRecoShowerArray->At(i);
2768 EdbTrackP* shower=(EdbTrackP*) eRecoShowerArray->At(i);
2769
2770 DoRun(shower);
2771
2772 } // (int i=0; i<eRecoShowerArrayN; i++)
2773
2774 if (gEDBDEBUGLEVEL >1) cout << "EdbShowerAlgESimple::DoRun(TObjArray* trackarray)...done." << endl;
2775 return;
2776}

◆ DumpNeuralNetworkWeight()

void EdbShowerAlgESimple::DumpNeuralNetworkWeight ( TString  weight,
Int_t  ANNType = 0 
)
inline
508 {
509 if (ANNType>=15) ANNType=14;
510 ANN_MLP_ARRAY[ANNType]->DumpWeights(weight);
511 return;
512 }

◆ FindClosestEfficiencyParametrization()

Int_t EdbShowerAlgESimple::FindClosestEfficiencyParametrization ( Double_t  TestAngle = 0.0,
Double_t  ReferenceEff = 1.0 
)

Returns the number of the closest EfficiencyParametrization to use. Supported are:
((in this order numbering, different from the ordering, fwm used for his calculations!:))
0:Neuch
1:All
2:MiddleFix
3:Low

4: UserEfficiency (can take input from eda, for example; or user measured input) // this cannot be used here,
because the weightfiles exist only for Neuch,All,MiddleFix,Low. UserEfficiency is just for different
angle<->efficiency assignment!

2934{
2941
2945
2946 if (gEDBDEBUGLEVEL >2) cout << "TestAngle " << TestAngle << endl;
2947 if (gEDBDEBUGLEVEL >2) cout << "ReferenceEff " << ReferenceEff << endl;
2948 if (gEDBDEBUGLEVEL >2) cout << "neuchmicro->Eval(TestAngle) "<< EffFunc_neuchmicro->Eval(TestAngle)<< endl;
2949 if (gEDBDEBUGLEVEL >2) cout << "all->Eval(TestAngle) "<< EffFunc_all->Eval(TestAngle)<< endl;
2950 if (gEDBDEBUGLEVEL >2) cout << "MiddleFix->Eval(TestAngle) "<< EffFunc_MiddleFix->Eval(TestAngle)<< endl;
2951 if (gEDBDEBUGLEVEL >2) cout << "LowEff->Eval(TestAngle) "<< EffFunc_LowEff->Eval(TestAngle)<< endl;
2952 if (gEDBDEBUGLEVEL >2) cout << "eEfficiencyParametrisation->Eval(TestAngle) "<< eEfficiencyParametrisation->Eval(TestAngle)<< endl;
2953
2954//Measure distance of angle to estimated angle:
2955 Double_t dist[9];
2956 dist[0]=TMath::Abs(ReferenceEff-EffFunc_neuchmicro->Eval(TestAngle));
2957 dist[1]=TMath::Abs(ReferenceEff-EffFunc_all->Eval(TestAngle));
2958 dist[2]=TMath::Abs(ReferenceEff-EffFunc_MiddleFix->Eval(TestAngle));
2959 dist[3]=TMath::Abs(ReferenceEff-EffFunc_LowEff->Eval(TestAngle));
2960 dist[4]=TMath::Abs(ReferenceEff-eEfficiencyParametrisation->Eval(TestAngle));
2961
2962 Double_t dist_min=1;
2963 Int_t best_i=-1;
2964 for (int i=0; i<4; i++) if (abs(dist[i])<dist_min) {
2965 dist_min=dist[i];
2966 best_i=i;
2967 }
2968
2969 if (gEDBDEBUGLEVEL >2) {
2970 cout << "Miminum distance = " << dist_min << ", best matching efficiency is = ";
2971 if (best_i==0) cout << " neuchmicro " <<endl;
2972 if (best_i==2) cout << " MiddleFix " <<endl;
2973 if (best_i==3) cout << " LowEff " <<endl;
2974 if (best_i==1) cout << " All " <<endl;
2975 //if (best_i==4) cout << " User Set eEfficiencyParametrisation " <<endl;
2976 }
2977
2978 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::FindClosestEfficiencyParametrization() ...Done." << endl;
2979 return best_i;
2980}
TF1 * EffFunc_MiddleFix
Definition: EdbShowerAlg.h:405
TF1 * EffFunc_neuchmicro
Definition: EdbShowerAlg.h:404
TF1 * EffFunc_all
Definition: EdbShowerAlg.h:401
TF1 * EffFunc_LowEff
Definition: EdbShowerAlg.h:406

◆ GetEffFunc_all()

TF1 * EdbShowerAlgESimple::GetEffFunc_all ( )
inline
475 {
476 return EffFunc_all;
477 }

◆ GetEffFunc_edefault()

TF1 * EdbShowerAlgESimple::GetEffFunc_edefault ( )
inline
478 {
479 return EffFunc_edefault;
480 }
TF1 * EffFunc_edefault
Definition: EdbShowerAlg.h:402

◆ GetEffFunc_elletroni()

TF1 * EdbShowerAlgESimple::GetEffFunc_elletroni ( )
inline
481 {
482 return EffFunc_elletroni;
483 }
TF1 * EffFunc_elletroni
Definition: EdbShowerAlg.h:403

◆ GetEffFunc_LowEff()

TF1 * EdbShowerAlgESimple::GetEffFunc_LowEff ( )
inline
490 {
491 return EffFunc_LowEff;
492 }

◆ GetEffFunc_MiddleFix()

TF1 * EdbShowerAlgESimple::GetEffFunc_MiddleFix ( )
inline
487 {
488 return EffFunc_MiddleFix;
489 }

◆ GetEffFunc_neuchmicro()

TF1 * EdbShowerAlgESimple::GetEffFunc_neuchmicro ( )
inline
484 {
485 return EffFunc_neuchmicro;
486 }

◆ GetEffFunc_UserEfficiency()

TF1 * EdbShowerAlgESimple::GetEffFunc_UserEfficiency ( )
inline
493 {
495 }
TF1 * EffFunc_UserEfficiency
Definition: EdbShowerAlg.h:407

◆ GetEfficiencyParametrisation()

TSpline3 * EdbShowerAlgESimple::GetEfficiencyParametrisation ( )
inline
543 {
545 }

◆ GetEnergy() [1/2]

Float_t EdbShowerAlgESimple::GetEnergy ( )
inline
464 {
465 return eEnergy;
466 }

◆ GetEnergy() [2/2]

Float_t EdbShowerAlgESimple::GetEnergy ( EdbTrackP track)
inline
467 {
468 return track->P();
469 }

◆ GetEnergyArray()

TArrayF * EdbShowerAlgESimple::GetEnergyArray ( )
inline
458 {
459 return eEnergyArray;
460 }

◆ GetEnergyArrayUnCorrected()

TArrayF * EdbShowerAlgESimple::GetEnergyArrayUnCorrected ( )
inline
461 {
463 }

◆ GetNeuralNetwork()

TMultiLayerPerceptron * EdbShowerAlgESimple::GetNeuralNetwork ( Int_t  ANNType = 0)
inline
497 {
498 if (ANNType>=15) ANNType=14;
499 return ANN_MLP_ARRAY[ANNType];
500 }

◆ GetNplIndexNr()

void EdbShowerAlgESimple::GetNplIndexNr ( Int_t  sh_npl,
Int_t &  check_Npl_index,
Int_t  ePlateNumberType 
)

— What is the sense of this?
— According to the shower it can have any different length between [1..57] plates.
— We cannot produce weights for each plate....
— So we have only some weights for different plates: acoording to the array set in plateBinning:
— 10,12,14,16,18,20,23,26,29,32,35,40,45,50
— So what if a shower has npl() of 19? take the weight for 18 or for 20 plates?
— therefore, we introduced ePlateNumberType. This will take the weight as like:
— In case ePlateNumberType > 0 the it will take he trained ANN file with Npl >= shower->Npl();
— in such a way that |ePlateNumberType| is the difference to the ePlateBinning[].
— Example: ePlateNumberType=1 -> 1st next availible weight will be taken: (npl==19->nplANN==20);
— Example: ePlateNumberType=2 -> 2nd next availible weight will be taken: (npl==19->nplANN==23);
— Example: ePlateNumberType=3 -> 3rd next availible weight will be taken: (npl==19->nplANN==26);
— Example: ePlateNumberType=-1-> 1st befo availible weight will be taken: (npl==19->nplANN==18);
— Example: ePlateNumberType=-2-> 2nd befo availible weight will be taken: (npl==19->nplANN==16);

2994{
3009
3010 if (gEDBDEBUGLEVEL >1) cout << "EdbShowerAlgESimple::GetNplIndexNr("<<check_Npl<<","<<check_Npl_index<<","<<ePlateNumberType<<")"<<endl;
3011
3012 // The Neural Network Training was done in a way that we had always Npl(Train)>=shower->Npl()
3013
3014//F. Brunet fix : use the max plate number NN for shower extending within more than 45 plates
3015 if(check_Npl >=45)check_Npl = 44;
3016
3017 // First we find closest (bigger) ePlateBinning index:
3018 int closestIndex=-1;
3019 for (int i=0; i<15; i++) {
3020 if (ePlateBinning[i]-check_Npl>0) {
3021 closestIndex=i;
3022 break;
3023 }
3024 }
3025 cout << "closestIndex= " << closestIndex << " ePlateBinning[closestIndex]= " << ePlateBinning[closestIndex] << endl;
3026 //-----------------------------------
3027 // for ePlateNumberType we have it already exact:
3028 ePlateNumber = closestIndex + ePlateNumberType;
3029 if (ePlateNumberType==1) ePlateNumber = closestIndex;
3030 cout << "ePlateNumber = closestIndex + ePlateNumberType = " << ePlateNumber << endl;
3031
3032 // check for over/underbound
3033 if (ePlateNumber<0) ePlateNumber=0;
3034 if (ePlateNumber>14) ePlateNumber=14;
3035 check_Npl_index=ePlateNumber;
3036 //-----------------------------------
3037// // This Plate Binnning is the _final_binning and will not be refined furthermore!
3038// ePlateBinning[0]=10;// ePlateBinning[1]=12;
3039// ePlateBinning[2]=14;// ePlateBinning[3]=16;
3040// ePlateBinning[4]=18;// ePlateBinning[5]=20;
3041// ePlateBinning[6]=23;// ePlateBinning[7]=26;
3042// ePlateBinning[8]=29;// ePlateBinning[9]=32;
3043// ePlateBinning[10]=35;// ePlateBinning[11]=40;
3044// ePlateBinning[12]=45;// ePlateBinning[13]=45;
3045// ePlateBinning[14]=45;
3046 //-----------------------------------
3047 cout << "All plates:"<<endl;
3048 cout << "Shower Npl:"<<endl;
3049 cout << "Available weights:"<<endl;
3050 cout << "Taken weight:"<<endl;
3051 cout << "[";
3052 for (int i=1; i<=57; i++) {
3053 cout << ".";
3054 }
3055 cout << "]" << endl;
3056 //-----------------------------------
3057 cout << "[";
3058 for (int i=1; i<=check_Npl; i++) {
3059 cout << ".";
3060 }
3061 cout << "]" << endl;
3062 //-----------------------------------
3063 cout << "[";
3064 for (int i=1; i<=57; i++) {
3065 Bool_t isExact=kFALSE;
3066 for (int j=0; j<15; j++) {
3067 if (i==ePlateBinning[j]) isExact=kTRUE;
3068 }
3069 if (isExact) {
3070 cout << "x";
3071 }
3072 else {
3073 cout << ".";
3074 }
3075 }
3076 cout << "]" << endl;
3077 //-----------------------------------
3078 cout << "[";
3079 for (int i=1; i<=57; i++) {
3080 if (i==ePlateBinning[ePlateNumber]) {
3081 cout << "X";
3082 }
3083 else {
3084 cout << ".";
3085 }
3086 }
3087 cout << "]" << endl;
3088 //-----------------------------------
3089
3090 if (gEDBDEBUGLEVEL >1) cout << "EdbShowerAlgESimple::check_Npl_index Shower:Npl()=" << check_Npl << "."<< endl;
3091 if (gEDBDEBUGLEVEL >1) cout << "EdbShowerAlgESimple::check_Npl_index ePlateNumber=" << ePlateNumber << "."<< endl;
3092 if (gEDBDEBUGLEVEL >1) cout << "EdbShowerAlgESimple::check_Npl_index ePlateBinning[=" << ePlateBinning[ePlateNumber] << "."<< endl;
3093
3094
3095 check_Npl=ePlateBinning[ePlateNumber];
3096 return;
3097}
Int_t ePlateNumber
Definition: EdbShowerAlg.h:343

◆ GetPara()

void EdbShowerAlgESimple::GetPara ( EdbTrackP track)

This parametrisation consists of these variables:
0) Axis TanTheta
1) NBT
2) BT_deltaR_mean
3) BT_deltaR_rms
4) BT_deltaT_mean
5) BT_deltaT_rms
6) longprofile[0]=eHisto_longprofile->GetBinContent(1) // number of basetracks in the SAME plate as the Initiator basetrack.
7) longprofile[1]=eHisto_longprofile->GetBinContent(2)
...
5+Npl()) longprofile[Npl()-1]=eHisto_longprofile->GetBinContent(Npl()) // number of basetracks in the LAST plate of the reconstructed shower.
==C== DUMMY routine to get the values deltarb and deltathetab filled:
==C== necessary since its an old relict from old times where this was saved only
==C== as root TTree.

--— IMPORTANT:: these areopen cut values for best combifinding of pair BT deltaR/Theta values --— IMPORTANT:: then you do NOT necessarily get back your values which you put in durign --— IMPORTANT:: your shower reconstruction cone ( deltaR/Theta cutvalues could be NO cutvalues --— IMPORTANT:: for some reconstruction algorithms for example, but we wanna have these values anyway. In Any Case: Frederics Cut looks only for best min_shower_deltar so we do also.

3302{
3303 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::GetPara Fill the Para structure with values from a track." << endl;
3304
3320
3321 EdbSegP* seg;
3322 Float_t shower_xb[5000];
3323 Float_t shower_yb[5000];
3324 Float_t shower_zb[5000];
3325 Float_t shower_txb[5000];
3326 Float_t shower_tyb[5000];
3327 Float_t shower_deltathetab[5000];
3328 Float_t shower_deltarb[5000];
3329 Float_t min_shower_deltathetab=99999; // Reset
3330 Float_t min_shower_deltar=99999; // Reset
3331 Float_t extrapo_diffz=0;
3332 Float_t extrapol_x=0;
3333 Float_t extrapol_y=0;
3334
3335 Float_t test_shower_deltathetab;
3336 Float_t test_shower_deltar;
3337 Float_t test_shower_deltax;
3338 Float_t test_shower_deltay;
3339
3340 Int_t TrackNSeg=track->N();
3341
3342 for (int ii=0; ii<TrackNSeg; ii++) {
3343 // the ii<N() => ii<TrackNSeg replacement is just an adaption from the
3344 // code snipplet taken from EdbShowRec
3345 seg=(EdbSegP*)track->GetSegment(ii);
3346 shower_xb[ii]=seg->X();
3347 shower_yb[ii]=seg->Y();
3348 shower_txb[ii]=seg->TX();
3349 shower_tyb[ii]=seg->TY();
3350 shower_zb[ii]=seg->Z();
3351 shower_deltathetab[ii]=0.5;
3352 shower_deltarb[ii]=200;
3353 }
3354
3355 //-------------------------------------
3356 for (int ii=0; ii<TrackNSeg; ii++) {
3357 seg=(EdbSegP*)track->GetSegment(ii);
3358 if (gEDBDEBUGLEVEL>2) {
3359 if (gEDBDEBUGLEVEL >2) cout << "====== --- DOING " << ii << endl;
3360 seg->PrintNice();
3361 }
3362 //-------------------------------------
3363 // InBT:
3364 if (ii==0) {
3365 shower_deltathetab[ii]=0.5;
3366 shower_deltarb[ii]=200;
3367 }
3368 // All other BTs:
3369 if (ii>0) {
3370
3371 // PUT HERE: calculation routine for shower_deltathetab, shower_deltarb
3372 // Exrapolate the BT [ii] to the position [jj] and then calc the
3373 // position and slope differences for the best matching next segment.
3374 // For the backward extrapolation of the shower_deltathetab and shower_deltarb
3375 // calulation for BaseTrack(ii), Basetrack(jj)->Z() hast to be smaller.
3376 min_shower_deltathetab=99999; // Reset
3377 min_shower_deltar=99999; // Reset
3378
3379 for (int jj=0; jj<track->N(); jj++) {
3380 if (ii==jj) continue;
3381
3382 // since we do not know if BTs are ordered by their Z positions:
3383 // and cannot cut directly on the number in the shower entry:
3384 if (shower_zb[ii]<shower_zb[jj]) continue;
3385
3386 extrapo_diffz=shower_zb[ii]-shower_zb[jj];
3387 if (TMath::Abs(extrapo_diffz)>4*1300+1.0) continue;
3388 if (TMath::Abs(extrapo_diffz)<1.0) continue; // remove same positions.
3389
3390 extrapol_x=shower_xb[ii]-shower_txb[ii]*extrapo_diffz; // minus, because its ii after jj.
3391 extrapol_y=shower_yb[ii]-shower_tyb[ii]*extrapo_diffz; // minus, because its ii after jj.
3392
3393 // Delta radius we need to extrapolate.
3394 test_shower_deltax=extrapol_x;//shower_txb[ii]*(shower_zb[ii]-shower_zb[jj])+shower_xb[ii];
3395 test_shower_deltay=extrapol_y;//shower_tyb[ii]*(shower_zb[ii]-shower_zb[jj])+shower_yb[ii];
3396 test_shower_deltax=test_shower_deltax-shower_xb[jj];
3397 test_shower_deltay=test_shower_deltay-shower_yb[jj];
3398 test_shower_deltar=TMath::Sqrt(test_shower_deltax*test_shower_deltax+test_shower_deltay*test_shower_deltay);
3399
3400 // Delta theta we do not need to extrapolate. (old version...)
3401 //test_shower_deltathetab=TMath::Sqrt(shower_txb[ii]*shower_txb[ii]+shower_tyb[ii]*shower_tyb[ii]);
3402 //test_shower_deltathetab=test_shower_deltathetab-TMath::Sqrt(shower_txb[jj]*shower_txb[jj]+shower_tyb[jj]*shower_tyb[jj]);
3403 //test_shower_deltathetab=TMath::Abs(test_shower_deltathetab);
3404 //----
3405 // As before in ShowRec this way of calculation is not equivalent as calculating
3406 // DeltaTheta domponentwise:
3407 // Code from libShower:
3408 // delta = sqrt((SX0-a->GetTXb(l2))*(SX0-a->GetTXb(l2))+((SY0-a->GetTYb(l2))*(SY0-a->GetTYb(l2))));
3409 test_shower_deltathetab=TMath::Sqrt(TMath::Power(shower_txb[ii]-shower_txb[jj],2)+TMath::Power(shower_tyb[ii]-shower_tyb[jj],2));
3410
3411 // Check if both dr,dt match parameter criteria and then just take these values.....
3412 // Maybe a change is necessary because it is not exactly the same as in the off. algorithm:
3413 if (test_shower_deltar<1000 && test_shower_deltathetab<2.0 ) { // open cut values
3414 // Make these values equal to the one in the "official algorithm"..... 150microns and 100mrad.
3415 //if (test_shower_deltar<150 && test_shower_deltathetab<0.15 ) { // frederics cut values
3422 if (test_shower_deltar<min_shower_deltar) {
3423 min_shower_deltathetab=test_shower_deltathetab;
3424 min_shower_deltar=test_shower_deltar;
3425 shower_deltathetab[ii]=min_shower_deltathetab;
3426 shower_deltarb[ii]=min_shower_deltar;
3427 }
3428 }
3429 }
3430 }
3431 //-------------------------------------
3432
3433 } // for (int ii=0;ii<N();ii++)
3434 //-------------------------------------
3435
3436 //-------------------------------------
3437 if (gEDBDEBUGLEVEL>2) {
3438 if (gEDBDEBUGLEVEL >2) cout << "EdbTrackP: N() nsegments= " << track->N() << endl;
3439 for (int ii=0; ii<track->N(); ii++) {
3440 if (gEDBDEBUGLEVEL >2) cout << "Shower: nentry= " << ii << " shower_zb[ii] = " << shower_zb[ii] << " shower_deltathetab[ii] = " << shower_deltathetab[ii] << " shower_deltarb[ii] = " << shower_deltarb[ii] <<endl;
3441 }
3442 }
3443 //-------------------------------------
3444
3445 // Profile starting with (arrayindex ==0): number of basetracks in the SAME plate as the Initiator basetrack.
3446 // Profile ending with (arrayindex ==56): // number of basetracks in the LAST plate of the reconstructed shower.
3447 Int_t longprofile[57];
3448
3449 //Int_t numberofilms=Npl(); // Historical reasons why there are more names for the same variable...
3450 //Int_t nbfilm=Npl(); // Historical reasons why there are more names for the same variable...
3451 Int_t nbfilm=track->Npl();
3452
3453 // Int_t numberoffilsforlongprofile=numberofilms+1; // it is one more going from [1..nbfilm] !!!// Not used in this routine.
3454
3455 Float_t Dr[57];
3456 Float_t X0[57];
3457 Float_t Y0[57];
3458 Float_t TX0,TY0;
3459 Float_t theta[57];
3460 Float_t dist;
3461 Float_t xb[5000];
3462 Float_t yb[5000];
3463 Float_t zb[5000];
3464 Float_t txb[5000];
3465 Float_t tyb[5000];
3466 Int_t nfilmb[5000];
3467 //Float_t deltathetab[5000]; // Not used in this routine.
3468 //Float_t deltarb[5000]; // Not used in this routine.
3469// Int_t sizeb=N();
3470 Int_t sizeb=track->N();
3471 //Int_t nentries_withisizeb=0;
3472
3473 //=C= =====================================================================
3474 eHisto_nbtk_av ->Reset();
3475 eHisto_nbtk ->Reset();
3476 eHisto_longprofile_av ->Reset();
3477 eHisto_deltaR_mean ->Reset();
3478 eHisto_deltaT_mean ->Reset();
3479 eHisto_deltaR_rms ->Reset();
3480 eHisto_deltaT_rms ->Reset();
3481 eHisto_longprofile ->Reset();
3482 eHisto_deltaR ->Reset();
3483 eHisto_deltaT ->Reset();
3484 //=C= =====================================================================
3485 for (int id=0; id<57; id++) {
3486 theta[id]= 0;
3487 X0[id] = id*1300.*track->GetSegment(0)->TX() + track->GetSegment(0)->X();
3488 Y0[id] = id*1300.*track->GetSegment(0)->TY() + track->GetSegment(0)->Y();
3489 longprofile[id]=-1;
3490 // this is important, cause it means that where is -1 there is no BT reconstructed anymore
3491 // this is different for example to holese where N=0 , so one can distinguish!
3492 }
3493
3494 if (gEDBDEBUGLEVEL>3) if (gEDBDEBUGLEVEL >2) cout << "eHisto_longprofile->GetBinWidth() " << eHisto_longprofile->GetBinWidth(1) << endl;
3495 if (gEDBDEBUGLEVEL>3) if (gEDBDEBUGLEVEL >2) cout << "eHisto_longprofile->GetBinCenter() " << eHisto_longprofile->GetBinCenter(1) << endl;
3496 if (gEDBDEBUGLEVEL>3) if (gEDBDEBUGLEVEL >2) cout << "eHisto_longprofile->GetNbinsX() " << eHisto_longprofile->GetNbinsX() << endl;
3497
3498 TX0 = track->GetSegment(0)->TX();
3499 TY0 = track->GetSegment(0)->TY();
3500 for (Int_t j = 0; j<57; ++j) {
3501 Dr[j] = 0.03*j*1300. +20.0;
3502 //if (gEDBDEBUGLEVEL >2) cout << " DEBUG j= " << j << " Dr[j]= " << Dr[j] << endl;
3503 } // old relict from FJ. Do not remove. //
3504 // it represents somehow conesize......(??)
3505
3506 //=C= =====================================================================
3507 for (Int_t ibtke = 0; ibtke < track->N(); ibtke++) {
3508 xb[ibtke]=track->GetSegment(ibtke)->X();
3509 yb[ibtke]=track->GetSegment(ibtke)->Y();
3510 zb[ibtke]=track->GetSegment(ibtke)->Z();
3511 txb[ibtke]=track->GetSegment(ibtke)->TX();
3512 tyb[ibtke]=track->GetSegment(ibtke)->TY();
3513 // abs() of filmPID with respect to filmPID of first BT, plus 1: (nfilmb(0):=1 per definition):
3514 // Of course PID() have to be read correctly (by ReadEdbPVrec) correctly.
3515 // Fedra should do it.
3516 nfilmb[ibtke]=TMath::Abs(track->GetSegment(ibtke)->PID()-track->GetSegment(0)->PID())+1;
3517
3518 if (gEDBDEBUGLEVEL>2) {
3519 if (gEDBDEBUGLEVEL >2) cout << "ibtke= " <<ibtke << " xb[ibtke]= " << xb[ibtke] << " nfilmb[ibtke]= " << nfilmb[ibtke] << endl;
3520 }
3521 }
3522 //=C= =====================================================================
3523 //=C= loop over the basetracks in the shower (boucle sur les btk)
3524 for (Int_t ibtke = 0; ibtke < track->N(); ibtke++) {
3525 dist = sqrt((xb[ibtke]- X0[nfilmb[ibtke]-1])*(xb[ibtke]- X0[nfilmb[ibtke]-1])+(yb[ibtke]- Y0[nfilmb[ibtke]-1])*(yb[ibtke]- Y0[nfilmb[ibtke]-1]));
3526
3527 // inside the cone
3528 //if (gEDBDEBUGLEVEL >2) cout << "ibtke= " <<ibtke << " dist = " << dist << " Dr[nfilmb[ibtke]-1] = " << Dr[nfilmb[ibtke]-1] << endl;
3529 //if (gEDBDEBUGLEVEL >2) cout << " nfilmb[ibtke] = " << nfilmb[ibtke] << " nbfilm = " << nbfilm << endl;
3530
3531 // In old times there was here an additional condition which checked the BTs for the ConeDistance.
3532 // Since at this point in buildparametrisation_FJ the shower is -already fully reconstructed -
3533 // (by either EdbShoAlgo, or ShowRec program or manual list) this intruduces an additional cut
3534 // which is not correct because conetube size is algorithm specifiv (see Version.log.txt #20052010)
3535 // Thats wy we drop it here....
3536 // if (dist<Dr[nfilmb[ibtke]-1]&&nfilmb[ibtke]<=nbfilm) { // original if line comdition
3537 if (dist>Dr[nfilmb[ibtke]-1]) {
3538 if (gEDBDEBUGLEVEL>2) {
3539 if (gEDBDEBUGLEVEL >2) cout << " WARNING , In old times this cut (dist>Dr[nfilmb[ibtke]-1]) (had also to be fulfilled!"<<endl;
3540 if (gEDBDEBUGLEVEL >2) cout << " For this specific shower it seems not the case....." << endl;
3541 if (gEDBDEBUGLEVEL >2) cout << " You might want to check this shower again manualy to make sure everything is correct....." << endl;
3542 }
3543 }
3544 if (nfilmb[ibtke]<=nbfilm) {
3545 //if (gEDBDEBUGLEVEL >2) cout << "DEBUG CUTCONDITION WITHOUT THE (?_?_? WRONG ?_?_?) CONE DIST CONDITION....." << endl;
3546 // if (gEDBDEBUGLEVEL >2) cout << yes, this additional cut is not necessary anymore, see above....
3547
3548 eHisto_longprofile ->Fill(nfilmb[ibtke]);
3549 eHisto_longprofile_av ->Fill(nfilmb[ibtke]);
3550
3551 Double_t DR=0; //Extrapolate the old stlye way:
3552 Double_t Dx=xb[ibtke]-(xb[0]+(zb[ibtke]-zb[0])*txb[0]);
3553 Double_t Dy=yb[ibtke]-(yb[0]+(zb[ibtke]-zb[0])*tyb[0]);
3554 DR=TMath::Sqrt(Dx*Dx+Dy*Dy);
3555 eHisto_transprofile_av->Fill(DR);
3556 eHisto_transprofile->Fill(DR);
3557
3558 theta[nfilmb[ibtke]]+= (TX0-txb[ibtke])*(TX0-txb[ibtke])+(TY0-tyb[ibtke])*(TY0-tyb[ibtke]);
3559 if (ibtke>0&&nfilmb[ibtke]<=nbfilm) {
3560 // eHisto_deltaR ->Fill(deltarb[ibtke]);
3561 // eHisto_deltaT ->Fill(deltathetab[ibtke]);
3562 // uses shower_deltarb,shower_deltathetab just calculated in dummy routine above.
3563 // The first BT is NOT used for this filling since the extrapolated values of
3564 // shower_deltarb and shower_deltathetab are set manually to 0.5 and 200, due to
3565 // historical reasons (these variables com from the e/pi separation stuff).
3566 eHisto_deltaR ->Fill(shower_deltarb[ibtke]);
3567 eHisto_deltaT ->Fill(shower_deltathetab[ibtke]);
3568 }
3569 }
3570 }//==C== END OF loop over the basetracks in the shower
3571 //if (gEDBDEBUGLEVEL >2) cout <<"---------------------------------------"<<endl;
3572 //=======================================================================================
3573 //==C== Fill NumberBT Histogram for all showers:
3574 eHisto_nbtk ->Fill(sizeb);
3575 eHisto_nbtk_av ->Fill(sizeb);
3576 //==C== Fill dR,dT Mean and RMS Histos for all showers:
3577 eHisto_deltaR_mean ->Fill(eHisto_deltaR->GetMean());
3578 eHisto_deltaT_mean ->Fill(eHisto_deltaT->GetMean());
3579 eHisto_deltaR_rms ->Fill(eHisto_deltaR->GetRMS());
3580 eHisto_deltaT_rms ->Fill(eHisto_deltaT->GetRMS());
3581 //=======================================================================================
3582
3583 // Fill the longprofile array: (NEW VERSION (arrayindex 0 is same plate as Initiator BT))
3584 for (Int_t i=1; i<=nbfilm; ++i) {
3585 // longprofile[i-1] = eHisto_longprofile->GetBinContent(i); /// OLD VERSION for (Int_t i=1; i<=nbfilm+1; ++i)
3586 longprofile[i-1] = (Int_t)(eHisto_longprofile->GetBinContent(i+1)); // NEW VERSION (arrayindex 0 is same plate as Initiator BT)
3587 //test+=longprofile[i-1] ;
3588 if (gEDBDEBUGLEVEL>1) {
3589 if (gEDBDEBUGLEVEL >2) cout << "i= " << i << " longprofile[i-1] "<< longprofile[i-1] << " eHisto_longprofile->GetBinContent(i) " << eHisto_longprofile->GetBinContent(i+1)<< endl;
3590 }
3591 if (i==nbfilm) {
3592 if (gEDBDEBUGLEVEL >2) cout << "i==nbfilm:" << endl;
3593 longprofile[i-1] = (Int_t)(eHisto_longprofile->GetBinContent(i+1));
3594 }
3595 //
3596 }
3597 // Rather strange but I have to put it explicetly here, otherwise last bin isnt filled...
3598 longprofile[nbfilm] = (Int_t)(eHisto_longprofile->GetBinContent(nbfilm+1));
3599 for (Int_t i=nbfilm+2; i<57; ++i) {
3600 longprofile[i-1] = -1;
3601 }
3602 //----------------------------------------------------------------------------------------
3603
3604 // - Inclusion. only for libShower, EdbShowerRec class!
3605 Float_t eShowerAxisAngle=track->GetSegment(0)->Theta();
3606 // - Inclusion. only for libShower, EdbShowerRec class; END;
3607
3608
3609 // Now set parametrisation values:
3610 eParaShowerAxisAngle=eShowerAxisAngle;
3611 eParanseg=track->N();
3616 for (int ii=0; ii<57; ii++) eParalongprofile[ii]=longprofile[ii];
3617 //for (int ii=0; ii<57; ii++) cout << " ii= " << ii << " eParalongprofile[ii]= " << eParalongprofile[ii] << " longprofile[ii]= " << longprofile[ii] << endl;
3618
3619 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::GetPara Fill the Para structure with values from a track.... done." << endl;
3620 return;
3621}
brick X0
Definition: RecDispMC.C:112
Float_t shower_yb[5000]
Definition: ShowRec.h:374
Float_t shower_tyb[5000]
Definition: ShowRec.h:377
Float_t shower_deltathetab[5000]
Definition: ShowRec.h:379
Float_t shower_zb[5000]
Definition: ShowRec.h:375
Float_t shower_deltarb[5000]
Definition: ShowRec.h:378
Float_t shower_xb[5000]
Definition: ShowRec.h:373
Float_t shower_txb[5000]
Definition: ShowRec.h:376
Definition: EdbSegP.h:21
Float_t TX() const
tangens = deltaX/deltaZ
Definition: EdbSegP.h:175
Float_t X() const
Definition: EdbSegP.h:173
Float_t Z() const
Definition: EdbSegP.h:153
Float_t Y() const
Definition: EdbSegP.h:174
void PrintNice() const
Definition: EdbSegP.cxx:392
Float_t TY() const
tangens = deltaY/deltaZ
Definition: EdbSegP.h:176
int sizeb
Definition: check_shower.C:38
UInt_t id
Definition: tlg2couples.C:117

◆ GetRecoShowerArrayN()

Int_t EdbShowerAlgESimple::GetRecoShowerArrayN ( ) const
inline
455 {
456 return eRecoShowerArrayN;
457 }

◆ GetSpecifications()

void EdbShowerAlgESimple::GetSpecifications ( )
3101 {
3102 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::GetSpecifications" << endl;
3103 return;
3104}

◆ GetSpecType()

Int_t EdbShowerAlgESimple::GetSpecType ( Int_t  SpecificationType)
inline
527 {
528 return eSpecificationType[SpecificationType];
529 }
Int_t eSpecificationType[7]
Definition: EdbShowerAlg.h:350

◆ Help()

void EdbShowerAlgESimple::Help ( )
3880{
3881 cout << "-------------------------------------------------------------------------------------------"<< endl<< endl;
3882 cout << "EdbShowerAlgESimple::Help:"<<endl;
3883 cout << "EdbShowerAlgESimple::Help:"<<endl;
3884 cout << "EdbShowerAlgESimple::Help: The easiest way: "<<endl;
3885 cout << "EdbShowerAlgESimple::Help: EdbShowerAlgESimple* ShowerAlgE1 = new EdbShowerAlgESimple(); "<<endl;
3886 cout << "EdbShowerAlgESimple::Help: ShowerAlgE1->DoRun(RecoShowerArray); // RecoShowerArray an TObjArray of EdbTrackP* "<<endl;
3887 cout << "EdbShowerAlgESimple::Help: ShowerAlgE1->WriteNewRootFile(); "<<endl;
3888 cout << "EdbShowerAlgESimple::Help: "<<endl;
3889 cout << "EdbShowerAlgESimple::Help: In shower_E.C // E.C you find more hints...."<<endl;
3890 cout << "EdbShowerAlgESimple::Help: Updates will follow........"<<endl;
3891 cout << "-------------------------------------------------------------------------------------------"<< endl<< endl;
3892}

◆ Init()

void EdbShowerAlgESimple::Init ( void  )
protected

This Plate Binnning is the _final_binning and will not be refined furthermore!

2574{
2576 ePlateBinning[0]=10;
2577 ePlateBinning[1]=12;
2578 ePlateBinning[2]=14;
2579 ePlateBinning[3]=16;
2580 ePlateBinning[4]=18;
2581 ePlateBinning[5]=20;
2582 ePlateBinning[6]=23;
2583 ePlateBinning[7]=26;
2584 ePlateBinning[8]=29;
2585 ePlateBinning[9]=32;
2586 ePlateBinning[10]=35;
2587 ePlateBinning[11]=40;
2588 ePlateBinning[12]=45;
2589 ePlateBinning[13]=45;
2590 ePlateBinning[14]=45; // 12,13,14 are the same, cause I have seen that it makes no sense to produce these,
2591 // since E,and sigmaE for 45 plates doesnt change anymore at all...
2592 // and how often !can! we scan more than 45 plates for a shower?
2593
2594 ePlateNumber=20; // 20 plates by default...
2595 eSpecificationType[0] = 1; // cp files by default...
2596 eSpecificationType[1] = 1; // ele by default...
2597 eSpecificationType[2] = 0; // neuch eff by default...
2598 eSpecificationType[3] = 0; // 0..20 GeV by default...
2599 eSpecificationType[4] = 0; // next before weight by default...
2600 eSpecificationType[5] = 5; // 20 plates by default...
2601
2603 eSpecificationTypeString[1]="GAMMA";
2604 eSpecificationTypeString[2]="Neuch";
2605 eSpecificationTypeString[3]="TypeA";
2606 eSpecificationTypeString[4]="Next";
2608
2610
2611 for (int i=0; i<15; i++) {
2613 }
2614
2615 EffFunc_all = new TF1("all","1.0-0.00000001*x",0,0.95);
2616 EffFunc_edefault = new TF1("default","0.94-0.216*x-0.767*x*x+1.865*x*x*x",0,0.95);
2617 EffFunc_elletroni = new TF1("elletroni","0.79+0.38*x-7.63*x*x+25.13*x*x*x-24.6*x*x*x*x",0,0.95);
2618 EffFunc_neuchmicro = new TF1("neuchmicro","0.94-0.955*x+1.80*x*x-0.95*x*x*x",0,0.95);
2619 EffFunc_MiddleFix = new TF1("MiddleFix","0.5*(0.888361-1.299*x+1.49198*x*x+1.64668*x*x*x-2.63758*x*x*x*x+0.79+0.38*x-7.63*x*x+25.13*x*x*x-24.6*x*x*x*x)",0,0.95);
2620 EffFunc_LowEff = new TF1("LowEff","0.85*0.5*(0.888361-1.299*x+1.49198*x*x+1.64668*x*x*x-2.63758*x*x*x*x+0.79+0.38*x-7.63*x*x+25.13*x*x*x-24.6*x*x*x*x)",0,0.95);
2621 EffFunc_UserEfficiency = new TF1("EffFunc_UserEfficiency","0.94-0.955*x+1.80*x*x-0.95*x*x*x",0,0.95); //new TSpline3();
2622
2623
2624 // Standard supposed efficiency of the showers/tracks that are to be evaluated:
2625 // The correct EffFunc_* is then choosen by taking the closes eff func to that one.
2626 // Initial parameters are the efficiency parameters for the "Neuch" parametrisation
2627 Double_t xarr[7]= {0,0.1,0.2,0.3,0.4,0.5,0.6};
2628 Double_t yarr[7]= {0.95,0.9,0.80,0.75,0.6,0.5,0.4};
2629 //if (NULL!=eEfficiencyParametrisation) delete eEfficiencyParametrisation;
2630 eEfficiencyParametrisation = new TSpline3("",xarr,yarr,5,0,0,0.65);
2631
2632 eHisto_nbtk_av = new TH1D("eHisto_nbtk_av","Average basetrack numbers",21,0.0,10.0);
2633 eHisto_longprofile_av = new TH1D("eHisto_longprofile_av","Basetracks per emulsion number",57,0.0,57.0);
2634 eHisto_transprofile_av = new TH1D("eHisto_transprofile_av","Basetracks in trans distance",8,0.0,800.0);
2635 eHisto_deltaR_mean = new TH1D("eHisto_deltaR_mean","Mean #deltar of all BTs in one shower",100,0.0,100.0);
2636 eHisto_deltaT_mean = new TH1D("eHisto_deltaT_mean","Mean #delta#theta of all BTs in one shower",100,0.0,0.1);
2637 eHisto_deltaR_rms = new TH1D("eHisto_deltaR_rms","RMS #deltar of all BTs in one shower",100,0.0,100.0);
2638 eHisto_deltaT_rms = new TH1D("eHisto_deltaT_rms","RMS #delta#theta of all BTs in one shower",100,0.0,0.1);
2639 eHisto_nbtk = new TH1D("eHisto_nbtk","Basetracks in the shower",50,0.0,100.0);
2640 eHisto_longprofile = new TH1D("eHisto_longprofile","Basetracks per emulsion number",57,0.0,57.0);
2641 eHisto_deltaR = new TH1D("eHisto_deltaR","Single #deltar of all BTs in Shower",100,0.0,150.0);
2642 eHisto_deltaT = new TH1D("eHisto_deltaT","Single #delta#theta of all BTs in Shower",150,0.0,0.15);
2643 eHisto_transprofile = new TH1D("eHisto_transprofile","Basetracks in trans distance",8,0.0,800.0);
2644
2645 // CreateANN
2646 CreateANN();
2647
2648 //First Time, we also Update()
2649 Update();
2650
2651 // Set Strings:
2652 InitStrings();
2653
2654 // Create EnergyArray
2655 eEnergyArray=new TArrayF(99999); // no way to adapt tarrayF on the fly...
2656 eEnergyArrayUnCorrected=new TArrayF(99999);
2657 eEnergyArraySigmaCorrected=new TArrayF(99999);
2658
2659
2660 // Read Energy Resolution Lookup tables:
2661 ReadTables();
2662
2663 cout << "EdbShowerAlgESimple::Init()...done."<< endl;
2664 return;
2665}
void CreateANN()
Definition: EdbShowerAlg.cxx:2691
void ReadTables()
Definition: EdbShowerAlg.cxx:3915
void InitStrings()
Definition: EdbShowerAlg.cxx:2670
Int_t ANN_nPlates_ARRAY[15]
Definition: EdbShowerAlg.h:366
TString eSpecificationTypeString[7]
Definition: EdbShowerAlg.h:351
void Update()
Definition: EdbShowerAlg.cxx:3189

◆ InitStrings()

void EdbShowerAlgESimple::InitStrings ( )
2671{
2674 eSpecificationTypeStringArray[1][0]="GAMMA";
2675 eSpecificationTypeStringArray[1][1]="ELECTRON";
2676 eSpecificationTypeStringArray[2][0]="Neuch";
2678 eSpecificationTypeStringArray[2][2]="MiddleFix";
2680 eSpecificationTypeStringArray[3][0]="TypeA";
2681 eSpecificationTypeStringArray[3][1]="TypeABCD";
2682 eSpecificationTypeStringArray[4][0]="Next";
2683 eSpecificationTypeStringArray[4][1]="Before";
2685 cout << "EdbShowerAlgESimple::InitStrings()...done."<< endl;
2686 return;
2687}
TString eSpecificationTypeStringArray[7][7]
Definition: EdbShowerAlg.h:352

◆ LoadSpecificationWeightFile()

void EdbShowerAlgESimple::LoadSpecificationWeightFile ( )
3134 {
3135 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::LoadSpecificationWeightFile from the following specifications:" << endl;
3137 TString weigth="NULL";
3138 return;
3139}
void PrintSpecifications()
Definition: EdbShowerAlg.cxx:3108

◆ Print()

void EdbShowerAlgESimple::Print ( )
3652 {
3653 cout << "-------------------------------------------------------------------------------------------"<< endl;
3654 cout << "EdbShowerAlgESimple::Print " << endl;
3656
3657 cout << "EdbShowerAlgESimple::Print Now correction factors and weightfilestrings:" << endl;
3658 for (int i=0; i<15; i++) {
3659 cout << "EdbShowerAlgESimple::Print eANN_MLP_CORR_0[" << i << "]="<<eANN_MLP_CORR_0[i] << endl;
3660 cout << "EdbShowerAlgESimple::Print eANN_MLP_CORR_1[" << i << "]="<<eANN_MLP_CORR_1[i] << endl;
3661 cout << "EdbShowerAlgESimple::Print ANN_WeightFile_ARRAY[" << i << "]="<<ANN_WeightFile_ARRAY[i] << endl;
3662 }
3663
3664 cout << "EdbShowerAlgESimple::Print ...done." << endl;
3665 cout << "-------------------------------------------------------------------------------------------"<< endl << endl;
3666 return;
3667}

◆ PrintEfficiencyParametrisation()

void EdbShowerAlgESimple::PrintEfficiencyParametrisation ( )
3899{
3900 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation for angles theta=0,0.1,..,0.6:" << endl;
3901 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation->Eval(0): "<< eEfficiencyParametrisation->Eval(0)<< endl;
3902 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation->Eval(0.1): "<< eEfficiencyParametrisation->Eval(0.1)<< endl;
3903 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation->Eval(0.2): "<< eEfficiencyParametrisation->Eval(0.2)<< endl;
3904 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation->Eval(0.3): "<< eEfficiencyParametrisation->Eval(0.3)<< endl;
3905 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation->Eval(0.4): "<< eEfficiencyParametrisation->Eval(0.4)<< endl;
3906 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation->Eval(0.5): "<< eEfficiencyParametrisation->Eval(0.5)<< endl;
3907 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation->Eval(0.6): "<< eEfficiencyParametrisation->Eval(0.6)<< endl;
3908 return;
3909}

◆ PrintSpecifications()

void EdbShowerAlgESimple::PrintSpecifications ( )
3108 {
3109 cout << "EdbShowerAlgESimple::PrintSpecifications" << endl;
3110 cout << "EdbShowerAlgESimple:: eSpecificationType[0] (Dataset: (linked tracks/full cp) LT/CP) = " << GetSpecType(0) << endl;
3111 cout << "EdbShowerAlgESimple:: eSpecificationType[1] (ShowerID: gamma/electron/pion weight) = " << GetSpecType(1) << endl;
3112 cout << "EdbShowerAlgESimple:: eSpecificationType[2] (ScanEff: Neuch/All/MiddleFix/LowEff) = " << GetSpecType(2) << endl;
3113 cout << "EdbShowerAlgESimple:: eSpecificationType[3] (E range: 0..20/0..40) = " << GetSpecType(3) << endl;
3114 cout << "EdbShowerAlgESimple:: eSpecificationType[4] (Npl weight: next before/after) = " <<GetSpecType(4) << endl;
3115 cout << "EdbShowerAlgESimple:: eSpecificationType[5] (Npl weight: 10,12,...,45) = " << GetSpecType(5) << endl;
3116 cout << "EdbShowerAlgESimple:: " << endl;
3117 cout << "EdbShowerAlgESimple:: eSpecificationTypeString[0] (take CP or linked_tracks) = " << eSpecificationTypeString[0].Data() << endl;
3118 cout << "EdbShowerAlgESimple:: eSpecificationTypeString[1] (gamma/electron/pion weight) = " << eSpecificationTypeString[1].Data() << endl;
3119 cout << "EdbShowerAlgESimple:: eSpecificationTypeString[2] (ScanEff: Neuch/All/MiddleFix/LowEff) = " << eSpecificationTypeString[2].Data() << endl;
3120 cout << "EdbShowerAlgESimple:: eSpecificationTypeString[3] (E range: 0..20/0..40) = " << eSpecificationTypeString[3].Data() << endl;
3121 cout << "EdbShowerAlgESimple:: eSpecificationTypeString[4] (Npl weight: next before/after) = " << eSpecificationTypeString[4].Data() << endl;
3122 cout << "EdbShowerAlgESimple:: eSpecificationTypeString[5] (Npl weight: 10,12,...,45) = " << eSpecificationTypeString[5].Data() << endl;
3123 cout << "EdbShowerAlgESimple:: " << endl;
3124 cout << "EdbShowerAlgESimple:: In case you want to change a specification then do for example:" << endl;
3125 cout << "EdbShowerAlgESimple:: EdbShowerAlgESimple->SetSpecification(2,3) " << endl;
3126 cout << "EdbShowerAlgESimple:: It will change the ScanEff from the actual value to _Low_." << endl;
3127 cout << "EdbShowerAlgESimple::" << endl;
3128 cout << "EdbShowerAlgESimple::PrintSpecifications...done." << endl;
3129 return;
3130}

◆ ReadCorrectionFactors()

void EdbShowerAlgESimple::ReadCorrectionFactors ( TString  weigthstring,
Float_t &  p0,
Float_t &  p1 
)

Read Linear Correction factors p0 p1 from this file.
Format: p0 p1

3625 {
3628 float pp0,pp1;
3629 const char* name=weigthstring.Data();
3630 FILE * pFile=NULL;
3631 if (gEDBDEBUGLEVEL >2) cout << " open file " << endl;
3632 pFile = fopen (name,"r");
3633 if (NULL==pFile) return;
3634 int dummy = fscanf (pFile, "%f %f", &pp0, &pp1);
3635 if (dummy!=2) {
3636 cout << "EdbShowerAlgESimple::ReadCorrectionFactors ERROR! Wrong formatted file! Could not read the two parameters! Set them to default values! Please check if the weightfiles in...: "<< endl;
3637 cout << weigthstring.Data() << endl;
3638 cout << "EdbShowerAlgESimple::ReadCorrectionFactors ... exist!"<< endl;
3639 p0=0.0;
3640 p1=1.0;
3641 }
3642 p0=pp0;
3643 p1=pp1;
3644 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::ReadCorrectionFactors p0: " << p0 << " p1 " << p1 << endl;
3645 if (gEDBDEBUGLEVEL >2) printf ("I have read: %f and %f , in total %d arguments.\n",pp0,pp1,dummy);
3646 fclose (pFile);
3647 return;
3648}
Int_t dummy
Definition: merge_Energy_SytematicSources_Electron.C:7
fclose(pFile)
pFile
Definition: merge_Energy_SytematicSources_Electron.C:25
const char * name
Definition: merge_Energy_SytematicSources_Electron.C:24
#define NULL
Definition: nidaqmx.h:84

◆ ReadTables()

void EdbShowerAlgESimple::ReadTables ( )
3916{
3917 cout << "EdbShowerAlgESimple::ReadTables"<<endl;
3918
3920
3921 cout << "EdbShowerAlgESimple::ReadTables...done."<<endl;
3922 return;
3923}
void ReadTables_Energy()
Definition: EdbShowerAlg.cxx:3927

◆ ReadTables_Energy()

void EdbShowerAlgESimple::ReadTables_Energy ( )
3928{
3929 cout << "EdbShowerAlgESimple::ReadTables_Energy"<<endl;
3930
3931 char name [100];
3932 FILE * fFile;
3933
3934 Int_t N_E=12;
3935 Int_t N_NPL=13;
3936 Int_t NPL[13]= {10,12,14,16,18,20,23,26,29,32,35,40,45};
3937 Double_t E[12] = {0.5,1.0,2.0, 4.0, 0.75, 1.5, 3.0, 6.0, 8.0, 16.0, 32.0 ,64.0};
3938 Int_t E_ASCEND[12] = {0,4,1,5,2,6,3,7,8,9,10,11};
3939 Int_t npl;
3940 Float_t energyresolution;
3941 Double_t E_Array_Ascending[12]= {0.5,0.75,1.0,1.5,2.0,3.0,4.0,6.0, 8.0, 16.0, 32.0 ,64.0};
3942 Double_t E_Resolution[12]= {0.5,0.75,1.0,1.5,2.0,3.0,4.0,6.0, 8.0, 16.0, 32.0 ,64.0}; // initialize with 100% resolution...
3943
3944 // First Create the ObjArray storing the Splines...
3945 eSplineArray_Energy_Stat_Electron = new TObjArray();
3946 eSplineArray_Energy_Stat_Gamma = new TObjArray();
3947 eSplineArray_Energy_Sys_Electron = new TObjArray();
3948 eSplineArray_Energy_Sys_Gamma = new TObjArray();
3949
3950 //We need to give the full path to look for the Statistics tables:
3951 TString basicstring = TString(gSystem->ExpandPathName("$FEDRA_ROOT"));
3952 TString addstring = ("/src/libShower/");
3953 TString tablestring;
3954
3955 // A) Table: Statistics: Electrons: "libShower_Energy_Statistics_Electron.txt"
3956 //sprintf(name,"tables/libShower_Energy_Statistics_Electron.txt");
3957 tablestring=basicstring+addstring+TString("tables/libShower_Energy_Statistics_Electron.txt");
3958 fFile = fopen (tablestring.Data(),"r");
3959 cout << "EdbShowerAlgESimple::ReadTables_Energy() name = " << tablestring.Data() << endl;
3960
3961 for (int i=0; i<N_NPL; i++) {
3962 // Read energy values into array:
3963 int narg=fscanf(fFile, "%d ",&npl);
3964 if (gEDBDEBUGLEVEL>2) cout << "NPL= " << npl << endl;
3965 for (int j=0; j<N_E-1; j++) {
3966 int narg=fscanf(fFile, "%f",&energyresolution);
3967 E_Resolution[j]=(Double_t)energyresolution;
3968 if (gEDBDEBUGLEVEL>2) cout << " energyresolution @ " << E_Array_Ascending[j] << " = " << energyresolution << endl;
3969 }
3970 narg=fscanf(fFile, "%f \n",&energyresolution);
3971 E_Resolution[N_E-1]=(Double_t)energyresolution;
3972 if (gEDBDEBUGLEVEL>2) cout << " energyresolution @ " << E_Array_Ascending[N_E-1] << " = " << energyresolution << endl;
3973
3974 TString splinename=TString(Form("Spline_Stat_Electron_Npl_%d",NPL[i]));
3975 TSpline3* Spline = new TSpline3(splinename,E_Array_Ascending,E_Resolution,10);
3976
3977 if (gEDBDEBUGLEVEL>2) cout << "Created Spline. Add Spline to ObjArray." << endl;
3979 }
3980 fclose (fFile);
3981
3982
3983 // B) Table: Statistics: Gamma: "libShower_Energy_Statistics_Gamma.txt"
3984 tablestring=basicstring+addstring+TString("tables/libShower_Energy_Statistics_Gamma.txt");
3985 fFile = fopen (tablestring.Data(),"r");
3986 cout << "EdbShowerAlgESimple::ReadTables_Energy() name = " << tablestring.Data() << endl;
3987
3988 for (int i=0; i<N_NPL; i++) {
3989 // Read energy values into array:
3990 int narg=fscanf(fFile, "%d ",&npl);
3991 if (gEDBDEBUGLEVEL>2) cout << "NPL= " << npl << endl;
3992 for (int j=0; j<N_E-1; j++) {
3993 int narg=fscanf(fFile, "%f",&energyresolution);
3994 E_Resolution[j]=(Double_t)energyresolution;
3995 if (gEDBDEBUGLEVEL>2) cout << " energyresolution @ " << E_Array_Ascending[j] << " = " << energyresolution << endl;
3996 }
3997 narg=fscanf(fFile, "%f \n",&energyresolution);
3998 E_Resolution[N_E-1]=(Double_t)energyresolution;
3999 if (gEDBDEBUGLEVEL>2) cout << " energyresolution @ " << E_Array_Ascending[N_E-1] << " = " << energyresolution << endl;
4000
4001 TString splinename=TString(Form("Spline_Stat_Electron_Npl_%d",NPL[i]));
4002 TSpline3* Spline = new TSpline3(splinename,E_Array_Ascending,E_Resolution,10);
4003
4004 if (gEDBDEBUGLEVEL>2) cout << "Created Spline. Add Spline to ObjArray." << endl;
4005 eSplineArray_Energy_Stat_Gamma->Add(Spline);
4006 }
4007 fclose (fFile);
4008
4009
4010
4011 // C) Table: Systematics: Electrons: "libShower_Energy_Systematics_Electron.txt"
4012 tablestring=basicstring+addstring+TString("tables/libShower_Energy_Systematics_Electron.txt");
4013 fFile = fopen (tablestring.Data(),"r");
4014 cout << "EdbShowerAlgESimple::ReadTables_Energy() name = " << tablestring.Data() << endl;
4015
4016 for (int i=0; i<N_NPL; i++) {
4017 // Read energy values into array:
4018 int narg=fscanf(fFile, "%d ",&npl);
4019 if (gEDBDEBUGLEVEL>2) cout << "NPL= " << npl << endl;
4020 for (int j=0; j<N_E-1; j++) {
4021 int narg=fscanf(fFile, "%f",&energyresolution);
4022 E_Resolution[j]=(Double_t)energyresolution;
4023 if (gEDBDEBUGLEVEL>2) cout << " energyresolution @ " << E_Array_Ascending[j] << " = " << energyresolution << endl;
4024 }
4025 narg=fscanf(fFile, "%f \n",&energyresolution);
4026 E_Resolution[N_E-1]=(Double_t)energyresolution;
4027 if (gEDBDEBUGLEVEL>2) cout << " energyresolution @ " << E_Array_Ascending[N_E-1] << " = " << energyresolution << endl;
4028
4029 TString splinename=TString(Form("Spline_Sysy_Electron_Npl_%d",NPL[i]));
4030 TSpline3* Spline = new TSpline3(splinename,E_Array_Ascending,E_Resolution,10);
4031
4032 if (gEDBDEBUGLEVEL>2) cout << "Created Spline. Add Spline to ObjArray." << endl;
4034 }
4035 fclose (fFile);
4036
4037
4038 // D) Table: Systematics: Gamma: "libShower_Energy_Systematics_Gamma.txt"
4039 tablestring=basicstring+addstring+TString("tables/libShower_Energy_Systematics_Gamma.txt");
4040 fFile = fopen (tablestring.Data(),"r");
4041 cout << "EdbShowerAlgESimple::ReadTables_Energy() name = " << tablestring.Data() << endl;
4042
4043 for (int i=0; i<N_NPL; i++) {
4044 // Read energy values into array:
4045 int narg=fscanf(fFile, "%d ",&npl);
4046 if (gEDBDEBUGLEVEL>2) cout << "NPL= " << npl << endl;
4047 for (int j=0; j<N_E-1; j++) {
4048 int narg=fscanf(fFile, "%f",&energyresolution);
4049 E_Resolution[j]=(Double_t)energyresolution;
4050 if (gEDBDEBUGLEVEL>2) cout << " energyresolution @ " << E_Array_Ascending[j] << " = " << energyresolution << endl;
4051 }
4052 narg=fscanf(fFile, "%f \n",&energyresolution);
4053 E_Resolution[N_E-1]=(Double_t)energyresolution;
4054 if (gEDBDEBUGLEVEL>2) cout << " energyresolution @ " << E_Array_Ascending[N_E-1] << " = " << energyresolution << endl;
4055
4056 TString splinename=TString(Form("Spline_Sysy_Electron_Npl_%d",NPL[i]));
4057 TSpline3* Spline = new TSpline3(splinename,E_Array_Ascending,E_Resolution,10);
4058
4059 if (gEDBDEBUGLEVEL>2) cout << "Created Spline. Add Spline to ObjArray." << endl;
4060 eSplineArray_Energy_Sys_Gamma->Add(Spline);
4061 }
4062 fclose (fFile);
4063
4064
4065 cout << "EdbShowerAlgESimple::ReadTables_Energy...done."<<endl;
4066 return;
4067}
TObjArray * eSplineArray_Energy_Sys_Electron
Definition: EdbShowerAlg.h:428
TObjArray * eSplineArray_Energy_Stat_Electron
Definition: EdbShowerAlg.h:426

◆ Set0()

void EdbShowerAlgESimple::Set0 ( )
protected

Reset Values

2530{
2537 ePlateNumber=20;
2540
2541 eWeightFileString="NULL";
2542 for (int i=0; i<7; i++) {
2543 eSpecificationType[i]=-1;
2545 for (int j=0; j<7; j++) {
2547 }
2548 }
2549 for (int i=0; i<15; i++) {
2551 eANN_MLP_CORR_0[i]=0.0;
2552 eANN_MLP_CORR_1[i]=1.0;
2553 }
2555
2556 // Reset "eEfficiencyParametrisation" to the "Neuch" eff.
2557 cout << "EdbShowerAlgESimple::Set0() Reset to the Neuch eff." << endl;
2558 cout << "EdbShowerAlgESimple::Set0() eEfficiencyParametrisation= " << eEfficiencyParametrisation << endl;
2559 //if (NULL!=eEfficiencyParametrisation) delete eEfficiencyParametrisation;
2560 Double_t xarr[6]= {0,0.1,0.2,0.3,0.4,0.5};
2561 Double_t yarr[6]= {0.95,0.9,0.80,0.75,0.6,0.5};
2562 eEfficiencyParametrisation = new TSpline3("",xarr,yarr,5,0,0,0.6);
2563 cout << "EdbShowerAlgESimple::Set0() Reset to the Neuch eff. ...done." << endl;
2564
2565
2566 cout << "EdbShowerAlgESimple::Set0()...done."<< endl;
2567 return;
2568}
Bool_t eSpecificationIsChanged
Definition: EdbShowerAlg.h:353

◆ SetCalibrationOffset()

void EdbShowerAlgESimple::SetCalibrationOffset ( Float_t  CalibrationOffset)
inline
514 {
515 eCalibrationOffset=CalibrationOffset;
516 }

◆ SetCalibrationSlope()

void EdbShowerAlgESimple::SetCalibrationSlope ( Float_t  CalibrationSlope)
inline
517 {
518 eCalibrationSlope=CalibrationSlope;
519 }

◆ SetEfficiencyParametrisation()

void EdbShowerAlgESimple::SetEfficiencyParametrisation ( TSpline3 *  EfficiencyParametrisation)
inline
539 {
540 eEfficiencyParametrisation=EfficiencyParametrisation;
541 EfficiencyParametrisation->Print();
542 }

◆ SetEfficiencyParametrisationAngles()

void EdbShowerAlgESimple::SetEfficiencyParametrisationAngles ( )
4074{
4075 //EffFunc_UserEfficiency->Set
4076 return;
4077}

◆ SetEfficiencyParametrisationValues()

void EdbShowerAlgESimple::SetEfficiencyParametrisationValues ( Double_t *  Angles,
Double_t *  EffValuesAtAngles 
)

Set User Efficiency Angles by hand!
This is interface routine which than can be used by eda...

4082{
4083 cout << "EdbShowerAlgESimple::SetEfficiencyParametrisationValues()" << endl;
4084 cout << "ATTENTION: Array has to consist of efficiencies at !seven! angles: 0, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6..." << endl;
4085 cout << "ATTENTION: If not, then it might crash! " << endl;
4088 //delete eEfficiencyParametrisation;
4089 eEfficiencyParametrisation = new TSpline3("",Angles,EffValuesAtAngles,5,0,0,0.6);
4090
4091 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation for angles theta=0,0.1,..,0.6:" << endl;
4092 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation->Eval(0): "<< eEfficiencyParametrisation->Eval(0)<< endl;
4093 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation->Eval(0.1): "<< eEfficiencyParametrisation->Eval(0.1)<< endl;
4094 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation->Eval(0.2): "<< eEfficiencyParametrisation->Eval(0.2)<< endl;
4095 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation->Eval(0.3): "<< eEfficiencyParametrisation->Eval(0.3)<< endl;
4096 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation->Eval(0.4): "<< eEfficiencyParametrisation->Eval(0.4)<< endl;
4097 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation->Eval(0.5): "<< eEfficiencyParametrisation->Eval(0.5)<< endl;
4098 cout << "EdbShowerAlgESimple::eEfficiencyParametrisation->Eval(0.6): "<< eEfficiencyParametrisation->Eval(0.6)<< endl;
4099
4100 cout << "EdbShowerAlgESimple::SetEfficiencyParametrisationValues() Do UnSetForceSpecificationReload" << endl;
4102 return;
4103}
void UnSetForceSpecificationReload()
Definition: EdbShowerAlg.h:533

◆ SetForceSpecificationReload()

void EdbShowerAlgESimple::SetForceSpecificationReload ( )
inline
530 {
532 }

◆ SetPlateNumber()

void EdbShowerAlgESimple::SetPlateNumber ( Int_t  PlateNumber)
inline
524 {
525 ePlateNumber=PlateNumber;
526 }

◆ SetPlateNumberType()

void EdbShowerAlgESimple::SetPlateNumberType ( Int_t  PlateNumberType)
inline
521 {
522 ePlateNumberType=PlateNumberType;
523 }

◆ SetRecoShowerArray()

void EdbShowerAlgESimple::SetRecoShowerArray ( TObjArray *  RecoShowerArray)
inline
448 {
450 eRecoShowerArrayN=eRecoShowerArray->GetEntriesFast();
451 }

◆ SetRecoShowerArrayN()

void EdbShowerAlgESimple::SetRecoShowerArrayN ( Int_t  RecoShowerArrayN)
inline
452 {
453 eRecoShowerArrayN=RecoShowerArrayN;
454 }

◆ SetSpecifications()

void EdbShowerAlgESimple::SetSpecifications ( Int_t  sp0,
Int_t  sp1,
Int_t  sp2,
Int_t  sp3,
Int_t  sp4,
Int_t  sp5 
)
3143 {
3144 eSpecificationType[0]=sp0;
3145 eSpecificationType[1]=sp1;
3146 eSpecificationType[2]=sp2;
3147 eSpecificationType[3]=sp3;
3148 eSpecificationType[4]=sp4;
3149 eSpecificationType[5]=sp5;
3153 Update();
3154 return;
3155}

◆ SetSpecificationType()

void EdbShowerAlgESimple::SetSpecificationType ( Int_t  SpecificationType,
Int_t  SpecificationTypeVal 
)
3161{
3162// if (GetSpecType(SpecificationType)==SpecificationTypeVal) {
3163// cout << "dbShowerAlgESimple::SetSpecificationType() The given value is the same value as before. So just set the eForceSpecificationReload to False."<<endl;
3164// eForceSpecificationReload=kTRUE;
3165// return;
3166// }
3167 if (gEDBDEBUGLEVEL >1) cout << "EdbShowerAlgESimple:: Change Specification (" << SpecificationType << ") from " << GetSpecType(SpecificationType) << " -> " << SpecificationTypeVal << " . Reprint the changed Specification String: " << endl;
3168
3169 eSpecificationType[SpecificationType]=SpecificationTypeVal;
3170 eSpecificationTypeString[SpecificationType] = eSpecificationTypeStringArray[SpecificationType][SpecificationTypeVal];
3171
3172 if (gEDBDEBUGLEVEL >1) cout << eSpecificationTypeString[SpecificationType].Data() << endl;
3173
3176
3177 cout << "EdbShowerAlgESimple::SetSpecificationType eSpecificationIsChanged: " << eSpecificationIsChanged << endl;
3178 cout << "EdbShowerAlgESimple::SetSpecificationType eForceSpecificationReload: " << eForceSpecificationReload << endl;
3179
3180 cout << "Reprint Specifications: " << endl;
3182
3183 Update();
3184 return;
3185}

◆ SetWeightFileString()

void EdbShowerAlgESimple::SetWeightFileString ( TString  weightstring)
2985{
2986 eWeightFileString=weightstring;
2987 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::SetWeightFileString() eWeightFileString = " << eWeightFileString << endl;
2988 return;
2989}

◆ TrainNeuralNetwork()

void EdbShowerAlgESimple::TrainNeuralNetwork ( TString  weight,
Int_t  ANNType = 0 
)
inline
502 {
503 if (ANNType>=15) ANNType=14;
504 ANN_MLP_ARRAY[ANNType]->Train(100);
505 return;
506 }

◆ UnSetForceSpecificationReload()

void EdbShowerAlgESimple::UnSetForceSpecificationReload ( )
inline
533 {
535 }

◆ Update()

void EdbShowerAlgESimple::Update ( )
3189 {
3190 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::Update Does the following things in the order:" << endl;
3191 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::Update * According to the switch: set the right ANN of the Array as generic one." << endl;
3192 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::Update * According to the switch: load the right weightfile as generic one." << endl;
3193 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::Update * According to the switch: set the right correction parameters...." << endl;
3194
3195 TString basicstring = TString(gSystem->ExpandPathName("$FEDRA_ROOT"));
3196 TString addstring = ("/src/libShower/weights/Energy/");
3197
3198 // (*) SpecType 0 :
3199 if (eSpecificationType[0]==1) {
3200 addstring+="volumeSpec_CP/";
3202 }
3203 else {
3204 addstring+="volumeSpec_LT/";
3206 cout << "EdbShowerAlgESimple::Update WARNING::eSpecificationTypeString[0]=LT NOT YET SUPPORTED!" << endl;
3207 }
3208
3209 // (*) SpecType 1 :
3210 if (eSpecificationType[1]==2) {
3211 addstring+="particleSpec_pion/";
3212 eSpecificationTypeString[1]="pion";
3213 cout << "EdbShowerAlgESimple::Update WARNING::eSpecificationTypeString[1]=pion NOT YET SUPPORTED!" << endl;
3214 }
3215 else if (eSpecificationType[1]==1) {
3216 addstring+="particleSpec_electron/";
3217 eSpecificationTypeString[1]="electron";
3218 }
3219 else {
3220 addstring+="particleSpec_gamma/";
3221 eSpecificationTypeString[1]="gamma";
3222 }
3223
3224 // (*) SpecType 2 :
3225 if (eSpecificationType[2]==3) {
3226 addstring+="efficiencySpec_LowEff/";
3227 eSpecificationTypeString[2]="LowEff";
3228 }
3229 else if (eSpecificationType[2]==2) {
3230 addstring+="efficiencySpec_MiddleFix/";
3231 eSpecificationTypeString[2]="MiddleFix";
3232 }
3233 else if (eSpecificationType[2]==1) {
3234 addstring+="efficiencySpec_All/";
3235 eSpecificationTypeString[2]="All";
3236 }
3237 else {
3238 addstring+="efficiencySpec_Neuch/";
3239 eSpecificationTypeString[2]="Neuch";
3240 }
3241
3242 // (*) SpecType 3 :
3243 if (eSpecificationType[3]==1) {
3244 addstring+="trainrangeSpec_extended/";
3245 eSpecificationTypeString[3]="extended";
3246 cout << "EdbShowerAlgESimple::Update WARNING::eSpecificationTypeString[3]=extended NOT YET SUPPORTED!" << endl;
3247 }
3248 else {
3249 addstring+="trainrangeSpec_normal/";
3250 eSpecificationTypeString[3]="normal";
3251 }
3252
3253 if (gEDBDEBUGLEVEL >2) cout << "EdbShowAlgE_Simple::Update " << addstring.Data() << endl;
3254
3255 if (NULL==ANN_MLP_ARRAY[0]) {
3256 CreateANN();
3257 }
3258
3259 if (gEDBDEBUGLEVEL >2) cout << "EdbShowAlgE_Simple::Update Now load the different ANN weightstrings:" << endl;
3260 // This was the loading part for the ANN and
3261 // this is the loading part for the Correction Factors, which we take from the
3262 // generic file EnergyCorrections_Npl_%d.txt
3263 for (int ll=0; ll<15; ll++) {
3264 TString weigthstring=basicstring+addstring+TString(Form("weights_Npl_%d.txt",ANN_nPlates_ARRAY[ll]));
3265 if (gEDBDEBUGLEVEL >2) cout << "weigthstring = " << weigthstring.Data() << endl;
3266 ANN_MLP_ARRAY[ll]->LoadWeights(weigthstring);
3267 ANN_WeightFile_ARRAY[ll]=weigthstring;
3268
3269 TString correctionsfactorstring=basicstring+addstring+TString(Form("EnergyCorrections_Npl_%d.txt",ANN_nPlates_ARRAY[ll]));
3270 if (gEDBDEBUGLEVEL >2) cout << "correctionsfactorstring = " << correctionsfactorstring.Data() << endl;
3271 Float_t p0,p1;
3272 p0=0.0;
3273 p1=1.0;
3274 if (ll<10) ReadCorrectionFactors(correctionsfactorstring,p0,p1);
3275// Correction files for ll>11 dont exist yet..
3276 eANN_MLP_CORR_0[ll]=p0;
3277 eANN_MLP_CORR_1[ll]=p1;
3278
3279 }
3280 if (gEDBDEBUGLEVEL >2) cout << "EdbShowerAlgESimple::Update * According to the switch: set the right ANN of the Array as generic one." << endl;
3281
3282 // (*) SpecType 4:
3283 // Get NplIndexNr:
3284 int check_Npl_index=0;
3285 GetNplIndexNr(ePlateNumber,check_Npl_index, eSpecificationType[4]);
3286
3287 if (check_Npl_index==0) eSpecificationTypeString[4]="Next";
3288 if (check_Npl_index==1) eSpecificationTypeString[4]="Before";
3289 eSpecificationTypeString[5]=TString(Form("%d",ePlateNumber));
3290
3291 // Set Generic ANN_Layout
3292 ANN_MLP=ANN_MLP_ARRAY[check_Npl_index];
3293
3294 if (gEDBDEBUGLEVEL >2) cout << "EdbShowAlgE_Simple::Update WARNING:: Weightfiles for _LT_ not produced yet!...." << endl;
3295 return;
3296}
void ReadCorrectionFactors(TString weigthstring, Float_t &p0, Float_t &p1)
Definition: EdbShowerAlg.cxx:3625

◆ WriteNewRootFile() [1/3]

void EdbShowerAlgESimple::WriteNewRootFile ( )
inline
576 {
577 WriteNewRootFile("Shower.root","treebranch");
578 }
void WriteNewRootFile()
Definition: EdbShowerAlg.h:576

◆ WriteNewRootFile() [2/3]

void EdbShowerAlgESimple::WriteNewRootFile ( TString  sourcefilename)
inline
579 {
580 WriteNewRootFile(sourcefilename,"treebranch");
581 }

◆ WriteNewRootFile() [3/3]

void EdbShowerAlgESimple::WriteNewRootFile ( TString  sourcefilename,
TString  treename 
)

eShowerTree->SetBranchAddress("output", &shower_output);

3673{
3674 cout << "-------------------------------------------------------------------------------------------"<< endl;
3675 cout << endl;
3676 cout << "EdbShowerAlgESimple::WriteNewRootFile -------------------------------------------"<< endl;
3677 cout << "EdbShowerAlgESimple::WriteNewRootFile This is NOT a good solution! But what shall" << endl;
3678 cout << "EdbShowerAlgESimple::WriteNewRootFile I do?" << endl;
3679 cout << "EdbShowerAlgESimple::WriteNewRootFile We rewrite completely the new tree and add " << endl;
3680 cout << "EdbShowerAlgESimple::WriteNewRootFile the new calculated energy values to it..." << endl;
3681 cout << "EdbShowerAlgESimple::WriteNewRootFile -------------------------------------------"<< endl;
3682 cout << endl;
3683
3684 //if (sourcefilename!="Shower.root") cout << "EdbShowerAlgESimple::WriteNewRootFile Attention: sourcefilename is different: "<< sourcefilename << endl;
3685 //if (treename!="treebranch") cout << "EdbShowerAlgESimple::WriteNewRootFile Attention: treename is different: "<< treename << endl;
3686
3687 cout << "EdbShowerAlgESimple::WriteNewRootFile Open sourcefilename:: "<< sourcefilename.Data() << endl;
3688 cout << "EdbShowerAlgESimple::WriteNewRootFile Open treename:: "<< treename.Data() << endl;
3689
3690 //- VARIABLES: shower_ "treebranch" reconstruction
3694 Float_t shower_xb[10000];
3695 Float_t shower_yb[10000];
3696 Float_t shower_zb[10000];
3697 Float_t shower_txb[10000];
3698 Float_t shower_tyb[10000];
3699 Float_t shower_deltarb[10000];
3700 Float_t shower_deltathetab[10000];
3701 Float_t shower_deltaxb[10000];
3702 Float_t shower_deltayb[10000];
3703 Int_t shower_nfilmb[10000];
3704 Float_t shower_chi2btkb[10000];
3705 Int_t shower_ntrace1simub[10000]; // MCEvt
3706 Int_t shower_ntrace2simub[10000]; // s->W()
3707 Float_t shower_ntrace3simub[10000]; // s->P()
3708 Int_t shower_ntrace4simub[10000]; // s->Flag()
3709 Float_t shower_tagprimary[10000];
3710 Int_t shower_idb[10000];
3711 Int_t shower_plateb[10000];
3712 Float_t shower_deltasigmathetab[58];
3713 Int_t shower_numberofilms;
3714 Float_t shower_purb; // purity of shower
3715 Int_t shower_eProb90;
3716 Int_t shower_eProb1;
3717 Int_t shower_Size; // number of BT in the shower
3718 Int_t shower_Size15; // number of BT in the shower (for 15 films crossed)
3719 Int_t shower_Size20; // number of BT in the shower (for 20 films crossed)
3720 Int_t shower_Size30; // number of BT in the shower (for 30 films crossed)
3721 Float_t shower_output; // Neural Network output for e/pi separation
3722 Float_t shower_output15; // Neural Network output for e/pi separation (for 15 films crossed)
3723 Float_t shower_output20; // Neural Network output for e/pi separation (for 20 films crossed)
3724 Float_t shower_output30; // Neural Network output for e/pi separation (for 30 films crossed)
3725 Float_t shower_output50; // Neural Network output for e/pi separation (for 50 films crossed)
3726 Float_t shower_purityb;
3727 Float_t shower_trackdensb;
3728 Float_t shower_E_MC;
3729 Float_t shower_EnergyCorrectedb;
3730 Float_t shower_EnergyUnCorrectedb;
3731 Float_t shower_EnergySigmaCorrectedb;
3732 Float_t shower_EnergySigmaUnCorrectedb;
3733
3734 Float_t shower_OldEnergyCorrectedb;
3735 Float_t shower_OldEnergyUnCorrectedb;
3736 Float_t shower_OldEnergySigmaCorrectedb;
3737 Float_t shower_OldEnergySigmaUnCorrectedb;
3738
3739 Int_t shower_mcDigitIndexTop[1000];
3740 Int_t shower_mcDigitIndexBottom[1000];
3741 //- Old Shower File:
3742 TFile* fileOld = new TFile(sourcefilename.Data(),"READ");
3743 // Set Addresses of treebranch tree:
3744 TTree* eShowerTree = (TTree*)fileOld->Get(treename);
3745 cout << "eShowerTree = " << eShowerTree << endl;
3746 eShowerTree->SetBranchAddress("number_eventb",&shower_number_eventb);
3747 eShowerTree->SetBranchAddress("sizeb",&shower_sizeb);
3748 eShowerTree->SetBranchAddress("sizeb15",&shower_sizeb15);
3749 eShowerTree->SetBranchAddress("sizeb20",&shower_sizeb20);
3750 eShowerTree->SetBranchAddress("sizeb30",&shower_sizeb30);
3752 eShowerTree->SetBranchAddress("output15", &shower_output15);
3753 eShowerTree->SetBranchAddress("output20", &shower_output20);
3754 eShowerTree->SetBranchAddress("output30", &shower_output30);
3755 eShowerTree->SetBranchAddress("output50", &shower_output50);
3756 eShowerTree->SetBranchAddress("eProb1", &shower_eProb1);
3757 eShowerTree->SetBranchAddress("eProb90", &shower_eProb90);
3758 eShowerTree->SetBranchAddress("isizeb",&shower_isizeb);
3759 eShowerTree->SetBranchAddress("xb",shower_xb);
3760 eShowerTree->SetBranchAddress("yb",shower_yb);
3761 eShowerTree->SetBranchAddress("zb",shower_zb);
3762 eShowerTree->SetBranchAddress("txb",shower_txb);
3763 eShowerTree->SetBranchAddress("tyb",shower_tyb);
3764 eShowerTree->SetBranchAddress("nfilmb",shower_nfilmb);
3765 eShowerTree->SetBranchAddress("ntrace1simub",shower_ntrace1simub); // s.eMCEvt
3766 eShowerTree->SetBranchAddress("ntrace2simub",shower_ntrace2simub); // s.eW
3767 eShowerTree->SetBranchAddress("ntrace3simub",shower_ntrace3simub); // s.eP
3768 eShowerTree->SetBranchAddress("ntrace4simub",shower_ntrace4simub); // s.eFlag
3769 eShowerTree->SetBranchAddress("chi2btkb",shower_chi2btkb);
3770 eShowerTree->SetBranchAddress("deltarb",shower_deltarb);
3771 eShowerTree->SetBranchAddress("deltathetab",shower_deltathetab);
3772 eShowerTree->SetBranchAddress("deltaxb",shower_deltaxb);
3773 eShowerTree->SetBranchAddress("deltayb",shower_deltayb);
3774 eShowerTree->SetBranchAddress("tagprimary",shower_tagprimary);
3775// eShowerTree->SetBranchAddress("energy_shot_particle",&shower_energy_shot_particle);
3776 eShowerTree->SetBranchAddress("E_MC",&shower_E_MC);
3777 eShowerTree->SetBranchAddress("showerID",&shower_showerID);
3778 eShowerTree->SetBranchAddress("idb",shower_idb);
3779 eShowerTree->SetBranchAddress("plateb",shower_plateb);
3780// eShowerTree->SetBranchAddress("deltasigmathetab",shower_deltasigmathetab);
3781// eShowerTree->SetBranchAddress("lenghtfilmb",&shower_numberofilms);
3782 eShowerTree->SetBranchAddress("purityb",&shower_purityb);
3783 eShowerTree->SetBranchAddress("Energy",&shower_OldEnergyCorrectedb);
3784 eShowerTree->SetBranchAddress("EnergyUnCorrected",&shower_OldEnergyUnCorrectedb);
3785 eShowerTree->SetBranchAddress("EnergySigma",&shower_OldEnergySigmaCorrectedb);
3786 eShowerTree->SetBranchAddress("EnergySigmaUnCorrected",&shower_OldEnergySigmaUnCorrectedb);
3787 eShowerTree->SetBranchAddress("mcDigitIndexTop",shower_mcDigitIndexTop);
3788 eShowerTree->SetBranchAddress("mcDigitIndexBottom",shower_mcDigitIndexBottom);
3789
3790
3791 Int_t nent=eShowerTree->GetEntries();
3792 cout << "EdbShowerAlgESimple::WriteNewRootFile eShowerTree="<< eShowerTree <<" ------"<< endl;
3793 cout << "EdbShowerAlgESimple::WriteNewRootFile eShowerTree->GetEntries()="<< nent <<" ------"<< endl;
3794
3795
3796 // Create the new File
3797 TFile* fileNew = new TFile("New.root","RECREATE");
3798 // Create the new Tree
3799 TTree* ShowerTreeNew=new TTree("treebranch","treebranch");
3800 cout << "EdbShowerAlgESimple::WriteNewRootFile ShowerTreeNew="<< ShowerTreeNew <<" ------"<< endl;
3801 // Create the new Branches:
3802 ShowerTreeNew->Branch("number_eventb",&shower_number_eventb,"number_eventb/I");
3803 ShowerTreeNew->Branch("sizeb",&shower_sizeb,"sizeb/I");
3804 ShowerTreeNew->Branch("sizeb15",&shower_sizeb15,"sizeb15/I");
3805 ShowerTreeNew->Branch("sizeb20",&shower_sizeb20,"sizeb20/I");
3806 ShowerTreeNew->Branch("sizeb30",&shower_sizeb30,"sizeb30/I");
3807 ShowerTreeNew->Branch("output15",&shower_output15,"output15/F");
3808 ShowerTreeNew->Branch("output20",&shower_output20,"output20/F");
3809 ShowerTreeNew->Branch("output30",&shower_output30,"output30/F");
3810 ShowerTreeNew->Branch("output50",&shower_output50,"output50/F");
3811 ShowerTreeNew->Branch("eProb90",&shower_eProb90,"eProb90/I");
3812 ShowerTreeNew->Branch("eProb1",&shower_eProb1,"eProb1/I");
3813 ShowerTreeNew->Branch("E_MC",&shower_E_MC,"E_MC/F");
3814 ShowerTreeNew->Branch("idb",shower_idb,"idb[sizeb]/I");
3815 ShowerTreeNew->Branch("plateb",shower_plateb,"plateb[sizeb]/I");
3816 ShowerTreeNew->Branch("showerID",&shower_showerID,"showerID/I");
3817 ShowerTreeNew->Branch("isizeb",&shower_isizeb,"isizeb/I");
3818 ShowerTreeNew->Branch("xb",shower_xb,"xb[sizeb]/F");
3819 ShowerTreeNew->Branch("yb",shower_yb,"yb[sizeb]/F");
3820 ShowerTreeNew->Branch("zb",shower_zb,"zb[sizeb]/F");
3821 ShowerTreeNew->Branch("txb",shower_txb,"txb[sizeb]/F");
3822 ShowerTreeNew->Branch("tyb",shower_tyb,"tyb[sizeb]/F");
3823 ShowerTreeNew->Branch("nfilmb",shower_nfilmb,"nfilmb[sizeb]/I");
3824 ShowerTreeNew->Branch("ntrace1simub",shower_ntrace1simub,"ntrace1simu[sizeb]/I");
3825 ShowerTreeNew->Branch("ntrace2simub",shower_ntrace2simub,"ntrace2simu[sizeb]/I");
3826 ShowerTreeNew->Branch("ntrace3simub",shower_ntrace3simub,"ntrace3simu[sizeb]/F");
3827 ShowerTreeNew->Branch("ntrace4simub",shower_ntrace4simub,"ntrace4simu[sizeb]/I");
3828 ShowerTreeNew->Branch("chi2btkb",shower_chi2btkb,"chi2btkb[sizeb]/F");
3829 ShowerTreeNew->Branch("deltarb",shower_deltarb,"deltarb[sizeb]/F");
3830 ShowerTreeNew->Branch("deltathetab",shower_deltathetab,"deltathetab[sizeb]/F");
3831 ShowerTreeNew->Branch("deltaxb",shower_deltaxb,"deltaxb[sizeb]/F");
3832 ShowerTreeNew->Branch("deltayb",shower_deltayb,"deltayb[sizeb]/F");
3833 ShowerTreeNew->Branch("tagprimary",shower_tagprimary,"tagprimary[sizeb]/F");
3834 ShowerTreeNew->Branch("purityb",&shower_purityb,"purityb/F");
3835 ShowerTreeNew->Branch("trackdensb",&shower_trackdensb,"trackdensb/F");
3836
3837 ShowerTreeNew->Branch("Energy",&shower_EnergyCorrectedb,"EnergyCorrectedb/F");
3838 ShowerTreeNew->Branch("EnergyUnCorrected",&shower_EnergyUnCorrectedb,"EnergyUnCorrectedb/F");
3839 ShowerTreeNew->Branch("EnergySigma",&shower_EnergySigmaCorrectedb,"EnergySigmaCorrectedb/F");
3840 ShowerTreeNew->Branch("EnergySigmaUnCorrected",&shower_EnergySigmaUnCorrectedb,"EnergySigmaUnCorrectedb/F");
3841
3842 ShowerTreeNew->Branch("OldEnergy",&shower_OldEnergyCorrectedb,"OldEnergy/F");
3843 ShowerTreeNew->Branch("OldEnergyUnCorrected",&shower_OldEnergyUnCorrectedb,"OldEnergyUnCorrected/F");
3844 ShowerTreeNew->Branch("OldEnergySigma",&shower_OldEnergySigmaCorrectedb,"OldEnergySigma/F");
3845 ShowerTreeNew->Branch("OldEnergySigmaUnCorrected",&shower_OldEnergySigmaUnCorrectedb,"OldEnergySigmaUnCorrected/F");
3846 ShowerTreeNew->Branch("mcDigitIndexTop",shower_mcDigitIndexTop,"mcDigitIndexTop[sizeb]/I");
3847 ShowerTreeNew->Branch("mcDigitIndexBottom",shower_mcDigitIndexBottom,"mcDigitIndexBottom[sizeb]/I");
3848
3849
3850 // Loop over Tree Entries (==different showers):
3851 for (int i=0; i<nent; ++i) {
3852 eShowerTree->GetEntry(i);
3853 //cout << "i = " << i << " " << eEnergyArrayCount << " shower_EnergyCorrectedb = " << eEnergyArray->At(i) << endl;
3854// eShowerTree->Show(i);
3855 shower_EnergyCorrectedb=eEnergyArray->At(i);
3856 shower_EnergyUnCorrectedb=eEnergyArrayUnCorrected->At(i);
3857 shower_EnergySigmaCorrectedb=eEnergyArraySigmaCorrected->At(i);
3858 shower_EnergySigmaUnCorrectedb=-999;
3859
3860 // Fill new Tree
3861 ShowerTreeNew->Fill();
3862 }
3863 ShowerTreeNew->Write();
3864 fileNew->Close();
3865 fileOld->Close();
3866
3867
3868 // This is realy not nice coding ....
3869 gSystem->Exec("mv -vf Shower.root Shower.Orig.root");
3870 gSystem->Exec("mv -vf New.root Shower.root");
3871 cout << "EdbShowerAlgESimple::WriteNewRootFile...done."<<endl;
3872 cout << "-------------------------------------------------------------------------------------------"<< endl<< endl;
3873 return;
3874}
Int_t shower_nfilmb[5000]
Definition: ShowRec.h:382
Int_t shower_showerID
Definition: ShowRec.h:370
Float_t shower_ntrace3simub[5000]
Definition: ShowRec.h:387
Float_t shower_deltaxb[5000]
Definition: ShowRec.h:380
Int_t shower_sizeb20
Definition: ShowRec.h:371
Float_t shower_purb
Definition: ShowRec.h:394
Int_t shower_numberofilms
Definition: ShowRec.h:393
Int_t shower_ntrace2simub[5000]
Definition: ShowRec.h:386
Float_t shower_energy_shot_particle
Definition: ShowRec.h:372
Float_t shower_trackdensb
Definition: ShowRec.h:395
Int_t shower_idb[5000]
Definition: ShowRec.h:390
Int_t shower_number_eventb
Definition: ShowRec.h:370
Int_t shower_sizeb15
Definition: ShowRec.h:371
Float_t shower_deltasigmathetab[58]
Definition: ShowRec.h:392
Float_t shower_tagprimary[5000]
Definition: ShowRec.h:389
Int_t shower_isizeb
Definition: ShowRec.h:370
Float_t shower_deltayb[5000]
Definition: ShowRec.h:381
Int_t shower_sizeb30
Definition: ShowRec.h:371
Int_t shower_ntrace4simub[5000]
Definition: ShowRec.h:388
Int_t shower_sizeb
Definition: ShowRec.h:370
Float_t shower_chi2btkb[5000]
Definition: ShowRec.h:384
Int_t shower_plateb[5000]
Definition: ShowRec.h:391
Int_t shower_ntrace1simub[5000]
Definition: ShowRec.h:385

Member Data Documentation

◆ ANN_Layout

TString EdbShowerAlgESimple::ANN_Layout
protected

◆ ANN_MLP

TMultiLayerPerceptron* EdbShowerAlgESimple::ANN_MLP
protected

◆ ANN_MLP_ARRAY

TMultiLayerPerceptron* EdbShowerAlgESimple::ANN_MLP_ARRAY[15]
protected

◆ ANN_n_InputNeurons

Int_t EdbShowerAlgESimple::ANN_n_InputNeurons
protected

◆ ANN_n_InputNeurons_ARRAY

Int_t EdbShowerAlgESimple::ANN_n_InputNeurons_ARRAY[15]
protected

◆ ANN_nPlates_ARRAY

Int_t EdbShowerAlgESimple::ANN_nPlates_ARRAY[15]
protected

◆ ANN_WeightFile_ARRAY

TString EdbShowerAlgESimple::ANN_WeightFile_ARRAY[15]
protected

◆ ANNTree

TTree* EdbShowerAlgESimple::ANNTree
protected

◆ eANN_MLP_CORR_0

Float_t EdbShowerAlgESimple::eANN_MLP_CORR_0[15]
protected

◆ eANN_MLP_CORR_1

Float_t EdbShowerAlgESimple::eANN_MLP_CORR_1[15]
protected

◆ eCalibrationOffset

Float_t EdbShowerAlgESimple::eCalibrationOffset
protected

◆ eCalibrationSlope

Float_t EdbShowerAlgESimple::eCalibrationSlope
protected

◆ eEfficiencyParametrisation

TSpline3* EdbShowerAlgESimple::eEfficiencyParametrisation
protected

◆ eEnergy

Float_t EdbShowerAlgESimple::eEnergy
protected

◆ eEnergyArray

TArrayF* EdbShowerAlgESimple::eEnergyArray
protected

◆ eEnergyArrayCount

Int_t EdbShowerAlgESimple::eEnergyArrayCount
protected

◆ eEnergyArraySigmaCorrected

TArrayF* EdbShowerAlgESimple::eEnergyArraySigmaCorrected
protected

◆ eEnergyArrayUnCorrected

TArrayF* EdbShowerAlgESimple::eEnergyArrayUnCorrected
protected

◆ eEnergyCorr

Float_t EdbShowerAlgESimple::eEnergyCorr
protected

◆ eEnergySigmaCorr

Float_t EdbShowerAlgESimple::eEnergySigmaCorr
protected

◆ eEnergyUnCorr

Float_t EdbShowerAlgESimple::eEnergyUnCorr
protected

◆ EffFunc_all

TF1* EdbShowerAlgESimple::EffFunc_all
protected

◆ EffFunc_edefault

TF1* EdbShowerAlgESimple::EffFunc_edefault
protected

◆ EffFunc_elletroni

TF1* EdbShowerAlgESimple::EffFunc_elletroni
protected

◆ EffFunc_LowEff

TF1* EdbShowerAlgESimple::EffFunc_LowEff
protected

◆ EffFunc_MiddleFix

TF1* EdbShowerAlgESimple::EffFunc_MiddleFix
protected

◆ EffFunc_neuchmicro

TF1* EdbShowerAlgESimple::EffFunc_neuchmicro
protected

◆ EffFunc_UserEfficiency

TF1* EdbShowerAlgESimple::EffFunc_UserEfficiency
protected

◆ eForceSpecificationReload

Bool_t EdbShowerAlgESimple::eForceSpecificationReload
protected

◆ eHisto_deltaR

TH1D* EdbShowerAlgESimple::eHisto_deltaR
protected

◆ eHisto_deltaR_mean

TH1D* EdbShowerAlgESimple::eHisto_deltaR_mean
protected

◆ eHisto_deltaR_rms

TH1D* EdbShowerAlgESimple::eHisto_deltaR_rms
protected

◆ eHisto_deltaT

TH1D* EdbShowerAlgESimple::eHisto_deltaT
protected

◆ eHisto_deltaT_mean

TH1D* EdbShowerAlgESimple::eHisto_deltaT_mean
protected

◆ eHisto_deltaT_rms

TH1D* EdbShowerAlgESimple::eHisto_deltaT_rms
protected

◆ eHisto_longprofile

TH1D* EdbShowerAlgESimple::eHisto_longprofile
protected

◆ eHisto_longprofile_av

TH1D* EdbShowerAlgESimple::eHisto_longprofile_av
protected

◆ eHisto_nbtk

TH1D* EdbShowerAlgESimple::eHisto_nbtk
protected

◆ eHisto_nbtk_av

TH1D* EdbShowerAlgESimple::eHisto_nbtk_av
protected

◆ eHisto_transprofile

TH1D* EdbShowerAlgESimple::eHisto_transprofile
protected

◆ eHisto_transprofile_av

TH1D* EdbShowerAlgESimple::eHisto_transprofile_av
protected

◆ eParaBT_deltaR_mean

Float_t EdbShowerAlgESimple::eParaBT_deltaR_mean
protected

◆ eParaBT_deltaR_rms

Float_t EdbShowerAlgESimple::eParaBT_deltaR_rms
protected

◆ eParaBT_deltaT_mean

Float_t EdbShowerAlgESimple::eParaBT_deltaT_mean
protected

◆ eParaBT_deltaT_rms

Float_t EdbShowerAlgESimple::eParaBT_deltaT_rms
protected

◆ eParalongprofile

Int_t EdbShowerAlgESimple::eParalongprofile[57]
protected

◆ eParaName

Int_t EdbShowerAlgESimple::eParaName
protected

◆ eParanseg

Int_t EdbShowerAlgESimple::eParanseg
protected

◆ eParaShowerAxisAngle

Float_t EdbShowerAlgESimple::eParaShowerAxisAngle
protected

◆ ePlateBinning

Int_t EdbShowerAlgESimple::ePlateBinning[15]
protected

◆ ePlateNumber

Int_t EdbShowerAlgESimple::ePlateNumber
protected

◆ ePlateNumberType

Int_t EdbShowerAlgESimple::ePlateNumberType
protected

◆ eRecoShowerArray

TObjArray* EdbShowerAlgESimple::eRecoShowerArray
protected

◆ eRecoShowerArrayN

Int_t EdbShowerAlgESimple::eRecoShowerArrayN
protected

◆ eSpecificationIsChanged

Bool_t EdbShowerAlgESimple::eSpecificationIsChanged
protected

◆ eSpecificationType

Int_t EdbShowerAlgESimple::eSpecificationType[7]
protected

◆ eSpecificationTypeString

TString EdbShowerAlgESimple::eSpecificationTypeString[7]
protected

◆ eSpecificationTypeStringArray

TString EdbShowerAlgESimple::eSpecificationTypeStringArray[7][7]
protected

◆ eSplineArray_Energy_Stat_Electron

TObjArray* EdbShowerAlgESimple::eSplineArray_Energy_Stat_Electron
protected

◆ eSplineArray_Energy_Stat_Gamma

TObjArray* EdbShowerAlgESimple::eSplineArray_Energy_Stat_Gamma
protected

◆ eSplineArray_Energy_Sys_Electron

TObjArray* EdbShowerAlgESimple::eSplineArray_Energy_Sys_Electron
protected

◆ eSplineArray_Energy_Sys_Gamma

TObjArray* EdbShowerAlgESimple::eSplineArray_Energy_Sys_Gamma
protected

◆ eSplineCurrent

TSpline3* EdbShowerAlgESimple::eSplineCurrent
protected

◆ eWeightFileString

TString EdbShowerAlgESimple::eWeightFileString
protected

◆ inANN

Double_t EdbShowerAlgESimple::inANN[70]
protected

◆ outANN

Double_t EdbShowerAlgESimple::outANN
protected

The documentation for this class was generated from the following files: