FEDRA emulsion software from the OPERA Collaboration
EdbSEQ Class Reference

#include <EdbSEQ.h>

Inheritance diagram for EdbSEQ:
Collaboration diagram for EdbSEQ:

Public Member Functions

void AddExcludeThetaRange (EdbSegP &s)
 
void CalculateDensityMT (EdbH1 &hEq)
 
double DNbt (double t)
 
double DNmt (double t)
 
void Draw ()
 
 EdbSEQ ()
 
void EqualizeMT (TObjArray &mti, TObjArray &mto, Double_t area)
 
void ExcludeThetaRange (TObjArray &mti, TObjArray &mto)
 
double FDNbt (double *x, double *par)
 
double FDNmt (double *x, double *par)
 
bool IsInsideThetaRange (EdbSegP *s)
 
void PreSelection (EdbPattern &pi, TObjArray &po)
 
void PrintLimits ()
 
void ResetExcludeThetaRange ()
 
void Set0 ()
 
void SetChiLimits (float cmin, float cmax)
 
void SetThetaLimits (float tmin, float tmax)
 
void SetWLimits (float wmin, float wmax)
 
void SetXLimits (float xmin, float xmax)
 
void SetYLimits (float ymin, float ymax)
 
TH1F * ThetaPlot (EdbPattern &p, const char *name="theta", const char *title="EdbSEQ theta distribution normalised to area")
 
TH1F * ThetaPlot (TObjArray &arr, const char *name="theta", const char *title="EdbSEQ theta distribution normalised to area")
 
float Wmt (EdbSegP &s)
 
virtual ~EdbSEQ ()
 
- Public Member Functions inherited from EdbSigma
double DAL (double t, double sxy, double sz, double dz)
 
double DALbt (double t)
 
double DALmt (double t)
 
double DAT (double t, double sxy, double dz)
 
double DATbt (double t)
 
double DATmt (double t)
 
double DP (double t, double sxy, double da, double dz)
 
double DPLbt (double t)
 
double DPLmt (double t)
 
double DPTbt (double t)
 
double DPTmt (double t)
 
void Draw ()
 
 EdbSigma ()
 RS agreement: 0 - in-view (default); 1 - in-zone; 2-in-plate; 3-in-brick. More...
 
double FDAL (double *x, double *par)
 
double FDALbt (double *x, double *par)
 
double FDAT (double *x, double *par)
 
double FDPLmt (double *x, double *par)
 
void Set0 ()
 
double SqSum (double a, double b)
 
virtual ~EdbSigma ()
 

Public Attributes

Double_t eArea
 effective area of the pattern to be equalized More...
 
TObjArray eExcludeThetaRange
 can be added EdbSegP with tx,ty, sigmaTX,sigmaTY to be excluded More...
 
EdbH1 eHEq
 
int eNP
 number of points for the functions calculation More...
 
Double_t eNsigma
 =4 More...
 
Double_t eS0mt
 = 270.*340.; // area unit for Nseg calculation More...
 
- Public Attributes inherited from EdbSigma
double eDZbase
 = 210. emulsion base More...
 
double eDZcell
 =1300. plate-to plate distance More...
 
double eDZem
 geometry defines the "in-view" mt/bt accuracy: More...
 
double eSaPlate
 "in-plate" basetrack angular tolrance (nonplanarity) More...
 
double eSaZone
 "in-zone" microtrack angular tolrance (incorrect shrinkage+offsets) More...
 
double eSxy
 grains measurement accuracy defines the "in view" bt and mt errors More...
 
double eSxyPlate
 additional terms defines "in-plate" resolution More...
 
double eSxyZone
 additional terms defines "in-zone" resolution More...
 
double eSz
 = 2. grains z-uncertainty "field depth" More...
 

Private Attributes

TVector2 * eChiLimits
 [min,max] More...
 
TVector2 * eThetaLimits
 [min,max] More...
 
TVector2 * eWLimits
 [min,max] More...
 
TVector2 * eXLimits
 [min,max] area limits for the preselection procedure More...
 
TVector2 * eYLimits
 [min,max] More...
 

Constructor & Destructor Documentation

◆ EdbSEQ()

EdbSEQ::EdbSEQ ( )
inline
32{ ((EdbSigma*)this)->Set0(); Set0();}
void Set0()
Definition: EdbSEQ.cxx:33
Definition: EdbSigma.h:8

◆ ~EdbSEQ()

EdbSEQ::~EdbSEQ ( )
virtual
24{
25 SafeDelete(eXLimits);
26 SafeDelete(eYLimits);
27 SafeDelete(eWLimits);
28 SafeDelete(eThetaLimits);
29 SafeDelete(eChiLimits);
30}
TVector2 * eXLimits
[min,max] area limits for the preselection procedure
Definition: EdbSEQ.h:25
TVector2 * eChiLimits
[min,max]
Definition: EdbSEQ.h:29
TVector2 * eWLimits
[min,max]
Definition: EdbSEQ.h:28
TVector2 * eYLimits
[min,max]
Definition: EdbSEQ.h:26
TVector2 * eThetaLimits
[min,max]
Definition: EdbSEQ.h:27

Member Function Documentation

◆ AddExcludeThetaRange()

void EdbSEQ::AddExcludeThetaRange ( EdbSegP s)
190{
191 eExcludeThetaRange.Add( new EdbSegP(s) );
192}
TObjArray eExcludeThetaRange
can be added EdbSegP with tx,ty, sigmaTX,sigmaTY to be excluded
Definition: EdbSEQ.h:21
Definition: EdbSegP.h:21
s
Definition: check_shower.C:55

◆ CalculateDensityMT()

void EdbSEQ::CalculateDensityMT ( EdbH1 hEq)
138{
139 Log(3,"EdbSEQ::CalculateDensityMT","eHEq.N(): %d", eHEq.N());
140 TF1 *fnmt = new TF1("fnmt",this,&EdbSEQ::FDNmt,eHEq.Xmin(),eHEq.Xmax(),0,"EdbSEQ","FDNmt");
141 fnmt->SetTitle("dNmt/dTheta in 1 view density function");
142 fnmt->SetNpx(eNP);
143 double *x=new double[eNP];
144 double *w=new double[eNP];
145 fnmt->CalcGaussLegendreSamplingPoints(eNP,x,w,1e-15);
146 int n = eHEq.N();
147 for(int i=0; i<n; i++) {
148 float tmi=i*eHEq.Xbin();
149 float tma=(i+1)*eHEq.Xbin();
150 double sbin = fnmt->IntegralFast(eNP,x,w,tmi,tma);
151 Log(3,"EdbSEQ::CalculateDensityMT","i,eNP,x,w,tmi,tma: %d %d %f %f %f %f %lf", i,eNP,x,w,tmi,tma, sbin );
152 if(sbin>0 && sbin<kMaxInt-2) eHEq.SetBin(i, int(sbin)+1 );
153 //Log(3,"EdbSEQ::CalculateDensityMT","Integral(%f,%f) = %f", tmi, tma, float(sbin) );
154 }
155 Log(3,"EdbSEQ::CalculateDensityMT","Integral(0,0.6) = %f", float(fnmt->IntegralFast(eNP,x,w,0,0.6)) );
156}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
float Xbin() const
Definition: EdbCell1.h:56
int N() const
Definition: EdbCell1.h:48
void SetBin(int ix, int n)
Definition: EdbCell1.h:65
float Xmax() const
Definition: EdbCell1.h:55
float Xmin() const
Definition: EdbCell1.h:54
int eNP
number of points for the functions calculation
Definition: EdbSEQ.h:19
EdbH1 eHEq
Definition: EdbSEQ.h:22
double FDNmt(double *x, double *par)
Definition: EdbSEQ.h:44
void w(int rid=2, int nviews=2)
Definition: test.C:27

◆ DNbt()

double EdbSEQ::DNbt ( double  t)

return dN/dTheta : the critical number of microtracks for the given theta range for a view area printf("t=%f\n",t);

98{
102 double dal0 = SqSum(DALbt(0),eSaPlate);
103 double dal = SqSum(DALbt(t),eSaPlate);
104 double dat = SqSum(DATbt(t),eSaPlate);
105 double dpl = DP( t, SqSum( DPLbt(t),eSxyPlate ), dal, eDZcell );
106 double dpt = DP( t, SqSum( DPTbt(t),eSxyPlate ), dat, eDZcell );
107 double nang = Sin(t)/dat * dal0/dal;
108 return nang * eArea/dpl/dpt /eNsigma/eNsigma/eNsigma/eNsigma;
109}
Double_t eArea
effective area of the pattern to be equalized
Definition: EdbSEQ.h:17
Double_t eNsigma
=4
Definition: EdbSEQ.h:16
double eDZcell
=1300. plate-to plate distance
Definition: EdbSigma.h:18
double DPLbt(double t)
Definition: EdbSigma.cxx:102
double eSaPlate
"in-plate" basetrack angular tolrance (nonplanarity)
Definition: EdbSigma.h:26
double DALbt(double t)
Definition: EdbSigma.cxx:88
double DATbt(double t)
Definition: EdbSigma.cxx:95
double DP(double t, double sxy, double da, double dz)
Definition: EdbSigma.cxx:51
double DPTbt(double t)
Definition: EdbSigma.cxx:110
double eSxyPlate
additional terms defines "in-plate" resolution
Definition: EdbSigma.h:25
double SqSum(double a, double b)
Definition: EdbSigma.h:36
TTree * t
Definition: check_shower.C:4

◆ DNmt()

double EdbSEQ::DNmt ( double  t)

return dN/dTheta : the critical number of microtracks for the given theta range for a view*nviews area

85{
88
89 double dpl = DP(t,DPLmt(t),DALmt(t),eDZbase);
90 double dpt = DP(t,DPTmt(t),DATmt(t),eDZbase);
91 double nang = Sin(t)/DATmt(t) * DALmt(0)/DALmt(t);
92 //return nang * eNviews*eS0mt/dpl/dpt /eNsigma/eNsigma/eNsigma/eNsigma;
93 return nang * eArea/dpl/dpt /eNsigma/eNsigma/eNsigma/eNsigma;
94}
double DATmt(double t)
Definition: EdbSigma.cxx:66
double DPLmt(double t)
Definition: EdbSigma.cxx:73
double DPTmt(double t)
Definition: EdbSigma.cxx:80
double eDZbase
= 210. emulsion base
Definition: EdbSigma.h:17
double DALmt(double t)
Definition: EdbSigma.cxx:59

◆ Draw()

void EdbSEQ::Draw ( )
113{
114 TCanvas *cn = new TCanvas("cn","Critical density",800,600);
115 cn->Divide(2,2);
116
117 TF1 *fnmt = new TF1("fnmt",this,&EdbSEQ::FDNmt,0.,1.,0,"EdbSEQ","FDNmt");
118 fnmt->SetTitle("dNmt/dTheta in 1 view density function");
119 fnmt->SetNpx(eNP);
120
121 cn->cd(1);
122 fnmt->Draw();
123 cn->cd(2);
124 fnmt->DrawIntegral();
125
126 TF1 *fnbt = new TF1("fnbt",this,&EdbSEQ::FDNbt,0.,1.,0,"EdbSEQ","FDNbt");
127 fnbt->SetTitle("dNbt/dTheta in 1 view density function");
128 fnbt->SetNpx(eNP);
129
130 cn->cd(3);
131 fnbt->Draw();
132 cn->cd(4);
133 fnbt->DrawIntegral();
134}
double FDNbt(double *x, double *par)
Definition: EdbSEQ.h:47
new TCanvas()

◆ EqualizeMT()

void EdbSEQ::EqualizeMT ( TObjArray &  mti,
TObjArray &  mto,
Double_t  area 
)

Input: mti - microtracks array Output: mto - selected microtracks array

221{
224
225 eArea=area;
226 float tbin=0.02;
227 int nbin= (int)(eThetaLimits->Y()/tbin) + 1;
228 eHEq.InitH1( nbin, 0, (nbin-1)*tbin );
229 CalculateDensityMT(eHEq); // define critical theta density
230 if(gEDBDEBUGLEVEL>2) eHEq.Print();
231
232 TClonesArray hw("EdbH1",eHEq.N()); // create W-hist per each theta bin
233 TArrayF thres(eHEq.N()); // W-thresholds to be defined
234 for(int i=0; i<eHEq.N(); i++) {
235 new(hw[i]) EdbH1();
236 ((EdbH1*)hw[i])->InitH1( 500, 5.*50, 17.*50. );
237 }
238
239 int nseg = mti.GetEntriesFast();
240 for(int i=0; i<nseg; i++) {
241 EdbSegP *s = (EdbSegP*)(mti.UncheckedAt(i));
242 int jtheta = eHEq.IX(s->Theta());
243 if(jtheta>eHEq.N()-1) continue;
244 float w = Wmt(*s);
245 ((EdbH1*)hw[jtheta])->Fill(w);
246 }
247
248 for(int i=0; i<eHEq.N(); i++) {
249 int eccess = ((EdbH1*)hw[i])->Integral() - eHEq.Bin(i);
250 if( eccess <= 0 ) continue; // the threshold remains 0
251
252 int nsum=0;
253 for(int j=0; j<((EdbH1*)hw[i])->N(); j++) {
254 nsum+=((EdbH1*)hw[i])->Bin(j);
255 if( nsum >= eccess ) { thres[i] = ((EdbH1*)hw[i])->X(j)-((EdbH1*)hw[i])->Xbin(); break; }
256 }
257 }
258
259 for(int i=0; i<nseg; i++) {
260 EdbSegP *s = (EdbSegP*)(mti.UncheckedAt(i));
261 int jtheta = eHEq.IX(s->Theta());
262 if(jtheta>eHEq.N()-1) continue;
263 float w = Wmt(*s);
264 if( w > thres[jtheta] )
265 if(!IsInsideThetaRange(s))
266 mto.Add(s);
267 }
268
269 Log(2,"EdbSEQ::EqualizeMT","selected %3.1f\% (%d of %d) segments in area of %7.1f cm²", 100.*mto.GetEntriesFast()/nseg, mto.GetEntriesFast(), nseg, (float)(eArea/1000./1000./10./10.) );
270}
fast 2-dim histogram class (used as a basis for EdbCell1)
Definition: EdbCell1.h:17
int Bin(int ix) const
Definition: EdbCell1.h:57
void Print()
Definition: EdbCell1.cpp:161
int IX(float x) const
Definition: EdbCell1.h:50
int InitH1(const EdbH1 &h)
Definition: EdbCell1.h:38
void CalculateDensityMT(EdbH1 &hEq)
Definition: EdbSEQ.cxx:137
float Wmt(EdbSegP &s)
Definition: EdbSEQ.h:60
bool IsInsideThetaRange(EdbSegP *s)
Definition: EdbSEQ.cxx:206
gEDBDEBUGLEVEL
Definition: energy.C:7

◆ ExcludeThetaRange()

void EdbSEQ::ExcludeThetaRange ( TObjArray &  mti,
TObjArray &  mto 
)

exclude 1 sigma in theta around all requested segments

196{
198 int nseg = mti.GetEntriesFast();
199 for(int i=0; i<nseg; i++) {
200 EdbSegP *s = (EdbSegP*)(mti.UncheckedAt(i));
201 if(!IsInsideThetaRange(s)) mto.Add(s);
202 }
203}

◆ FDNbt()

double EdbSEQ::FDNbt ( double *  x,
double *  par 
)
inline
47{return DNbt(*x);}
double DNbt(double t)
Definition: EdbSEQ.cxx:97

◆ FDNmt()

double EdbSEQ::FDNmt ( double *  x,
double *  par 
)
inline
44{return DNmt(*x);}
double DNmt(double t)
Definition: EdbSEQ.cxx:84

◆ IsInsideThetaRange()

bool EdbSEQ::IsInsideThetaRange ( EdbSegP s)

exclude 1 sigma in theta around all requested segments

207{
209 int nrange = eExcludeThetaRange.GetEntriesFast();
210 for(int j=0; j<nrange; j++) {
211 EdbSegP *sex = (EdbSegP*)(eExcludeThetaRange.UncheckedAt(j));
212 float dtx = s->TX() - sex->TX();
213 float dty = s->TY() - sex->TY();
214 if( (dtx*dtx/sex->STX()/sex->STX()+dty*dty/sex->STY()/sex->STY())<1. ) return 1;
215 }
216 return 0;
217}
Float_t STX() const
Definition: EdbSegP.h:164
Float_t TX() const
tangens = deltaX/deltaZ
Definition: EdbSegP.h:175
Float_t TY() const
tangens = deltaY/deltaZ
Definition: EdbSegP.h:176
Float_t STY() const
Definition: EdbSegP.h:165

◆ PreSelection()

void EdbSEQ::PreSelection ( EdbPattern pi,
TObjArray &  po 
)
160{
161 int nseg = pi.N();
162 for(int i=0; i<nseg; i++) {
163 EdbSegP *s = pi.GetSegment(i);
164 if(eChiLimits) {
165 if( s->Chi2() < eChiLimits->X() ) continue;
166 if( s->Chi2() > eChiLimits->Y() ) continue;
167 }
168 if(eXLimits) {
169 if( s->X() < eXLimits->X() ) continue;
170 if( s->X() > eXLimits->Y() ) continue;
171 }
172 if(eYLimits) {
173 if( s->Y() < eYLimits->X() ) continue;
174 if( s->Y() > eYLimits->Y() ) continue;
175 }
176 if(eWLimits) {
177 if( s->W() < eWLimits->X() ) continue;
178 if( s->W() > eWLimits->Y() ) continue;
179 }
180 if(eThetaLimits) {
181 if( s->Theta() < eThetaLimits->X() ) continue;
182 if( s->Theta() > eThetaLimits->Y() ) continue;
183 }
184 po.Add(s);
185 }
186}
Int_t N() const
Definition: EdbPattern.h:86
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:66
TProfile * po
Definition: testChi2Ordering.C:29

◆ PrintLimits()

void EdbSEQ::PrintLimits ( )
50{
51 if(eXLimits) printf("XLimits: %f %f\n",eXLimits->X(),eXLimits->Y());
52 if(eYLimits) printf("YLimits: %f %f\n",eYLimits->X(),eYLimits->Y());
53 if(eWLimits) printf("WLimits: %f %f\n",eWLimits->X(),eWLimits->Y());
54 if(eChiLimits) printf("ChiLimits: %f %f\n",eChiLimits->X(),eChiLimits->Y());
55 if(eThetaLimits) printf("ThetaLimits: %f %f\n",eThetaLimits->X(),eThetaLimits->Y());
56}

◆ ResetExcludeThetaRange()

void EdbSEQ::ResetExcludeThetaRange ( )
inline
55{eExcludeThetaRange.Delete();}

◆ Set0()

void EdbSEQ::Set0 ( )
34{
35 eS0mt = 270.*340.; // area unit for Nseg calculation in microns
36 //eNviews = 1; // number of views used for total area calculation
37 eArea = eS0mt; // effective area of the pattern to be equalize
38 eNsigma=4;
39 eNP = 100;
40 eExcludeThetaRange.SetOwner(1);
41 eXLimits =0;
42 eYLimits =0;
44 eWLimits =0;
45 eChiLimits =0;
46}
Double_t eS0mt
= 270.*340.; // area unit for Nseg calculation
Definition: EdbSEQ.h:15

◆ SetChiLimits()

void EdbSEQ::SetChiLimits ( float  cmin,
float  cmax 
)
inline
41{SafeDelete(eChiLimits); eChiLimits=new TVector2(cmin,cmax);}

◆ SetThetaLimits()

void EdbSEQ::SetThetaLimits ( float  tmin,
float  tmax 
)
inline
40{SafeDelete(eThetaLimits); eThetaLimits=new TVector2(tmin,tmax);}

◆ SetWLimits()

void EdbSEQ::SetWLimits ( float  wmin,
float  wmax 
)
inline
39{SafeDelete(eWLimits); eWLimits=new TVector2(wmin,wmax);}

◆ SetXLimits()

void EdbSEQ::SetXLimits ( float  xmin,
float  xmax 
)
inline
37{SafeDelete(eXLimits); eXLimits=new TVector2(xmin,xmax);}

◆ SetYLimits()

void EdbSEQ::SetYLimits ( float  ymin,
float  ymax 
)
inline
38{SafeDelete(eYLimits); eYLimits=new TVector2(ymin,ymax);}

◆ ThetaPlot() [1/2]

TH1F * EdbSEQ::ThetaPlot ( EdbPattern p,
const char *  name = "theta",
const char *  title = "EdbSEQ theta distribution normalised to area" 
)
60{
61 int n = arr.N(); if(!n) return 0;
62 TObjArray a(n);
63 for(int i=0; i<n; i++) a.Add(arr.GetSegment(i));
64 return ThetaPlot(a, name, title);
65}
void a()
Definition: check_aligned.C:59
TH1F * ThetaPlot(TObjArray &arr, const char *name="theta", const char *title="EdbSEQ theta distribution normalised to area")
Definition: EdbSEQ.cxx:68
const char * name
Definition: merge_Energy_SytematicSources_Electron.C:24

◆ ThetaPlot() [2/2]

TH1F * EdbSEQ::ThetaPlot ( TObjArray &  arr,
const char *  name = "theta",
const char *  title = "EdbSEQ theta distribution normalised to area" 
)
69{
70 int n = arr.GetEntries(); if(!n) return 0;
71 float tbin=0.01;
72 int nbin= (int)(eThetaLimits->Y()/tbin) + 1;
73 TH1F *h = new TH1F( name, Form("%s normalised to area of %8.1f mm²",title,eArea/1000./1000.), nbin, 0., (nbin-1)*tbin );
74 for(int i=0; i<n; i++) {
75 EdbSegP *s = (EdbSegP*)(arr.At(i));
76 h->Fill(s->Theta());
77 }
78 h->GetXaxis()->SetTitle("theta");
79 if(eArea>0) h->Scale(1000000./eArea);
80 return h;
81}

◆ Wmt()

float EdbSEQ::Wmt ( EdbSegP s)
inline
60{ return s.W()*50. + s.DZ(); }

Member Data Documentation

◆ eArea

Double_t EdbSEQ::eArea

effective area of the pattern to be equalized

◆ eChiLimits

TVector2* EdbSEQ::eChiLimits
private

[min,max]

◆ eExcludeThetaRange

TObjArray EdbSEQ::eExcludeThetaRange

can be added EdbSegP with tx,ty, sigmaTX,sigmaTY to be excluded

◆ eHEq

EdbH1 EdbSEQ::eHEq

◆ eNP

int EdbSEQ::eNP

number of points for the functions calculation

◆ eNsigma

Double_t EdbSEQ::eNsigma

=4

◆ eS0mt

Double_t EdbSEQ::eS0mt

= 270.*340.; // area unit for Nseg calculation

◆ eThetaLimits

TVector2* EdbSEQ::eThetaLimits
private

[min,max]

◆ eWLimits

TVector2* EdbSEQ::eWLimits
private

[min,max]

◆ eXLimits

TVector2* EdbSEQ::eXLimits
private

[min,max] area limits for the preselection procedure

◆ eYLimits

TVector2* EdbSEQ::eYLimits
private

[min,max]


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