FEDRA emulsion software from the OPERA Collaboration
EdbDataStore Class Reference

#include <EdbDataStore.h>

Inheritance diagram for EdbDataStore:
Collaboration diagram for EdbDataStore:

Public Member Functions

void AddPattern (EdbPattern *pat)
 
EdbSegPAddSegment (EdbSegP *seg, EdbSegmentCut *cut=0, int Plt0=0, int Plt1=100)
 add methods More...
 
void AddTrack (EdbTrackP *tr)
 
void AddVertex (EdbVertex *v)
 
void Clear (bool hard=false)
 clear methods More...
 
void ClearGeom ()
 
void ClearRaw (bool hard=false)
 
void ClearSeg (bool hard=false)
 
void ClearTracks (bool hard=false)
 
void ClearVTX ()
 
void DoEfficiency (TF1 *eff_seg, TF1 *eff_mtk)
 
void DoSmearing (EdbScanCond *cond_btk, EdbScanCond *cond_mtk=0)
 methods for simulation: More...
 
 EdbDataStore ()
 
EdbLayerFindLayer (int plate, int side=0)
 
EdbTrackPFindLongTrk (int nsmin=8)
 
EdbPatternFindPattern (int plate, int side=0)
 
EdbVertexFindPrimVtx ()
 
EdbTrackPFindTrack (int id)
 find methods More...
 
EdbVertexFindVertex (int id)
 
long Gen_mtk_BG (long NBG, int Plate, int Side, TH2 *pdf_Ang, TH2 *pdf_WT=0)
 
EdbPatternGetPattern (int n, bool btk=true)
 
EdbPatternsVolumeGetPV (bool btk=true)
 
EdbPatternGetRawPat (int n)
 
EdbPatternGetSegPat (int n)
 
EdbTrackPGetTrack (int n)
 get methods More...
 
EdbVertexGetVertex (int n)
 
void LoadMCVertices (TObjArray *vtx)
 restore MC info methods More...
 
void MakePattern (double z, int plate, int side)
 
int Nplt ()
 
int Nt ()
 count methods: More...
 
int Nv ()
 
void PrintBrief ()
 print methods More...
 
void PrintPatterns ()
 
void PrintTracks (int vlev=0)
 
void Restore_PatFromGeom (int np0=0, int np1=1000)
 
void Restore_PIDFromID ()
 
void Restore_SegFromTrx (EdbSegmentCut *cut=0, int Plt0=0, int Plt1=1000)
 
void Restore_TrxFromVtx ()
 
void SavePlateToRaw (const char *fname, int PID, Option_t *option="RECREATE")
 
void SaveToRaw (const char *dir="./", const EdbID &idset="0.0.0.0", Option_t *option="RECREATE", bool doaff=true)
 save methods: More...
 
void SetOwnTracks (bool own=true)
 setown methods: More...
 
void SetOwnTrkSegs ()
 
void SetOwnVertices (bool own=true)
 
void TransferGeometry (EdbDataStore *ds)
 
void TransferTo (EdbDataStore *ds, char level, EdbSegmentCut *cut=0, int FromPlate=0, int ToPlate=57)
 transfer methods More...
 
 ~EdbDataStore ()
 

Static Public Member Functions

static void TransferSegs (EdbPatternsVolume *pv0, EdbPatternsVolume *pv1, EdbSegmentCut *cut=0, int FromPlate=0, int ToPlate=57)
 

Public Attributes

EdbBrickPeBrick
 
EdbPatternsVolume eRawPV
 geometry More...
 
EdbPatternsVolume eSegPV
 
TObjArray eTracks
 
TObjArray eVTX
 

Constructor & Destructor Documentation

◆ EdbDataStore()

EdbDataStore::EdbDataStore ( )

11 {
12 eRawPV.Set0();
13 eSegPV.Set0();
14};
EdbPatternsVolume eRawPV
geometry
Definition: EdbDataStore.h:82
EdbPatternsVolume eSegPV
Definition: EdbDataStore.h:83
void Set0()
Definition: EdbPattern.cxx:1585

◆ ~EdbDataStore()

EdbDataStore::~EdbDataStore ( )

16{ Clear();};
void Clear(bool hard=false)
clear methods
Definition: EdbDataStore.h:30

Member Function Documentation

◆ AddPattern()

void EdbDataStore::AddPattern ( EdbPattern pat)

257 {
258 assert(pat);
259 assert(FindLayer(pat->Plate(),pat->Side()));
260 if(pat->Side()==0)eSegPV.AddPattern(pat);
261 else eRawPV.AddPattern(pat);
262}
EdbLayer * FindLayer(int plate, int side=0)
Definition: EdbDataStore.cxx:206
Int_t Side() const
Definition: EdbPattern.h:328
Int_t Plate() const
Definition: EdbPattern.h:327
void AddPattern(EdbPattern *pat)
Definition: EdbPattern.cxx:1707

◆ AddSegment()

EdbSegP * EdbDataStore::AddSegment ( EdbSegP s,
EdbSegmentCut cut = 0,
int  Plt0 = 0,
int  Plt1 = 100 
)

add methods


126 {
127 if(s->Plate()<Plt0 || s->Plate()>Plt1)return 0;
128 if(cut){
129 float par[]={s->X(),s->Y(),s->TX(),s->TY(),s->W()};
130 if(!cut->PassCut(par))return 0;
131 }
132 if(gEDBDEBUGLEVEL>5){
133 s->PrintNice();
134 printf("plt#%d side#%d\n",s->Plate(),s->Side());
135 }
136 EdbPattern* p=FindPattern(s->Plate(),s->Side());
137 assert(p);
138 EdbSegP* seg=p->AddSegment(*s);
139 //p->SetID(s->Plate());
140 seg->SetPID(p->ID());
141 return seg;
142}
EdbPattern * FindPattern(int plate, int side=0)
Definition: EdbDataStore.cxx:220
Definition: EdbPattern.h:273
Definition: EdbSegP.h:21
void SetPID(int pid)
Definition: EdbSegP.h:129
TCut cut
Definition: check_shower.C:6
s
Definition: check_shower.C:55
gEDBDEBUGLEVEL
Definition: energy.C:7
p
Definition: testBGReduction_AllMethods.C:8

◆ AddTrack()

void EdbDataStore::AddTrack ( EdbTrackP tr)
inline
62{assert(tr!=0); eTracks.Add(tr);}
TTree * tr
Definition: Shower_E_FromShowerRoot.C:5
TObjArray eTracks
Definition: EdbDataStore.h:84

◆ AddVertex()

void EdbDataStore::AddVertex ( EdbVertex v)
inline
63{assert(v!=0); eVTX.Add(v);}
TObjArray eVTX
Definition: EdbDataStore.h:85

◆ Clear()

void EdbDataStore::Clear ( bool  hard = false)
inline

clear methods

30{ClearTracks(hard);ClearVTX();ClearRaw(hard);ClearSeg(hard);}
void ClearRaw(bool hard=false)
Definition: EdbDataStore.cxx:319
void ClearSeg(bool hard=false)
Definition: EdbDataStore.cxx:293
void ClearVTX()
Definition: EdbDataStore.cxx:285
void ClearTracks(bool hard=false)
Definition: EdbDataStore.cxx:289

◆ ClearGeom()

void EdbDataStore::ClearGeom ( )

276 {
277 if(eBrick==0)return;
278 for(int np=0;np<eBrick->Npl();++np){
279 eBrick->GetPlate(np)->Clear();
280 delete eBrick->GetPlate(np);
281 }
282 eBrick->Clear();
283}
void Clear()
Definition: EdbBrick.h:54
EdbPlateP * GetPlate(int i)
Definition: EdbBrick.h:53
int Npl() const
Definition: EdbBrick.h:51
EdbBrickP * eBrick
Definition: EdbDataStore.h:81

◆ ClearRaw()

void EdbDataStore::ClearRaw ( bool  hard = false)

319 {
321
322
323 int N=pv->Npatterns();
324 if(!N)return;
325 TObjArray* Sgs=0;
326 for(int k=0;k<N;k++) {
327 Sgs=pv->GetPattern(k)->GetSegments();
328 if(Sgs)Sgs->Delete();
329// delete pv->GetPattern(k);
330 }
331 if(hard){
332 pv->ePatterns->SetOwner(1);
333 pv->ePatterns->Clear();
334 pv->Clear();
335 pv->DropCell();
336 pv->DropCouples();
337 }
338}
Definition: EdbPattern.h:334
Int_t Npatterns() const
Definition: EdbPattern.h:366
int DropCouples()
Definition: EdbPattern.cxx:1592
void DropCell()
Definition: EdbPattern.cxx:1760
EdbPattern * GetPattern(int id) const
Definition: EdbPattern.cxx:1721
TObjArray * ePatterns
collection of patterns
Definition: EdbPattern.h:342
TClonesArray * GetSegments() const
Definition: EdbPattern.h:69

◆ ClearSeg()

void EdbDataStore::ClearSeg ( bool  hard = false)

293 {
295// int N=pv->Npatterns();
296// if(!N)return;
297// TObjArray* Sgs=0;
298// for(int k=0;k<N;k++) {
299// Sgs=pv->GetPattern(k)->GetSegments();
300// if(Sgs)Sgs->Delete();
301// }
302 int N=pv->Npatterns();
303 if(!N)return;
304 TObjArray* Sgs=0;
305 for(int k=0;k<N;k++) {
306 Sgs=pv->GetPattern(k)->GetSegments();
307 if(Sgs)Sgs->Delete();
308// delete pv->GetPattern(k);
309 }
310 if(hard){
311 pv->ePatterns->SetOwner(1);
312 pv->DropCell();
313 pv->DropCouples();
314 pv->ePatterns->Clear();
315 pv->Clear();
316 }
317}

◆ ClearTracks()

void EdbDataStore::ClearTracks ( bool  hard = false)

289 {
290 eTracks.Clear();
291}

◆ ClearVTX()

void EdbDataStore::ClearVTX ( )

285 {
286 eVTX.Clear();
287}

◆ DoEfficiency()

void EdbDataStore::DoEfficiency ( TF1 *  eff_seg,
TF1 *  eff_mtk 
)

386 {
387 EdbSegP* s=0;
388 double eff, ran, th;
389 if(eff_seg){
390 for(int np=0;np<eSegPV.Npatterns();++np){
391 EdbPattern* pat=GetSegPat(np);
392 for(int ns=0; ns< pat->N();++ns){
393 s=pat->GetSegment(ns);
394 th=s->Theta();
395 eff=(th<=eff_seg->GetXmax())?eff_seg->Eval(th):0;
396 ran=gRandom->Rndm();
397 if(ran>eff)s->SetW(-s->W());
398 }
399 }
400 }
401 if(eff_mtk){
402 for(int np=0;np<eRawPV.Npatterns();++np){
403 int n1=0;
404 EdbPattern* pat=GetRawPat(np);
405 for(int ns=0; ns< pat->N();++ns){
406 s=pat->GetSegment(ns);
407 th=s->Theta();
408 Log(5,"DoEfficiency",Form("mtk#%d/%d theta=%2.4f [%2.4f %2.4f]\n",ns,pat->N(),th,s->TX(),s->TY()));
409 eff=(th<=eff_mtk->GetXmax())?eff_mtk->Eval(th):0;
410 ran=gRandom->Rndm();
411 Log(5,"DoEfficiency",Form("eff=%2.4f, ran=%2.4f => %s\n",eff,ran,eff>ran?"survive":"KILL"));
412 if(ran>eff)s->SetW(-s->W());
413 else ++n1;
414 }
415 Log(2,"DoEfficiency",Form("Plate #%d side#%d - efficiency %d/%d\n",pat->Plate(),pat->Side(),n1,pat->N()));
416 }
417 }
418}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
EdbPattern * GetRawPat(int n)
Definition: EdbDataStore.h:48
EdbPattern * GetSegPat(int n)
Definition: EdbDataStore.h:47
Int_t N() const
Definition: EdbPattern.h:86
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:66

◆ DoSmearing()

void EdbDataStore::DoSmearing ( EdbScanCond cond_btk,
EdbScanCond cond_mtk = 0 
)

methods for simulation:


372 {
373 EdbPVGen pvg;
374 if(cond_btk){
375 pvg.SetVolume(&eSegPV);
376 pvg.SetScanCond(cond_btk);
377 pvg.SmearSegments();
378 }
379 if(cond_mtk){
380 pvg.SetVolume(&eRawPV);
381 pvg.SetScanCond(cond_mtk);
382 pvg.SmearSegments();
383 }
384}
Definition: EdbPVGen.h:18
void SmearSegments()
Definition: EdbPVGen.cxx:190
void SetVolume(EdbPatternsVolume *pv)
Definition: EdbPVGen.h:34
void SetScanCond(EdbScanCond *scan)
Definition: EdbPVGen.h:35

◆ FindLayer()

EdbLayer * EdbDataStore::FindLayer ( int  plate,
int  side = 0 
)

206 {
207 assert(plate>=0);
208 assert(side>=0 && side<3);
209
210 EdbPlateP* plt=0;
211 for(int np=0;np<eBrick->Npl();++np){
212 plt=eBrick->GetPlate(np);
213 if(plt->ID()==plate){
214 return plt->GetLayer(side);
215 }
216 }
217 return 0;
218}
int ID() const
Definition: EdbLayer.h:73
Definition: EdbBrick.h:14
EdbLayer * GetLayer(int i)
Definition: EdbBrick.h:29
Int_t plate
Definition: merge_Energy_SytematicSources_Electron.C:1

◆ FindLongTrk()

EdbTrackP * EdbDataStore::FindLongTrk ( int  nsmin = 8)

find longest track

236 {
238 int nsmax=0;
239 int imax=0;
240 EdbTrackP* trk=0;
241 for(int nt=0;nt<Nt();++nt){
242 trk=GetTrack(nt);
243 if(trk->PDG()!=22 && trk->P()>1.&& trk->N()>nsmax){
244 imax=nt;
245 nsmax=trk->N();
246 }
247 }
248 if(nsmax<nsmin)return 0;
249 return GetTrack(imax);
250}
int Nt()
count methods:
Definition: EdbDataStore.h:37
EdbTrackP * GetTrack(int n)
get methods
Definition: EdbDataStore.h:45
Float_t P() const
Definition: EdbSegP.h:152
Definition: EdbPattern.h:113
Int_t PDG() const
Definition: EdbPattern.h:150
Int_t N() const
Definition: EdbPattern.h:177

◆ FindPattern()

EdbPattern * EdbDataStore::FindPattern ( int  plate,
int  side = 0 
)

220 {
221 Log(3,"EdbDataStore::FindPattern",Form("search plt#%d side=%d",plate,side));
222 assert(plate>=0);
223 assert(side>=0 && side<3);
224
225 EdbPatternsVolume* pv=0;
226 pv=(side==0)?&eSegPV:&eRawPV;
227 EdbPattern* pat=0;
228 for(int np=0;np<pv->Npatterns();++np){
229 pat=pv->GetPattern(np);
230 //Log(6,"pat",Form("#%d/%d plt#%d side#%d",np,pv->Npatterns(),pat->Plate(),pat->Side()));
231 if(pat->Plate()==plate && pat->Side()==side)return pat;
232 }
233 return 0;
234}

◆ FindPrimVtx()

EdbVertex * EdbDataStore::FindPrimVtx ( )

252 {
253 EdbTrackP* trk=FindLongTrk(8);
254 return trk->VertexS();
255}
EdbTrackP * FindLongTrk(int nsmin=8)
Definition: EdbDataStore.cxx:236
EdbVertex * VertexS()
Definition: EdbPattern.cxx:1191

◆ FindTrack()

EdbTrackP * EdbDataStore::FindTrack ( int  id)

find methods


187 {
188 EdbTrackP* trk=0;
189 for(int nt=0;nt<eTracks.GetEntries();++nt){
190 trk=GetTrack(nt);
191// printf("id%d/%d=%d\n",nt,eTracks.GetEntries(),trk->ID());
192 if(trk->ID()==id)return trk;
193 }
194 return 0;
195}
Int_t ID() const
Definition: EdbSegP.h:147

◆ FindVertex()

EdbVertex * EdbDataStore::FindVertex ( int  id)

197 {
198 EdbVertex* vtx=0;
199 for(int nv=0;nv<eVTX.GetEntries();++nv){
200 vtx=GetVertex(nv);
201 if(vtx->ID()==id)return vtx;
202 }
203 return 0;
204}
EdbVertex * GetVertex(int n)
Definition: EdbDataStore.h:46
Definition: EdbVertex.h:69
Int_t ID() const
Definition: EdbVertex.h:126

◆ Gen_mtk_BG()

long EdbDataStore::Gen_mtk_BG ( long  NBG,
int  Plate,
int  Side,
TH2 *  pdf_Ang,
TH2 *  pdf_WT = 0 
)

Generate background microtracks, using distributions: 1)from pdf_Ang(TY:TX) - should be normalized so that its maximum bin content = 1 2)from pdb_WT (W:tan(Theta)) - each row with fixed X(i.e. slope) should have maximum content = 1

Side==0 - this is basetrack. So we need to project each segment to RAW patterns

420 {
424 EdbPattern* pat=FindPattern(Plate,Side);
425 EdbPattern *p1=0,*p2=0;
426 if(Side==0){
428 p1=FindPattern(Plate,1);
429 p2=FindPattern(Plate,2);
430 }
431 assert(pat!=0);
432 EdbLayer * lay=FindLayer (Plate,Side);
433 assert(lay!=0);
434 float x,y,z,tx,ty,w=16,p;
435 z=lay->Z();
436 int Ntot=0;
437 EdbSegP *seg=0,*s1=0;
438 for(long n=0; n<NBG; ++n){
439 tx=gRandom->Uniform(pdf_Ang->GetXaxis()->GetXmin(),pdf_Ang->GetXaxis()->GetXmax());
440 ty=gRandom->Uniform(pdf_Ang->GetYaxis()->GetXmin(),pdf_Ang->GetYaxis()->GetXmax());
441 p=pdf_Ang->Interpolate(tx,ty);
442 Log(5,"Gen_mtk_BG",Form("tt=[%4.2f %4.2f] p=%6.4f\n",tx,ty,p));
443 if(gRandom->Uniform()>p){n--; continue;}
444 Ntot++;
445 x=gRandom->Uniform(lay->Xmin(),lay->Xmax());
446 y=gRandom->Uniform(lay->Ymin(),lay->Ymax());
447 seg=pat->AddSegment(-1,x,y,tx,ty,w,0);
448 if(p1){
449 s1=p1->AddSegment(*seg);
450 s1->PropagateTo(p1->Z());
451 }
452 if(p2){
453 s1=p2->AddSegment(*seg);
454 s1->PropagateTo(p2->Z());
455 }
456 }
457 Log(0,"Gen_mtk_BG",Form("Generated %ld/%ld\n",Ntot,NBG));
458 return Ntot;
459}
Definition: EdbLayer.h:39
float Ymin() const
Definition: EdbLayer.h:87
float Ymax() const
Definition: EdbLayer.h:88
float Xmax() const
Definition: EdbLayer.h:86
float Xmin() const
Definition: EdbLayer.h:85
float Z() const
Definition: EdbLayer.h:77
void PropagateTo(float z)
Definition: EdbSegP.cxx:293
Float_t Z() const
Definition: EdbPattern.h:84
EdbSegP * AddSegment(int i, EdbSegP &s)
Definition: EdbPattern.cxx:72
Definition: Side.h:11
EdbSegP * s1
Definition: tlg2couples.C:29
long NBG
Definition: mc2raw.cxx:39
void w(int rid=2, int nviews=2)
Definition: test.C:27

◆ GetPattern()

EdbPattern * EdbDataStore::GetPattern ( int  n,
bool  btk = true 
)
inline
49{return GetPV(btk)->GetPattern(n);}
EdbPatternsVolume * GetPV(bool btk=true)
Definition: EdbDataStore.h:50

◆ GetPV()

EdbPatternsVolume * EdbDataStore::GetPV ( bool  btk = true)
inline
50{return btk?(&eSegPV):(&eRawPV);}

◆ GetRawPat()

EdbPattern * EdbDataStore::GetRawPat ( int  n)
inline
48{return (n<eRawPV.Npatterns())?eRawPV.GetPattern(n):0;}

◆ GetSegPat()

EdbPattern * EdbDataStore::GetSegPat ( int  n)
inline
47{return (n<eSegPV.Npatterns())?eSegPV.GetPattern(n):0;}

◆ GetTrack()

EdbTrackP * EdbDataStore::GetTrack ( int  n)
inline

get methods

45{return (n<eTracks.GetEntries())?(EdbTrackP*)eTracks.At(n):0;}

◆ GetVertex()

EdbVertex * EdbDataStore::GetVertex ( int  n)
inline
46{return (n<eVTX.GetEntries())?(EdbVertex*)eVTX.At(n):0;}

◆ LoadMCVertices()

void EdbDataStore::LoadMCVertices ( TObjArray *  vtx)

restore MC info methods


18 {
19 Log(2,"EdbDataStore::LoadMCVertices",Form("Reading %d vertices",vtx->GetEntries()));
20 for(int nv=0; nv<vtx->GetEntries();++nv){
21 AddVertex((EdbVertex*)vtx->At(nv));
22 }
23}
void AddVertex(EdbVertex *v)
Definition: EdbDataStore.h:63

◆ MakePattern()

void EdbDataStore::MakePattern ( double  z,
int  plate,
int  side 
)

264 {
265 Log(3,"EdbDataStore::MakePattern()",Form("plt#%d (side %d) z=%2.1f",plate,side,z));
266 EdbPattern* p=new EdbPattern(0,0,z,0);
267 p->SetID(plate);
268 p->SetScanID(EdbID(0,plate,0,0));
269 p->SetSide(side);
270 p->SetPID(side?(plate*10+side):plate);
271 if(side==0)eSegPV.AddPattern(p);
272 else eRawPV.AddPattern(p);
273
274}
Definition: EdbID.h:7

◆ Nplt()

int EdbDataStore::Nplt ( )
inline
39{return eRawPV.Npatterns()/2;}

◆ Nt()

int EdbDataStore::Nt ( )
inline

count methods:

37{return eTracks.GetEntries();}

◆ Nv()

int EdbDataStore::Nv ( )
inline
38{return eVTX.GetEntries();}

◆ PrintBrief()

void EdbDataStore::PrintBrief ( )

print methods


340 {
341 printf("EdbDataStore. We have here:\n");
342 printf(" --%d vertcies\n --%d tracks\n --%d Plates (geometry)\n",eVTX.GetEntries(),eTracks.GetEntries(),(eBrick)?eBrick->Npl():0);
343 printf(" --%d BT Patterns\n --%d MT Patterns\n",eSegPV.Npatterns(),eRawPV.Npatterns());
344}

◆ PrintPatterns()

void EdbDataStore::PrintPatterns ( )

362 {
363 printf("EdbDataStore. We have %d BT patterns:\n",eSegPV.Npatterns());
364 EdbPattern* pat=0;
365 for(int np=0;np<eSegPV.Npatterns();++np){
366 pat=GetSegPat(np);
367 printf(" (pat #%3d) [X,Y,Z]=[% 7.2f % 7.2f % 7.2f]\t Nseg=%d\n",
368 pat->Plate(),pat->X(),pat->Y(),pat->Z(),pat->N());
369 }
370}
Float_t Y() const
Definition: EdbPattern.h:83
Float_t X() const
Definition: EdbPattern.h:82

◆ PrintTracks()

void EdbDataStore::PrintTracks ( int  vlev = 0)

346 {
347 printf("EdbDataStore. We have %d tracks\n",eTracks.GetEntries());
348 if(vlev<1)return;
349 EdbTrackP* trk=0;
350 for(int nt=0;nt<eTracks.GetEntries();++nt){
351 trk=GetTrack(nt);
352 printf(" (trk #%3d) PDG=% 4d P=% 4.2f [X,Y,Z]=[% 7.2f % 7.2f % 7.2f]\t[TX,TY]=[% 2.2f % 2.2f] Nseg=%d\n",
353 trk->ID(),trk->PDG(),trk->P(),trk->X(),trk->Y(),trk->Z(),trk->TX(),trk->TY(),trk->N());
354 if(vlev<2)continue;
355 for(int ns=0;ns<trk->N();++ns){
356 EdbSegP* s=trk->GetSegment(ns);
357 printf(" seg#%d Plate=%d pid=%d flag=%d\n",s->ID(),s->Plate(),s->PID(),s->Flag());
358 }
359 }
360}
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
Float_t TY() const
tangens = deltaY/deltaZ
Definition: EdbSegP.h:176
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:195

◆ Restore_PatFromGeom()

void EdbDataStore::Restore_PatFromGeom ( int  np0 = 0,
int  np1 = 1000 
)

fill patterns array from eBrick plates

144 {
146 assert(eBrick);
147 EdbPlateP* plt=0;
148 for(int np=0;np<eBrick->Npl();++np){
149 if(np<np0 || np>np1)continue;
150 plt=eBrick->GetPlate(np);
151 assert(plt);
152 for(int side=0;side<3;++side)
153 MakePattern(plt->Z()+plt->GetLayer(side)->Z(),np,side);
154 }
155}
void MakePattern(double z, int plate, int side)
Definition: EdbDataStore.cxx:264

◆ Restore_PIDFromID()

void EdbDataStore::Restore_PIDFromID ( )

restore patterns' "plate" and "side" from their ID

105 {
107 EdbPattern* p;
108 int plate, side;
109 for(int np=0;np<eSegPV.Npatterns();np++){
110 p=eSegPV.GetPattern(np);
111 plate=p->ID();
112 p->SetPID(p->ID());
113 p->SetScanID(EdbID(0,plate,0,0));
114 p->SetSide(0);
115 }
116 for(int np=0;np<eRawPV.Npatterns();np++){
117 p=eRawPV.GetPattern(np);
118 plate=p->ID()/2;
119 p->SetPID(p->ID());
120 side=1+p->ID()%2;
121 p->SetScanID(EdbID(0,plate,0,0));
122 p->SetSide(side);
123 }
124};

◆ Restore_SegFromTrx()

void EdbDataStore::Restore_SegFromTrx ( EdbSegmentCut cut = 0,
int  Plt0 = 0,
int  Plt1 = 1000 
)

restore segments from tracks

171 {
173 EdbTrackP* t=0;
174 EdbSegP* s=0;
175 for(int nt=0;nt<eTracks.GetEntries();++nt){
176 t=GetTrack(nt);
177 if(t->N()==0)continue;
178 for(int ns=0;ns<t->N();++ns){
179 s=t->GetSegment(ns);
180 s->SetMC(t->MCEvt(),t->MCTrack());
181// s->SetFlag(1);
182 AddSegment(s,cut,Plt0,Plt1);
183 }
184 }
185}
EdbSegP * AddSegment(EdbSegP *seg, EdbSegmentCut *cut=0, int Plt0=0, int Plt1=100)
add methods
Definition: EdbDataStore.cxx:126
TTree * t
Definition: check_shower.C:4

◆ Restore_TrxFromVtx()

void EdbDataStore::Restore_TrxFromVtx ( )

restore tracks from vertex

157 {
159 ClearTracks();
160 EdbVertex* v;
161 for(int nv=0;nv<eVTX.GetEntries();nv++){
162 v=GetVertex(nv);
163 for(int nt=0;nt<v->N();nt++){
164 EdbTrackP* trk=v->GetTrack(nt);
165 trk->AddVTA(v->GetVTa(nt));
166 if(v->GetVTa(nt)->Zpos()==1)AddTrack(trk);
167 }
168 }
169}
void AddTrack(EdbTrackP *tr)
Definition: EdbDataStore.h:62
void AddVTA(EdbVTA *vta)
Definition: EdbPattern.cxx:450
Int_t Zpos() const
Definition: EdbVertex.h:47
EdbVTA * GetVTa(int i)
Definition: EdbVertex.h:139
Int_t N() const
Definition: EdbVertex.h:121
EdbTrackP * GetTrack(int i)
Definition: EdbVertex.h:141

◆ SavePlateToRaw()

void EdbDataStore::SavePlateToRaw ( const char *  fname,
int  PID,
Option_t *  option = "RECREATE" 
)

Convert this event to RAW data format!

Get plate geometry:


maximum segments per view. When this is reached - create a new view.

761 {
763 EdbSegP* seg=0;
764 EdbSegment* seg1=0;
765 EdbPattern* ptop=FindPattern(PID,2);
766 EdbPattern* pbot=FindPattern(PID,1);
767 if(ptop->N()==0 && pbot->N()==0)return;
768
772 if(strcmp(option,"UPDATE")==0){
773 //test if file exists
774 FILE *file= fopen(fname, "r");
775 if (file) fclose(file);
776 else{
777 option="RECREATE";}
778 }
779 EdbRun run(fname,option);
780 Log(1,"EdbDataStore::SavePlateToRaw",Form("Open file \"%s\" opt=\"%s\"",fname,option));
781 Log(2,"EdbDataStore::SavePlateToRaw","created run");
782 Log(2,"EdbDataStore::SavePlateToRaw",Form("cycle through %d TOP segs",ptop->N()));
783 int Vid=0;
784 const int NsegPerView=1000;
785 int Ns=0;
786 run.GetView()->SetNframes(1,0);
787 run.GetView()->GetHeader()->SetViewID(++Vid);
788 for(int ns=0;ns<ptop->N();++ns){
789 if(Ns==NsegPerView){
790 Ns=0;
791 run.AddView();
792 run.GetView()->Clear();
793 run.GetView()->GetHeader()->SetViewID(++Vid);
794 }
795 seg=ptop->GetSegment(ns);
796 if(seg->W()<0)continue;
797 seg1=run.GetView()->AddSegment(p->Xp(*seg),p->Yp(*seg),seg->Z(),p->TXp(*seg),p->TYp(*seg),seg->DZ(),0,(int)(seg->W()),seg->ID());
798 seg1->SetSigma(seg->P(),1);
799 ++Ns;
800 }
801 run.AddView();
802 run.GetView()->Clear();
803
804 Log(2,"EdbDataStore::SavePlateToRaw",Form("cycle through %d BOT segs",pbot->N()));
805 run.GetView()->SetNframes(0,1);
806 run.GetView()->GetHeader()->SetViewID(++Vid);
807 Ns=0;
808 for(int ns=0;ns<pbot->N();++ns){
809 if(Ns==NsegPerView){
810 Ns=0;
811 run.AddView();
812 run.GetView()->Clear();
813 run.GetView()->GetHeader()->SetViewID(++Vid);
814 }
815 seg=pbot->GetSegment(ns);
816 if(seg->W()<0)continue;
817 seg1=run.GetView()->AddSegment(p->Xp(*seg),p->Yp(*seg),seg->Z(),p->TXp(*seg),p->TYp(*seg),seg->DZ(),1,(int)(seg->W()),seg->ID());
818 seg1->SetSigma(seg->P(),1);
819 ++Ns;
820 }
821 run.AddView();
822 run.GetView()->Clear();
823 Log(2,"EdbDataStore::SavePlateToRaw","Close plate");
824 run.Save();
825 run.Close();
826};
Definition: EdbRun.h:75
void Close()
Definition: EdbRun.cxx:439
EdbView * GetView() const
Definition: EdbRun.h:110
void Save()
Definition: EdbRun.cxx:424
void AddView()
Definition: EdbRun.cxx:305
Float_t DZ() const
Definition: EdbSegP.h:154
Float_t W() const
Definition: EdbSegP.h:151
segment of the track
Definition: EdbSegment.h:63
void SetSigma(float sx, float sy)
Definition: EdbSegment.h:85
void SetViewID(int id)
Definition: EdbView.h:99
EdbViewHeader * GetHeader() const
Definition: EdbView.h:163
EdbSegment * AddSegment(float x, float y, float z, float tx, float ty, float dz=0, int side=0, int puls=0, int id=-1)
Definition: EdbView.h:231
void Clear()
Definition: EdbView.cxx:79
void SetNframes(int top, int bot)
Definition: EdbView.h:184
EdbRun * run
Definition: check_raw.C:38
TFile * file
Definition: write_pvr.C:3
const char * fname
Definition: mc2raw.cxx:41
fclose(pFile)
Definition: Flexmotn.h:102

◆ SaveToRaw()

void EdbDataStore::SaveToRaw ( const char *  dir = "./",
const EdbID idset = "0.0.0.0",
Option_t *  option = "RECREATE",
bool  doaff = true 
)

save methods:


Save the raw data to FEDRA brick structure

generate plate filename

save plate.

save affine tranformations

save affine parameters

828 {
830 EdbScanProc sp;
831
832 sp.eProcDirClient=dir;
835
836 EdbID PL=idset;
837 EdbID PL0=idset;
838 EdbPlateP* p=0;
839
840 float z,z0;
841 for(int i=0;i<eBrick->Npl();i++){
842 p=eBrick->GetPlate(i);
843 z=p->Z();
844 PL.ePlate=i+1;
845 sp.CheckPlateDir(PL);
846 TString fnm;
847 sp.MakeFileName(fnm,PL,"raw.root",true);
848 SavePlateToRaw(fnm.Data(),i,option);
849 if(doaff && i>0){
851 sp.MakeAffName(fnm,PL,PL0);
852 FILE* afF=fopen(fnm.Data(),"w");
853 fprintf(afF,"ZLAYER 0 %4f 0 0\n",z0-z);
854 fclose(afF);
855 }
856 PL0=PL;
857 z0=z;
858 }
859}
brick z0
Definition: RecDispMC.C:106
void SavePlateToRaw(const char *fname, int PID, Option_t *option="RECREATE")
Definition: EdbDataStore.cxx:761
Int_t eBrick
Definition: EdbID.h:10
Int_t ePlate
Definition: EdbID.h:11
scanned data processing
Definition: EdbScanProc.h:12
bool CheckBrickDir(EdbID id, bool create=true)
Definition: EdbScanProc.cxx:1780
bool CheckAFFDir(int brick, bool create=true)
Definition: EdbScanProc.cxx:1771
void MakeAffName(TString &s, int id1[4], int id2[4], const char *suffix="aff.par")
Definition: EdbScanProc.cxx:1851
void MakeFileName(TString &s, int id[4], const char *suffix, bool inplate=true)
Definition: EdbScanProc.cxx:1819
bool CheckPlateDir(EdbID id, bool create=true)
Definition: EdbScanProc.cxx:1790
TString eProcDirClient
directory path for root data
Definition: EdbScanProc.h:14
EdbID idset
Definition: emrec.cpp:35

◆ SetOwnTracks()

void EdbDataStore::SetOwnTracks ( bool  own = true)
inline

setown methods:

41{eTracks.SetOwner(own);}

◆ SetOwnTrkSegs()

void EdbDataStore::SetOwnTrkSegs ( )
inline
43{for(int nt=0;nt<Nt();++nt)GetTrack(nt)->SetOwner();}
void SetOwner()
Definition: EdbPattern.h:139

◆ SetOwnVertices()

void EdbDataStore::SetOwnVertices ( bool  own = true)
inline
42{eVTX.SetOwner(own);}

◆ TransferGeometry()

void EdbDataStore::TransferGeometry ( EdbDataStore ds)

25 {
26 Log(2,"transfer geometry","");
27 ds->eBrick->Clear();
28 ds->eBrick->Copy(*eBrick);
29}
virtual void Clear()
Definition: EdbDisplayBase.h:142
EdbDisplay * ds
Definition: check_vertex.C:16

◆ TransferSegs()

void EdbDataStore::TransferSegs ( EdbPatternsVolume pv0,
EdbPatternsVolume pv1,
EdbSegmentCut cut = 0,
int  FromPlate = 0,
int  ToPlate = 57 
)
static

find corresponding DEST pattern

62 {
63 Log(1,"EdbDataStore::TransferSegs()","Transfer segments");
64 assert(pv0);
65 assert(pv1);
66 EdbPattern* pat0=0;
67 EdbPattern* pat1=0;
68 EdbSegP* seg=0;
69 int npat=pv0->Npatterns();
70
71 if(pv1->Npatterns()==0)pv0->PassProperties(*pv1);
72
73 int nseg=0;
74 float p[5];
75 for(int np=0;np<npat;++np){
76
77 pat0=pv0->GetPattern(np);
78 if(pat0->Plate()<FromPlate)continue;
79 if(pat0->Plate()>ToPlate)continue;
80 pat1=pv1->GetPattern(np);
81 assert(pat1);
82 assert(pat1->Z()==pat0->Z());
83
84 nseg=pat0->N();
85 for(int ns=0;ns<nseg;++ns){
86 seg=pat0->GetSegment(ns);
87 p[0]=seg->X();
88 p[1]=seg->Y();
89 p[2]=seg->TX();
90 p[3]=seg->TY();
91 p[4]=seg->W();
92 if(cut && !cut->PassCut(p))continue;
93 seg=pat1->AddSegment(*seg);
94 seg->SetTrack(0);
95 };
96// printf("dest.pat#%d has #%d segments\n",pat1->ID(),pat1->N());
97 pat1->SetID(np);
98 pat1->SetSegmentsPID();
99
100 }
101
102}
Int_t npat
Definition: Xi2HatStartScript.C:33
void SetSegmentsPID()
Definition: EdbPattern.cxx:1455
void SetID(int id)
Definition: EdbPattern.h:309
void PassProperties(EdbPatternsVolume &pvol)
Definition: EdbPattern.cxx:1729
void SetTrack(int trid)
Definition: EdbSegP.h:131

◆ TransferTo()

void EdbDataStore::TransferTo ( EdbDataStore ds,
char  level,
EdbSegmentCut cut = 0,
int  FromPlate = 0,
int  ToPlate = 57 
)

transfer methods


31 {
32 assert(level>=0 && level<16);
33
34 if(level&0x1){
35 Log(2,"EdbDataStore::TransferTo","transfer MTK segments\n");
36 TransferSegs(&eRawPV,&(ds->eRawPV),cut,FromPlate,ToPlate);
37 Log(2,"EdbDataStore::TransferTo","done");
38 }
39
40 if(level&0x2){
41 Log(2,"EdbDataStore::TransferTo","transfer BTK segments\n");
42 TransferSegs(&eSegPV,&(ds->eSegPV),cut,FromPlate,ToPlate);
43 Log(2,"EdbDataStore::TransferTo","done");
44 }
45
46 if(level&0x4){
47 Log(2,"EdbDataStore::TransferTo","transfer Tracks\n");
48 for(int nt=0;nt<eTracks.GetEntries();++nt){
49 if(GetTrack(nt)->N())ds->AddTrack(GetTrack(nt));
50 }
51 ds->SetOwnTracks(0);
52 Log(2,"EdbDataStore::TransferTo","done");
53 }
54 if(level&0x8){
55 Log(2,"EdbDataStore::TransferTo","transfer Vertices\n");
56 for(int nv=0;nv<eVTX.GetEntries();++nv)ds->AddVertex(GetVertex(nv));
57 ds->SetOwnVertices(0);
58 Log(2,"EdbDataStore::TransferTo","done");
59 }
60}
static void TransferSegs(EdbPatternsVolume *pv0, EdbPatternsVolume *pv1, EdbSegmentCut *cut=0, int FromPlate=0, int ToPlate=57)
Definition: EdbDataStore.cxx:62

Member Data Documentation

◆ eBrick

EdbBrickP* EdbDataStore::eBrick

◆ eRawPV

EdbPatternsVolume EdbDataStore::eRawPV

geometry

◆ eSegPV

EdbPatternsVolume EdbDataStore::eSegPV

◆ eTracks

TObjArray EdbDataStore::eTracks

◆ eVTX

TObjArray EdbDataStore::eVTX

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