FEDRA emulsion software from the OPERA Collaboration
EdbGA Class Reference

grains analysys More...

#include <EdbGA.h>

Inheritance diagram for EdbGA:
Collaboration diagram for EdbGA:

Public Member Functions

void CheckViewGrains (const char *options="")
 
void CheckViewGrains (int vid, const char *options="")
 
 EdbGA ()
 
 EdbGA (char *fname, float bx, float by, float bz)
 
void GetClustPFile (const char *file)
 
TTree * GetTree (void)
 
void GrainStat (TClonesArray *clusters, float &x0, float &y0, float &z0)
 
int GrainStat2 (TClonesArray *clusters, float &x0, float &y0, float &z0, float &vol, float &amin, float &amax, float &zmin, float &zmax, int &fmin, int &fmax)
 
void InitTree (const char *file="grain_chains.root")
 
int MakeGrainsTree (TClonesArray *clust, TIndexCell &chains, const char *options="")
 
void SelectGrains (const char *selection, const char *outfile="grains_chains_selection.root")
 
void SelectGrains (TCut c1, const char *outfile="grains_chains_selection.root")
 
void SetBin (float bx, float by, float bz)
 
void SetRun (char *fname)
 
void VerticalChains (TClonesArray *clusters, TIndexCell &chains)
 
void VerticalChainsA (TClonesArray *clusters)
 
virtual ~EdbGA ()
 

Private Attributes

Float_t eBinX
 
Float_t eBinY
 
Float_t eBinZ
 
TClonesArray * eClusters
 
TTree * eGrains
 
TFile * eGrainsFile
 
EdbRuneRun
 
Int_t eVid
 view under prcessing More...
 

Detailed Description

grains analysys

Constructor & Destructor Documentation

◆ EdbGA() [1/2]

EdbGA::EdbGA ( )
inline
34{eRun=0;}
EdbRun * eRun
Definition: EdbGA.h:21

◆ EdbGA() [2/2]

EdbGA::EdbGA ( char *  fname,
float  bx,
float  by,
float  bz 
)

◆ ~EdbGA()

EdbGA::~EdbGA ( )
virtual

39{
40 if(eGrainsFile) {
41 eGrainsFile->Write();
42 eGrainsFile->Purge();
43 eGrainsFile->Close();
44 }
45 if(eRun) delete eRun;
46}
TFile * eGrainsFile
Definition: EdbGA.h:29

Member Function Documentation

◆ CheckViewGrains() [1/2]

void EdbGA::CheckViewGrains ( const char *  options = "")

94{
95 int nview = eRun->GetEntries();
96 for(int i=0; i< nview ; i++) CheckViewGrains(i,options);
97}
void CheckViewGrains(const char *options="")
Definition: EdbGA.cxx:93
int GetEntries() const
Definition: EdbRun.h:136

◆ CheckViewGrains() [2/2]

void EdbGA::CheckViewGrains ( int  vid,
const char *  options = "" 
)

101{
102 eVid = vid;
103 TClonesArray *clusters=0;
104
105 clusters = eRun->GetEntryClusters(vid);
106
107 printf("%d clusters are read for view %d \n", clusters->GetEntriesFast(), vid );
108
109 TIndexCell chains;
110 VerticalChains(clusters,chains);
111 MakeGrainsTree(clusters,chains,options);
112}
Int_t eVid
view under prcessing
Definition: EdbGA.h:27
void VerticalChains(TClonesArray *clusters, TIndexCell &chains)
Definition: EdbGA.cxx:151
int MakeGrainsTree(TClonesArray *clust, TIndexCell &chains, const char *options="")
Definition: EdbGA.cxx:201
TClonesArray * GetEntryClusters(int entry) const
Definition: EdbRun.cxx:507
sort collection with attributes
Definition: TIndexCell.h:19

◆ GetClustPFile()

void EdbGA::GetClustPFile ( const char *  file)

50{
51 TFile *fcl = new TFile(file);
52 TTree *clust = (TTree *)fcl->Get("clust");
53 EdbClustP *clp = 0;
54 clust->SetBranchAddress("cl",&clp);
55
56 int nentr = (int)(clust->GetEntries());
57 TClonesArray *clarr = new TClonesArray("EdbCluster");
58
59 int ncl =0;
60 for(int i=0; i<nentr; i++) {
61 clust->GetEntry(i);
62 //if( clp->GetZ()>-50.) continue;
63 //if( clp->Xp()<500.) continue;
64 //if( clp->Xp()>600.) continue;
65 //if( clp->Yp()<500.) continue;
66 //if( clp->Yp()>600.) continue;
67
68 new((*clarr)[ncl])
69 EdbCluster( clp->Xcg(),clp->Ycg(),clp->Z(),
70 clp->GetArea(), clp->Peak(), clp->GetFrame(), 0);
71 ncl++;
72 }
73 fcl->Close();
74 delete fcl;
75 printf("%d clusters are readed from %s\n",ncl,file);
76
77
78 SetBin(4,4,1);
79 TIndexCell chains;
80 printf("chains\n");
81 VerticalChains(clarr,chains);
82// for(int i=0; i<10; i++) {
83// printf("chainsA\n");
84// VerticalChainsA(clarr);
85// }
86 printf("\n tree: \n\n");
87 InitTree("grains_tree.root");
88 MakeGrainsTree(clarr,chains);
89 eGrainsFile->Close();
90}
cluster reconstruction
Definition: EdbIP.h:58
float Xcg() const
Definition: EdbIP.h:71
float Ycg() const
Definition: EdbIP.h:72
float Peak() const
Definition: EdbIP.h:73
Definition: EdbCluster.h:19
Int_t GetFrame() const
Definition: EdbCluster.h:56
virtual Float_t Z() const
Definition: EdbCluster.h:80
Float_t GetArea() const
Definition: EdbCluster.h:54
void SetBin(float bx, float by, float bz)
Definition: EdbGA.h:41
void InitTree(const char *file="grain_chains.root")
Definition: EdbGA.cxx:176
TFile * file
Definition: write_pvr.C:3

◆ GetTree()

TTree * EdbGA::GetTree ( void  )
inline
57{return eGrains;}
TTree * eGrains
Definition: EdbGA.h:30

◆ GrainStat()

void EdbGA::GrainStat ( TClonesArray *  clusters,
float &  x0,
float &  y0,
float &  z0 
)

291 {
292 double xs=0,ys=0,zs=0;
293 EdbCluster *cl=0;
294 int ncl = clusters->GetEntriesFast();
295 for( int i=0; i<ncl; i++ ) {
296 cl = (EdbCluster *)clusters->At(i);
297 xs += cl->X();
298 ys += cl->Y();
299 zs += cl->Z();
300 }
301 x0 = xs/ncl;
302 y0 = ys/ncl;
303 z0 = zs/ncl;
304 }
brick z0
Definition: RecDispMC.C:106
virtual Float_t Y() const
Definition: EdbCluster.h:79
virtual Float_t X() const
Definition: EdbCluster.h:78

◆ GrainStat2()

int EdbGA::GrainStat2 ( TClonesArray *  clusters,
float &  x0,
float &  y0,
float &  z0,
float &  vol,
float &  amin,
float &  amax,
float &  zmin,
float &  zmax,
int &  fmin,
int &  fmax 
)

309 {
310 float x,y,z,area,xs=0,ys=0,zs=0,areasum=0;
311 int frame;
312 EdbCluster *cl=0;
313 amin = 10000000;
314 amax = 0;
315 zmin = 1.e30;
316 zmax = -1.e30;
317 fmin = 10000000;
318 fmax = 0;
319
320 int ncl = clusters->GetEntriesFast();
321 for( int i=0; i<ncl; i++ )
322 {
323 cl = (EdbCluster *) clusters->At(i);
324 x = cl->X();
325 y = cl->Y();
326 z = cl->Z();
327 frame = cl->GetFrame() ;
328 area = cl->GetArea() ;
329 if(area < amin) amin = area ;
330 if(area > amax) amax = area ;
331 if(frame < fmin ) { fmin = frame ; zmax = z ;}
332 if(frame > fmax ) { fmax = frame ; zmin = z ;}
333 /* xs += x*area;
334 ys += y*area;
335 zs += z*area;
336 */ xs += x;
337 ys += y;
338 zs += z;
339 areasum += area;
340 }
341 /* x0 = xs/areasum;
342 y0 = ys/areasum;
343 z0 = zs/areasum;
344 */ x0 = xs/ncl;
345 y0 = ys/ncl;
346 z0 = zs/ncl;
347
348 float zlen = zmax - zmin;
349 if(zlen && ncl>1) vol = areasum*zlen/(ncl-1) ; // should be (areasum/ncl) * (zlen*ncl/(ncl-1))
350 else vol = -999;
351
352 return 0;
353}
EdbPatternsVolume * vol
Definition: RecDispNU.C:116

◆ InitTree()

void EdbGA::InitTree ( const char *  file = "grain_chains.root")

177{
178 eGrainsFile= new TFile(file,"RECREATE");
179 eGrains= new TTree("grains","grains");
180
181 Int_t ncl=0,fmin,fmax;
182 eClusters=new TClonesArray("EdbCluster");
183 Float_t x0,y0,z0,vol,amin,amax, zmin, zmax;
184
185 eGrains->Branch("vid",&eVid,"eVid/I");
186 eGrains->Branch("ncl",&ncl,"ncl/I");
187 eGrains->Branch("clusters",&eClusters);
188 eGrains->Branch("x0",&x0,"x0/F");
189 eGrains->Branch("y0",&y0,"y0/F");
190 eGrains->Branch("z0",&z0,"z0/F");
191 eGrains->Branch("vol",&vol,"vol/F");
192 eGrains->Branch("amin",&amin,"amin/F");
193 eGrains->Branch("amax",&amax,"amax/F");
194 eGrains->Branch("zmin",&zmin,"zmin/F");
195 eGrains->Branch("zmax",&zmax,"zmax/F");
196 eGrains->Branch("fmin",&fmin,"fmin/I");
197 eGrains->Branch("fmax",&fmax,"fmax/I");
198}
TClonesArray * eClusters
Definition: EdbGA.h:31

◆ MakeGrainsTree()

int EdbGA::MakeGrainsTree ( TClonesArray *  clust,
TIndexCell chains,
const char *  option = "" 
)

process the options

202{
203
205 int min_ncl=1; // minimum number of clusters of each grains
206 if (! sscanf(option,"MIN_NCL") ) sscanf(option,"MIN_NCL=%d",&min_ncl);
207
208 // set branch addresses
209 Int_t ncl=0,fmin, fmax;
210 Float_t x0,y0,z0,vol,amin,amax, zmin, zmax;
211 eGrains->SetBranchAddress("vid",&eVid);
212 eGrains->SetBranchAddress("ncl",&ncl);
213 eGrains->SetBranchAddress("clusters",&eClusters);
214 eGrains->SetBranchAddress("x0",&x0);
215 eGrains->SetBranchAddress("y0",&y0);
216 eGrains->SetBranchAddress("z0",&z0);
217 eGrains->SetBranchAddress("vol",&vol);
218 eGrains->SetBranchAddress("amin",&amin);
219 eGrains->SetBranchAddress("amax",&amax);
220 eGrains->SetBranchAddress("zmin",&zmin);
221 eGrains->SetBranchAddress("zmax",&zmax);
222 eGrains->SetBranchAddress("fmin",&fmin);
223 eGrains->SetBranchAddress("fmax",&fmax);
224
225 int ngr=0;
226 int entry=0;
227
228 TIndexCell *cx;
229 TIndexCell *cy;
230 TIndexCell *cz;
231
232 float zlvl;
233 EdbCluster *cl;
234
235 int nix,niy,niz,nie;
236
237 nix = chains.GetEntries();
238 for(int ix=0; ix<nix; ix++) {
239 cx = chains.At(ix);
240 niy = cx->GetEntries();
241 for(int iy=0; iy<niy; iy++) {
242 cy = cx->At(iy);
243 zlvl = -9999999.;
244 niz = cy->GetEntries();
245 for(int iz=0; iz<niz; iz++) {
246 cz = cy->At(iz);
247 nie = cz->GetEntries();
248 for(int ie=0; ie<nie; ie++) {
249 entry = cz->At(ie)->Value();
250 cl = (EdbCluster *)clust->UncheckedAt(entry);
251
252 if( cl->Z()-zlvl > 2.*eBinZ ) {
253 if( ncl>0 ) {
254 if(ncl>=min_ncl) {
255 //GrainStat(eClusters,x0,y0,z0);
256 GrainStat2(eClusters,x0,y0,z0,vol,amin,amax, zmin, zmax, fmin, fmax);
257 eGrains->Fill();
258 ngr++;
259 }
260 eClusters->Clear();
261 ncl=0;
262 }
263 }
264 new((*eClusters)[ncl]) EdbCluster( *cl );
265 ncl++;
266 zlvl=cl->Z();
267 }
268 }
269
270 if( ncl>0 ) {
271 if(ncl>=min_ncl) {
272 //GrainStat(eClusters,x0,y0,z0);
273 GrainStat2(eClusters,x0,y0,z0,vol,amin,amax, zmin, zmax, fmin, fmax);
274 eGrains->Fill();
275 ngr++;
276 }
277 eClusters->Clear();
278 ncl=0;
279 }
280 }
281 }
282
283 eGrains->Write();
284 printf("%d grains are selected\n",ngr);
285
286 return ngr;
287 }
TLegendEntry * entry
Definition: Canv_SYSTEMATICS_ALLCOMBINED__RMSEnergy__vs__Energy__ELECTRON.C:130
Float_t eBinZ
Definition: EdbGA.h:25
int GrainStat2(TClonesArray *clusters, float &x0, float &y0, float &z0, float &vol, float &amin, float &amax, float &zmin, float &zmax, int &fmin, int &fmax)
Definition: EdbGA.cxx:306
Long_t Value() const
Definition: TIndexCell.h:79
Int_t GetEntries() const
Definition: TIndexCell.h:82
TIndexCell const * At(Int_t narg, Int_t vind[]) const
Definition: TIndexCell.cpp:519

◆ SelectGrains() [1/2]

void EdbGA::SelectGrains ( const char *  selection,
const char *  outfile = "grains_chains_selection.root" 
)

save only the selected grains

361{
363
364 // set branch addresses
365 Int_t ncl,fmin, fmax;
366 Float_t x0,y0,z0,vol,amin,amax, zmin, zmax;
367 eGrains->SetBranchAddress("vid",&eVid);
368 eGrains->SetBranchAddress("ncl",&ncl);
369 eGrains->SetBranchAddress("clusters",&eClusters);
370 eGrains->SetBranchAddress("x0",&x0);
371 eGrains->SetBranchAddress("y0",&y0);
372 eGrains->SetBranchAddress("z0",&z0);
373 eGrains->SetBranchAddress("vol",&vol);
374 eGrains->SetBranchAddress("amin",&amin);
375 eGrains->SetBranchAddress("amax",&amax);
376 eGrains->SetBranchAddress("zmin",&zmin);
377 eGrains->SetBranchAddress("zmax",&zmax);
378 eGrains->SetBranchAddress("fmin",&fmin);
379 eGrains->SetBranchAddress("fmax",&fmax);
380
381 //select the event list
382 printf("Total: %d grains\n", (int)eGrains->GetEntries());
383 eGrains->Draw(">>lst",selection);
384 TEventList* lst = (TEventList*) gDirectory->Get("lst");
385
386 // create a new tree
387 TFile* nfile = new TFile(outfile,"recreate");
388 TTree* tree2 = (TTree*) eGrains->CloneTree(0);
389 tree2->CopyAddresses(eGrains);
390
391 int nentries = lst->GetN(); printf("Select %d grains",nentries); printf(" with \"%s\"\n", selection);
392 printf("Create the new tree... ");
393 for(int i=0;i<nentries;i++)
394 {
395 eGrains->GetEntry(lst->GetEntry(i));
396 tree2->Fill();
397 if((i%5000)==0) printf("\b\b\b\b%3d%%",(int)((double)i/double(nentries)*100.)) ;
398 }
399 printf("\b\b\b\b100%%\n");
400 tree2->Write("grains",TObject::kOverwrite);
401 nfile->Close();
402}
int nentries
Definition: check_shower.C:40

◆ SelectGrains() [2/2]

void EdbGA::SelectGrains ( TCut  c1,
const char *  outfile = "grains_chains_selection.root" 
)

356{
357 SelectGrains(c1.GetTitle(), outfile);
358}
void SelectGrains(TCut c1, const char *outfile="grains_chains_selection.root")
Definition: EdbGA.cxx:355
TCanvas * c1
Definition: energy.C:13

◆ SetBin()

void EdbGA::SetBin ( float  bx,
float  by,
float  bz 
)
inline
42 { eBinX=bx; eBinY=by; eBinZ=bz; }
Float_t eBinX
Definition: EdbGA.h:23
Float_t eBinY
Definition: EdbGA.h:24

◆ SetRun()

void EdbGA::SetRun ( char *  fname)
inline
39 { if(eRun) delete eRun;
40 eRun=new EdbRun(fname); }
Definition: EdbRun.h:75
const char * fname
Definition: mc2raw.cxx:41

◆ VerticalChains()

void EdbGA::VerticalChains ( TClonesArray *  clusters,
TIndexCell chains 
)

chains: "x:y:z:entry"

152{
154
155 Long_t v[5];
156 float xb=eBinX, yb=eBinY, zb=eBinZ; // bins
157
158 EdbCluster *cl;
159
160 int ncl = clusters->GetEntriesFast();
161 for( int i=0; i<ncl; i++ ) {
162 cl = (EdbCluster*)(clusters->UncheckedAt(i));
163 v[0] = (Long_t)(cl->X()/xb);
164 v[1] = (Long_t)(cl->Y()/yb);
165 v[2] = (Long_t)(cl->Z()/zb);
166 v[3] = (Long_t)i;
167 chains.Add(4,v);
168 }
169
170 chains.Sort();
171 //chains.PrintStat();
172 //chains.PrintPopulation(3);
173}
void Sort(Int_t upto=kMaxInt)
Definition: TIndexCell.cpp:539
Int_t Add(Int_t narg, Long_t varg[])
Definition: TIndexCell.cpp:602

◆ VerticalChainsA()

void EdbGA::VerticalChainsA ( TClonesArray *  clusters)

116{
117 float eDX = 1301, eDY=1025;
118 eBinX=4; eBinY=4;
119 const int nx = (int)(eDX/eBinX);
120 const int ny = (int)(eDY/eBinY);
121 const int nz = 50;
122
123 int *v = new int[nx*nx*ny];
124
125 int ix,iy,iz,iv;
126
127 EdbCluster *cl;
128
129 for(iz=0; iz<nz; iz++)
130 for(iy=0; iy<ny; iy++)
131 for(ix=0; ix<nx; ix++) {
132 iv = ix + nx*iy + nx*ny*iz;
133 v[iv] = -1;
134 }
135
136 int ncl = clusters->GetEntriesFast();
137 for( int i=0; i<ncl; i++ ) {
138 cl = (EdbCluster*)(clusters->UncheckedAt(i));
139 ix = (int)(cl->X()/eBinX);
140 iy = (int)(cl->Y()/eBinY);
141 iz = (int)(cl->GetFrame());
142
143 iv = ix + nx*iy + nx*ny*iz;
144 v[iv] = i;
145 }
146
147 delete[] v;
148}

Member Data Documentation

◆ eBinX

Float_t EdbGA::eBinX
private

◆ eBinY

Float_t EdbGA::eBinY
private

◆ eBinZ

Float_t EdbGA::eBinZ
private

◆ eClusters

TClonesArray* EdbGA::eClusters
private

◆ eGrains

TTree* EdbGA::eGrains
private

◆ eGrainsFile

TFile* EdbGA::eGrainsFile
private

◆ eRun

EdbRun* EdbGA::eRun
private

◆ eVid

Int_t EdbGA::eVid
private

view under prcessing


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