FEDRA emulsion software from the OPERA Collaboration
EdbViewGen Class Reference

#include <EdbViewGen.h>

Inheritance diagram for EdbViewGen:
Collaboration diagram for EdbViewGen:

Public Member Functions

 EdbViewGen ()
 
 EdbViewGen (EdbViewDef &vd)
 
int GenAlfaGrains ()
 
int GenDeltaGrains ()
 
int GenFogGrains (int ngr, TObjArray &grains, int side=0)
 
int GenFrames (EdbView &v)
 
int GenGrainClusters (EdbView &v, EdbCluster &g)
 
int GenGrains (EdbView &v)
 
int GenSegGrains (EdbSegment &s)
 
float GrainPathMip (float lambda=2.)
 
 ~EdbViewGen ()
 
- Public Member Functions inherited from EdbViewDef
 EdbViewDef ()
 
void Print ()
 
void SetDef ()
 
 ~EdbViewDef ()
 

Additional Inherited Members

- Public Attributes inherited from EdbViewDef
Float_t eClaSX
 
Float_t eClaSY
 
Float_t eClaSZ
 
Float_t eClaSZvar
 
Float_t eDZdead
 
Float_t eFogDens
 
Float_t eFogGrainArea
 
Float_t eGrainArea
 
Float_t eGrainSX
 
Float_t eGrainSY
 
Float_t eGrainSZ
 
Int_t eNframes
 
Float_t eTXopt
 
Float_t eTYopt
 
Float_t eX0
 
Float_t eX0opt
 
Float_t eXmax
 
Float_t eXmin
 
Float_t eY0
 
Float_t eY0opt
 
Float_t eYmax
 
Float_t eYmin
 
Float_t eZdead
 
Float_t eZmax
 
Float_t eZmin
 
Float_t eZxy
 

Constructor & Destructor Documentation

◆ EdbViewGen() [1/2]

EdbViewGen::EdbViewGen ( )

◆ EdbViewGen() [2/2]

EdbViewGen::EdbViewGen ( EdbViewDef vd)
inline
18: EdbViewDef(vd) {}
void vd(int trmin=2, float amin=.0)
Definition: check_vertex.C:217

◆ ~EdbViewGen()

EdbViewGen::~EdbViewGen ( )
inline
19{}

Member Function Documentation

◆ GenAlfaGrains()

int EdbViewGen::GenAlfaGrains ( )
inline
26{return 0;} // bg from alfa-particles

◆ GenDeltaGrains()

int EdbViewGen::GenDeltaGrains ( )
inline
27{return 0;} // bg from delta-electrons

◆ GenFogGrains()

int EdbViewGen::GenFogGrains ( int  ngr,
TObjArray &  grains,
int  side = 0 
)
213{
214 float x,y,z;
215 int area;
216 int n=0;
217
218 printf("GenFogGrains: X:(%f %f) Y:(%f %f) \n", eXmin,eXmax, eYmin,eYmax );
219 while( n<ngr ) {
220 x = eXmin + gRandom->Rndm()*(eXmax-eXmin);
221 y = eYmin + gRandom->Rndm()*(eYmax-eYmin);
222 z = eZmin + gRandom->Rndm()*(eZmax-eZmin);
223
224 area = (int)(eFogGrainArea); // + gRandom->Poisson(4) - 2);
225 grains.Add( new EdbCluster(x,y,z, area, 0, 0, side, -1) );
226 n++;
227 }
228 return n;
229 }
Definition: EdbCluster.h:19
Float_t eYmin
Definition: EdbViewDef.h:21
Float_t eFogGrainArea
Definition: EdbViewDef.h:33
Float_t eXmax
Definition: EdbViewDef.h:20
Float_t eZmax
Definition: EdbViewDef.h:17
Float_t eZmin
Definition: EdbViewDef.h:17
Float_t eYmax
Definition: EdbViewDef.h:21
Float_t eXmin
Definition: EdbViewDef.h:20
grains()
Definition: grains.C:3

◆ GenFrames()

int EdbViewGen::GenFrames ( EdbView v)
28{
30 float step=(eZmax-eZmin-1.)/(eNframes-1);
31 for(int i=0; i<eNframes; i++)
32 v.AddFrame( i, eZmin + 0.5 + i*step );
33
34 return eNframes;
35 }
Int_t eNframes
Definition: EdbViewDef.h:16
void AddFrame(int id, float z, int ncl=0, int npix=0)
Definition: EdbView.cxx:268
void SetNframes(int top, int bot)
Definition: EdbView.h:184

◆ GenGrainClusters()

int EdbViewGen::GenGrainClusters ( EdbView v,
EdbCluster g 
)
86{
87 // generate clusters for the input grain g
88
89 int ncl=0;
90 if( g.GetArea()==0 ) return ncl;
91 int nfr= v.GetNframes();
92 float z;
93 float dz = 0.5 * gRandom->Gaus(eGrainSZ,eGrainSZ/3.);
94 //eClaSZ*TMath::Sqrt(g.GetArea())/3; //ToDo
95 float dx,dy;
96 EdbCluster *c;
97 for(int i=0; i<nfr; i++) {
98 z= v.GetFrame(i)->GetZ();
99 if( TMath::Abs(z - g.GetZ()) < dz ) {
100 c= new EdbCluster(g);
101 gRandom->Rannor(dx,dy);
102 c->SetZ(z);
103 c->SetX( g.X()+ dx*eClaSX ); // TODO: navesti nauku
104 c->SetY( g.Y()+ dy*eClaSY );
105 c->SetFrame(i);
106 v.AddCluster(c);
107 ncl++;
108 }
109 }
110 //printf("ncl = %d\n",ncl);
111 return ncl;
112}
brick dz
Definition: RecDispMC.C:107
virtual Float_t Y() const
Definition: EdbCluster.h:79
virtual Float_t X() const
Definition: EdbCluster.h:78
void SetFrame(int f)
Definition: EdbCluster.h:47
virtual void SetX(float x)
Definition: EdbCluster.h:81
virtual void SetY(float y)
Definition: EdbCluster.h:82
Float_t GetArea() const
Definition: EdbCluster.h:54
virtual void SetZ(float z)
Definition: EdbCluster.h:83
Float_t GetZ() const
Definition: EdbCluster.h:53
float GetZ() const
Definition: EdbFrame.h:42
Float_t eClaSX
Definition: EdbViewDef.h:28
Float_t eGrainSZ
Definition: EdbViewDef.h:30
Float_t eClaSY
Definition: EdbViewDef.h:28
EdbFrame * GetFrame(int i) const
Definition: EdbView.h:221
Int_t GetNframes() const
Definition: EdbView.h:206
EdbCluster * AddCluster(EdbCluster *c)
Definition: EdbView.h:225

◆ GenGrains()

int EdbViewGen::GenGrains ( EdbView v)
39{
40 TObjArray grains;
41 int nseg = v.Nsegments();
42 //printf("nsegments = %d\n", nseg);
43 EdbSegment *s = 0;
44 int ngrseg=0;
45 int ngr;
46 for(int iseg=0; iseg<nseg; iseg++) {
47 s = v.GetSegment(iseg);
48 ngr = GenSegGrains( *s );
49 printf("ngr = %d\n", ngr);
50 for(int igr=0; igr<ngr; igr++) grains.Add(s->GetElements()->At(igr));
51 ngrseg+=ngr;
52 }
53 //printf("ngrseg = %d\n", ngrseg);
54
55 int nfog = (int)(eFogDens*(eXmax-eXmin)*(eYmax-eYmin)*(eZmax-eZmin)/1000.);
56 printf("ngrfog = %d\n", nfog);
57 GenFogGrains(nfog,grains);
58
59 EdbCluster *g;
60 ngr = grains.GetEntries();
61 printf("%d grains generated\n",ngr);
62
63 // remove grains from the dead layer
64 int icrem=0;
65 for(int i=0; i<ngr; i++) {
66 g = (EdbCluster*)grains.At(i);
67 if( TMath::Abs(g->GetZ()-eZdead)<eDZdead/2. ) { g->SetArea(0); icrem++;}
68 }
69 printf("%d grains removed in dead layer\n",icrem);
70
71 for(int i=0; i<ngr; i++) {
72 g = (EdbCluster*)grains.At(i);
73 GenGrainClusters(v,*g);
74 //v.AddCluster(new EdbCluster(*g));
75 }
76
77 for(int iseg=0; iseg<nseg; iseg++)
78 if( v.GetSegment(iseg)->GetElements() )
79 v.GetSegment(iseg)->GetElements()->Delete();
80
81 return ngr;
82}
void SetArea(float a)
Definition: EdbCluster.h:45
segment of the track
Definition: EdbSegment.h:63
TObjArray * GetElements() const
Definition: EdbSegment.h:117
Float_t eFogDens
Definition: EdbViewDef.h:32
Float_t eDZdead
Definition: EdbViewDef.h:37
Float_t eZdead
Definition: EdbViewDef.h:36
int GenFogGrains(int ngr, TObjArray &grains, int side=0)
Definition: EdbViewGen.cxx:212
int GenSegGrains(EdbSegment &s)
Definition: EdbViewGen.cxx:115
int GenGrainClusters(EdbView &v, EdbCluster &g)
Definition: EdbViewGen.cxx:85
EdbSegment * GetSegment(int i) const
Definition: EdbView.h:219
Int_t Nsegments() const
Definition: EdbView.h:216
s
Definition: check_shower.C:55

◆ GenSegGrains()

int EdbViewGen::GenSegGrains ( EdbSegment s)
116{
117 // input: segments of the view, assumed that puls of the
118 // segment is correspond to the number of grains to be generated
119 // output: array of grains belong to the segment
120 // ToDo: dz only positive?
121
122 using namespace TMath;
123
124 // could be switched on if serve...
125 //
126 // int ngr = s.GetPuls();
127 // TArrayF arr(ngr);
128 // gRandom->RndmArray(ngr, arr.GetArray());
129 // for(int i=0; i<ngr; i++) printf("%f ",arr[i]);
130 // printf("\n");
131 // TArrayI ind(ngr);
132 // Sort( n, arr.GetArray(), ing.GetArray(), 0 );
133
134 float x = s.GetX0();
135 float y = s.GetY0();
136 float z = s.GetZ0();
137 float dz = s.GetDz();
138 float zmax=z+dz;
139 float cs = Cos(ATan(Sqrt(s.GetTx()*s.GetTx()+s.GetTy()*s.GetTy())));
140 int area;
141 float step;
142 int ncl=0;
143 int nmax=10000;
144
145 float lambda = dz/s.GetPuls();
146
147 for(int i=0; i<nmax; i++ ) {
148 step = GrainPathMip(lambda)*cs;
149 z += step;
150 //printf("Zstep: %5.3f \n",step);
151 x += s.GetTx()*step;
152 y += s.GetTy()*step;
153 if( z > zmax ) break;
154 if( x < eXmin ) continue;
155 if( x > eXmax ) continue;
156 if( y < eYmin ) continue;
157 if( y > eYmax ) continue;
158 area = (int)(eGrainArea); // + gRandom->Poisson(4) - 2);
159
160 double sx,sy,sz, r=gRandom->Gaus(0.,eGrainSX); // the smearing of the physical grain center position
161 gRandom->Sphere(sx,sy,sz,r);
162 x+=sx; y+=sy; z+=sz;
163 if(x<eXmin) continue;
164 if(x>eXmax) continue;
165 if(y<eYmin) continue;
166 if(y<eYmin) continue;
167 s.AddElement( new EdbCluster( x+sx, y+sy, z+sz, area, 0, 0, s.GetSide(), s.GetID()) );
168 ncl++;
169 }
170 s.SetPuls(ncl);
171 return ncl;
172 }
Float_t eGrainArea
Definition: EdbViewDef.h:34
Float_t eGrainSX
Definition: EdbViewDef.h:30
float GrainPathMip(float lambda=2.)
Definition: EdbViewGen.cxx:175
void r(int rid=2)
Definition: test.C:201

◆ GrainPathMip()

float EdbViewGen::GrainPathMip ( float  lambda = 2.)
176{
177 //=========================================
178 //
179 // Some emulsion matter properties:
180 // density: 2.4 g/cm^3
181 // <A> = 18.2
182 // <Z> = 8.9
183 // X0 = 5.5 cm
184 // nuclear collision lenght: 33 cm
185 //
186 // <dE/dX> = 37 KeV/100 microns
187 //
188 // visible grain density: 30/100 microns
189 // grain size: 0.6 micron
190 //
191 // assuming the ionization energy about 4 eV (to check!)
192 // we expect about 100 ionizations/micron
193 //
194 //=========================================
195
196 // lambda // mean path length
197
198 int nstep=100;
199 double prob0 = 1./nstep; //prob of grain after single ionization
200 float grsiz=gRandom->Uniform(0.9,1.1);
201 float x=0;
202 for(int j=0; j<100000; j++) {
203 x+=gRandom->Exp(lambda/nstep);
204 if( x>grsiz&&gRandom->Rndm()<prob0 ) {
205 break;
206 }
207 }
208 return x;
209}

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