FEDRA emulsion software from the OPERA Collaboration
EdbScanSet Class Reference

#include <EdbScanSet.h>

Inheritance diagram for EdbScanSet:
Collaboration diagram for EdbScanSet:

Public Member Functions

Bool_t AddID (EdbID *id, Int_t step)
 
Int_t AssembleBrickFromPC ()
 
EdbBrickPBrick ()
 
const EdbBrickPBrick () const
 
void Copy (EdbScanSet &sc)
 
 EdbScanSet (EdbID id)
 
 EdbScanSet (int brickid=0)
 
EdbIDFindNextPlateID (Int_t p, Bool_t direction)
 
EdbIDFindPlateID (Int_t p)
 
Bool_t GetAffP2P (Int_t p1, Int_t p2, EdbAffine2D &aff)
 
Float_t GetDZP2P (Int_t p1, Int_t p2)
 
EdbIDGetID (Int_t i)
 
EdbLayerGetLayer (Int_t p, Int_t side)
 
EdbPlatePGetPlate (Int_t p)
 
void MakeNominalSet (EdbID id, int from_plate=57, int to_plate=1, float z0=0, float dz=-1300, float shr=1, float dzbase=210., float dzem=45.)
 
Int_t MakeParFiles (Int_t piece=0, const char *dir=".")
 
void MakePIDList ()
 
void Print ()
 
void ReadGeom (const char *fname)
 
Int_t ReadIDS (const char *file)
 
void RemovePlate (int pid)
 
Bool_t SetAsReferencePlate (Int_t pid)
 
void SetID (EdbID id)
 
Int_t ShiftBrickZ (Float_t z)
 
Int_t TransformBrick (EdbAffine2D aff)
 
void TransformBrick (EdbScanSet &ss)
 
Int_t TransformSidesIntoBrickRS ()
 
void UpdateBrickWithP2P (EdbLayer &la, int plate1, int plate2)
 
void UpdateGap (float dz, int pid1, int pid2)
 
void UpdateIDS (int brick, int ma, int mi)
 
bool ValidSide (int side) const
 
void WriteGeom (const char *fname)
 
Int_t WriteIDS (const char *file=0)
 
float Zlayer (int plate, int side)
 
virtual ~EdbScanSet ()
 

Public Attributes

EdbBrickP eB
 all layers of the brick defined here More...
 
EdbID eID
 id of the scanset itself More...
 
TList eIDS
 list of the identifiers to be processed More...
 
TObjArray ePC
 Plate Couples. Each couple represented as EdbPlateP where 2 sides in reality corresponds to 2 plates. More...
 
TIndexCell ePID
 correspondance between index in eB and the plate id More...
 
Int_t eReferencePlate
 

Constructor & Destructor Documentation

◆ EdbScanSet() [1/2]

EdbScanSet::EdbScanSet ( EdbID  id)

◆ EdbScanSet() [2/2]

EdbScanSet::EdbScanSet ( int  brickid = 0)
30{
31 eB.SetID(id);
34}
Int_t eBrick
Definition: EdbID.h:10
void SetID(int id)
Definition: EdbLayer.h:94
EdbBrickP eB
all layers of the brick defined here
Definition: EdbScanSet.h:13
EdbID eID
id of the scanset itself
Definition: EdbScanSet.h:19
Int_t eReferencePlate
Definition: EdbScanSet.h:18
UInt_t id
Definition: tlg2couples.C:117

◆ ~EdbScanSet()

virtual EdbScanSet::~EdbScanSet ( )
inlinevirtual
24{}

Member Function Documentation

◆ AddID()

Bool_t EdbScanSet::AddID ( EdbID id,
Int_t  step 
)
477{
478 if (FindPlateID(id->GetPlate())) return kFALSE;
479 for (EdbID *it = (EdbID*)eIDS.First(); it ; it = (EdbID*)eIDS.After(it)) {
480 if (step < 0 && id->GetPlate() > it->GetPlate()) {
481 eIDS.AddBefore(it, id);
482 return kTRUE;
483 }
484 if (step > 0 && id->GetPlate() < it->GetPlate()) {
485 eIDS.AddBefore(it, id);
486 return kTRUE;
487 }
488 }
489 eIDS.Add(id);
490 return kTRUE;
491}
Definition: EdbID.h:7
EdbPlateP * GetPlate(Int_t p)
Definition: EdbScanSet.h:56
TList eIDS
list of the identifiers to be processed
Definition: EdbScanSet.h:16
EdbID * FindPlateID(Int_t p)
Definition: EdbScanSet.cxx:465

◆ AssembleBrickFromPC()

Int_t EdbScanSet::AssembleBrickFromPC ( )

input: ePC - couples of layers represented as "plates" output: brick with a first plate as a reference one

131{
134
135 ePID.Delete();
136 eB.Clear();
137 EdbPlateP *plate=0;
138
139 Float_t dz0=214,dz1=45,dz2=45; //TODO!
140 EdbAffine2D aff;
141 Int_t npc = ePC.GetEntries(); if(npc<1) return 0;
142 EdbPlateP *pc=0; //the couple of plates
143 Float_t z;
144 printf("\nAssembleBrickFromPC\n");
145 for(Int_t i=-1; i<npc; i++) {
146 plate = new EdbPlateP();
147
148 if(i==-1) {
149 pc = (EdbPlateP *)(ePC.At(0));
150 plate->SetID(pc->GetLayer(1)->ID());
151 z=0;
152 }
153 else {
154 pc = (EdbPlateP *)(ePC.At(i));
155 plate->SetID(pc->GetLayer(2)->ID());
156 z = eB.GetPlate(eB.Npl()-1)->Z() - pc->GetLayer(1)->Z();
157 aff.Reset();
158 aff.Transform( pc->GetLayer(1)->GetAffineXY() );
159 aff.Invert();
160 aff.Transform( eB.GetPlate(eB.Npl()-1)->GetAffineXY() );
161 plate->GetAffineXY()->Transform(&aff);
162 }
163 plate->SetZlayer(z, z - dz0/2 + dz1, z+dz0/2+dz2);
164 plate->GetLayer(0)->SetZlayer(0,-dz0/2,dz0/2); // internal plate coord
165 plate->GetLayer(2)->SetZlayer(-dz0/2,-dz0/2-dz2,-dz0/2);
166 plate->GetLayer(1)->SetZlayer( dz0/2, dz0/2, dz0/2+dz1);
167
168 eB.AddPlate(plate);
169 }
170 MakePIDList();
171
173
174 return eB.Npl();
175}
Definition: EdbAffine.h:17
void Invert()
Definition: EdbAffine.cxx:103
void Reset()
Definition: EdbAffine.cxx:72
void Transform(const EdbAffine2D *a)
Definition: EdbAffine.cxx:93
void Clear()
Definition: EdbBrick.h:54
EdbPlateP * GetPlate(int i)
Definition: EdbBrick.h:53
void AddPlate(EdbPlateP *pl)
Definition: EdbBrick.h:52
int Npl() const
Definition: EdbBrick.h:51
int ID() const
Definition: EdbLayer.h:73
EdbAffine2D * GetAffineXY()
Definition: EdbLayer.h:119
float Z() const
Definition: EdbLayer.h:77
Definition: EdbBrick.h:14
EdbLayer * GetLayer(int i)
Definition: EdbBrick.h:29
void MakePIDList()
Definition: EdbScanSet.cxx:178
TIndexCell ePID
correspondance between index in eB and the plate id
Definition: EdbScanSet.h:14
TObjArray ePC
Plate Couples. Each couple represented as EdbPlateP where 2 sides in reality corresponds to 2 plates.
Definition: EdbScanSet.h:15
void Delete()
Definition: TIndexCell.cpp:49
Int_t plate
Definition: merge_Energy_SytematicSources_Electron.C:1
Double_t Z
Definition: tlg2couples.C:104

◆ Brick() [1/2]

EdbBrickP & EdbScanSet::Brick ( )
inline
50{return eB;}

◆ Brick() [2/2]

const EdbBrickP & EdbScanSet::Brick ( ) const
inline
49{return eB;}

◆ Copy()

void EdbScanSet::Copy ( EdbScanSet sc)

copy from ss to this

52{
54 eB.Copy(sc.eB);
55 eID = sc.eID;
56
57 int npc = sc.ePC.GetEntries();
58 for(int i=0; i<npc; i++) {
59 EdbPlateP *p = new EdbPlateP();
60 p->Copy( *((EdbPlateP *)sc.ePC.At(i)) );
61 ePC.Add(p);
62 }
63
64 TIter next( &(sc.eIDS) );
65 TObject *obj=0;
66 while((obj = next()))
67 eIDS.Add( new EdbID( *((EdbID*)obj) ) );
68
71}
void Copy(EdbBrickP &b)
Definition: EdbBrick.cxx:85
p
Definition: testBGReduction_AllMethods.C:8

◆ FindNextPlateID()

EdbID * EdbScanSet::FindNextPlateID ( Int_t  p,
Bool_t  direction 
)

< find next valid identifier in respect to plate in positive or negative direction (+-1)

441{
443 Log(2,"FindNextPlateID","plate = %d direction = %d", plate, positive_direction);
444 EdbID *id = 0;
445 Int_t plmin=kMaxInt, plmax=kMinInt;
446 for (Int_t i = 0; i < eIDS.GetEntries(); i++) {
447 id = (EdbID *)eIDS.At(i);
448 if (id->ePlate>plmax) plmax=id->ePlate;
449 if (id->ePlate<plmin) plmin=id->ePlate;
450 }
451 int step = positive_direction? 1 : -1;
452 int first = plate +step;
453 int last = positive_direction? plmax:plmin; last += step;
454 printf("first=%d, last = %d step = %d\n",first,last, step);
455 if( positive_direction && first>=last) return 0;
456 if( !positive_direction && first<=last) return 0;
457 for(int i=first; i!=last; i+=step) {
458 id = FindPlateID(i);
459 if(id) return id;
460 }
461 return 0;
462}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75

◆ FindPlateID()

EdbID * EdbScanSet::FindPlateID ( Int_t  p)
466{
467 EdbID *id;
468 for (Int_t i = 0; i < eIDS.GetEntries(); i++) {
469 id = (EdbID *)eIDS.At(i);
470 if (id->ePlate == plate) return id;
471 }
472 return 0;
473}

◆ GetAffP2P()

Bool_t EdbScanSet::GetAffP2P ( Int_t  p1,
Int_t  p2,
EdbAffine2D aff 
)

input: p1 - id of the "from-plate" p2 - id of the "to-plate" output: aff - the transformation - applied to "from" segment became "to" segment

248{
252 Int_t ip1,ip2;
253 if(ePID.Find(p1)) ip1 = ePID.Find(p1)->At(0)->Value(); else return false;
254 if(ePID.Find(p2)) ip2 = ePID.Find(p2)->At(0)->Value(); else return false;
255 aff.Reset();
256 aff.Transform( eB.GetPlate(ip2)->GetAffineXY() );
257 aff.Invert();
258 aff.Transform( eB.GetPlate(ip1)->GetAffineXY() );
259 return true;
260}
Long_t Value() const
Definition: TIndexCell.h:79
TIndexCell const * At(Int_t narg, Int_t vind[]) const
Definition: TIndexCell.cpp:519
TIndexCell * Find(Int_t narg, Long_t varg[]) const
Definition: TIndexCell.cpp:613

◆ GetDZP2P()

Float_t EdbScanSet::GetDZP2P ( Int_t  p1,
Int_t  p2 
)

input: p1 - id of the "from-plate" p2 - id of the "to-plate" return: dz - segment of from-plate propagated to dz become in to-plate

264{
268 Int_t ip1,ip2;
269 if(ePID.Find(p1)) ip1 = ePID.Find(p1)->At(0)->Value(); else return false;
270 if(ePID.Find(p2)) ip2 = ePID.Find(p2)->At(0)->Value(); else return false;
271 return eB.GetPlate(ip2)->Z() - eB.GetPlate(ip1)->Z();
272}

◆ GetID()

EdbID * EdbScanSet::GetID ( Int_t  i)
inline
46{return (EdbID *)(eIDS.At(i));}

◆ GetLayer()

EdbLayer * EdbScanSet::GetLayer ( Int_t  p,
Int_t  side 
)
inline
53{ if(GetPlate(p)) if(ValidSide(side)) return GetPlate(p)->GetLayer(side); return 0; }
bool ValidSide(int side) const
Definition: EdbScanSet.h:52

◆ GetPlate()

EdbPlateP * EdbScanSet::GetPlate ( Int_t  p)
inline
56 {
57 if (ePID.Find(p))
58 return eB.GetPlate(ePID.Find(p)->At(0)->Value());
59 else return 0;
60 }

◆ MakeNominalSet()

void EdbScanSet::MakeNominalSet ( EdbID  id,
int  from_plate = 57,
int  to_plate = 1,
float  z0 = 0,
float  dz = -1300,
float  shr = 1,
float  dzbase = 210.,
float  dzem = 45. 
)
76{
77 eID=id;
78 Log(2,"EdbScanSet::MakeNominalSet","%s",eID.AsString());
79 eReferencePlate = from_plate;
80 eB.SetID(id.eBrick);
81 Int_t step= (from_plate>to_plate)?-1:1;
82 Int_t p=from_plate;
83 Float_t z=z0;
84 do {
85 EdbPlateP *plate = new EdbPlateP();
86 plate->SetID(p);
87 plate->SetZ(z);
88 plate->SetPlateLayout(dzbase,dzem,dzem);
89
90// plate->SetZlayer(z,z-dzbase/2-dzem,z+dzbase/2+dzem);
91// plate->GetLayer(0)->SetZlayer(0,-dzbase/2,dzbase/2);
92// plate->GetLayer(1)->SetZlayer(dzbase/2,dzbase/2,dzbase/2+dzem);
93// plate->GetLayer(2)->SetZlayer(-dzbase/2,-dzbase/2-dzem,-dzbase/2);
94
95 plate->GetLayer(1)->SetShrinkage(shr);
96 plate->GetLayer(2)->SetShrinkage(shr);
98 eIDS.Add(new EdbID(id.eBrick,p,id.eMajor,id.eMinor));
99 p += step;
100 z += dz;
101 } while( p!=(to_plate+step) );
102 MakePIDList();
103 if(gEDBDEBUGLEVEL>2) Print();
104}
brick z0
Definition: RecDispMC.C:106
brick dz
Definition: RecDispMC.C:107
char * AsString() const
Definition: EdbID.cxx:26
void Print()
Definition: EdbScanSet.cxx:315
gEDBDEBUGLEVEL
Definition: energy.C:7

◆ MakeParFiles()

Int_t EdbScanSet::MakeParFiles ( Int_t  piece = 0,
const char *  dir = "." 
)

prepare par files for the dataset for "recset-like" processing

108{
110
111 TString name;
112 Log(2,"MakeParFiles","for brick %d with %d plates:\n", eB.ID(), eB.Npl());
113 EdbPlateP *p=0;
114 for(Int_t i=0; i<eB.Npl(); i++) {
115 p = eB.GetPlate(i);
116 name.Form("%s/%2.2d_%3.3d.par",dir,p->ID(),piece);
117 Log(3,"MakeParFiles","%s\n",name.Data());
118 FILE *f = fopen(name.Data(),"w");
119 EdbAffine2D *a=p->GetAffineXY();
120 fprintf(f,"INCLUDE default.par\n");
121 fprintf(f,"ZLAYER 0 %15.3f 0. 0.\n", p->Z());
122 fprintf(f,"AFFXY 0 %15.6f %9.6f %9.6f %9.6f %15.6f %15.6f\n",
123 a->A11(),a->A12(),a->A21(),a->A22(),a->B1(),a->B2());
124 fclose(f);
125 }
126 return eB.Npl();
127}
FILE * f
Definition: RecDispMC.C:150
void a()
Definition: check_aligned.C:59
fclose(pFile)
const char * name
Definition: merge_Energy_SytematicSources_Electron.C:24

◆ MakePIDList()

void EdbScanSet::MakePIDList ( )
179{
180 ePID.Delete();
181 Long_t v[2];
182 for(Int_t i=0; i<eB.Npl(); i++) {
183 v[0]= eB.GetPlate(i)->ID();
184 v[1]= i;
185 ePID.Add(2,v);
186 }
187 ePID.Sort();
188}
void Sort(Int_t upto=kMaxInt)
Definition: TIndexCell.cpp:539
Int_t Add(Int_t narg, Long_t varg[])
Definition: TIndexCell.cpp:602

◆ Print()

void EdbScanSet::Print ( )
316{
317 EdbPlateP *p=0;
318 printf("EdbScanSet: %s\n", eID.AsString() );
319// int npc = ePC.GetEntries();
320// printf("%d couples\n",npc);
321// for(Int_t i=0; i<npc; i++) {
322// p = (EdbPlateP *)(ePC.At(i));
323// p->GetLayer(1)->Print();
324// p->GetLayer(2)->Print();
325// }
326 printf("Brick %d with %d plates:\n", eB.ID(), eB.Npl());
327 for(Int_t i=0; i<eB.Npl(); i++) {
328 p = eB.GetPlate(i);
329 if(!p) continue;
330 EdbAffine2D *a=p->GetAffineXY();
331 if(!a) continue;
332 printf("%3d %12.2f %9.6f %9.6f %9.6f %9.6f %15.6f %15.6f\n",
333 p->ID(), p->Z(), a->A11(),a->A12(),a->A21(),a->A22(),a->B1(),a->B2());
334 }
335 printf("for this brick %d identifiers are defined:\n", eIDS.GetEntries());
336 for(Int_t i=0; i<eIDS.GetEntries(); i++)
337 {
338 EdbID *id = (EdbID*)eIDS.At(i);
339 TString str = id->AsString();
340 EdbPlateP *p = GetPlate(id->ePlate);
341 if(p) str += Form(" shr1/base/shr2 = %.3f / %.1f / %.3f", p->GetLayer(1)->Shr(),p->GetLayer(0)->DZ(), p->GetLayer(2)->Shr());
342 printf("%s\n",str.Data() );
343 }
344}

◆ ReadGeom()

void EdbScanSet::ReadGeom ( const char *  fname)

todo

381{
383
384 char line[256];
385 FILE *f = fopen (file, "r");
386 if(!f) {Log(1,"EdbScanSet::ReadGeom","Can not open file %s", file); return; }
387
388 int b,p,mi,ma;
389 int n=0;
390
391 if(fgets(line, 256, f) != NULL) if( sscanf(line, "%d", &n) != 1 ) return;
392 printf("%d\n",n);
393 if(!n) return;
394
395 ePID.Delete();
396 eB.Clear();
397
398 EdbID *id1=0, *id2=0;
399 float dz=0;
400 for(int i=0; i<n; i++) {
401 if(i==0) {
402 if( fgets(line, 256, f) == NULL ) return;
403 if( sscanf(line, "%d.%d.%d.%d", &b, &p, &mi, &ma) != 4 ) break;
404 id1 = new EdbID(b,p,mi,ma);
405 if(n==1) eIDS.Add(id1);
406 }
407 if( n>1 && i<n-1) {
408 if( fgets(line, 256, f) == NULL ) return;
409 if( sscanf(line, "%f", &dz) != 1 ) break;
410 if( fgets(line, 256, f) == NULL ) return;
411 if( sscanf(line, "%d.%d.%d.%d", &b, &p, &mi, &ma) != 4 ) break;
412 id2 = new EdbID(b,p,mi,ma);
413 eIDS.Add(id2);
414 }
415 }
416 fclose(f);
417}
TFile * file
Definition: write_pvr.C:3
#define NULL
Definition: nidaqmx.h:84

◆ ReadIDS()

Int_t EdbScanSet::ReadIDS ( const char *  file)
348{
349 char line[256];
350 FILE *f = fopen (file, "rt");
351 Int_t b,p,mi,ma;
352 Int_t cnt=0;
353 while(fgets(line, 256, f) != NULL)
354 {
355 if( sscanf(line, "%d, %d, %d, %d", &b, &p, &mi, &ma) != 4 ) break;
356 eIDS.Add(new EdbID(b,p,mi,ma));
357 cnt++;
358 }
359 fclose(f);
360 return cnt;
361}

◆ RemovePlate()

void EdbScanSet::RemovePlate ( int  pid)

Remove both from eB and from eIDS all data corresponding to plate id note that gaps and affine transformations may become inconsistent and had to be updated after this operation

495{
500 EdbID *id = FindPlateID(pid);
501 if(id) {
502 eIDS.Remove(id);
503 MakePIDList();
504 }
505}
void RemovePlate(int pid)
Definition: EdbBrick.cxx:104
int pid[1000]
Definition: m2track.cpp:13

◆ SetAsReferencePlate()

Bool_t EdbScanSet::SetAsReferencePlate ( Int_t  pid)
192{
194 if(!p) return kFALSE;
196 Float_t z = p->Z();
197 EdbAffine2D aff( *(p->GetAffineXY()) );
198 aff.Invert();
199 TransformBrick(aff);
200 ShiftBrickZ(-z);
201 return kTRUE;
202}
Int_t ShiftBrickZ(Float_t z)
Definition: EdbScanSet.cxx:239
Int_t TransformBrick(EdbAffine2D aff)
Definition: EdbScanSet.cxx:216

◆ SetID()

void EdbScanSet::SetID ( EdbID  id)
inline
47{eID=id;}

◆ ShiftBrickZ()

Int_t EdbScanSet::ShiftBrickZ ( Float_t  z)
240{
241 Int_t npl=eB.Npl();
242 if(npl>0) for(Int_t i=0; i<npl; i++) eB.GetPlate(i)->ShiftZ(z);
243 return npl;
244}
void ShiftZ(float dz)
Definition: EdbLayer.cxx:121

◆ TransformBrick() [1/2]

Int_t EdbScanSet::TransformBrick ( EdbAffine2D  aff)
217{
218 Int_t npl=eB.Npl();
219 if(npl>0)
220 for(Int_t i=0; i<npl; i++) eB.GetPlate(i)->GetAffineXY()->Transform(aff);
221 return npl;
222}

◆ TransformBrick() [2/2]

void EdbScanSet::TransformBrick ( EdbScanSet ss)
226{
227 int n = eIDS.GetEntries();
228 for(int i=0; i<n; i++) {
229 EdbID *id = GetID(i); if(!id) continue;
230 EdbPlateP *p = GetPlate(id->ePlate); if(!p) continue;
231 EdbPlateP *p1 = ss.GetPlate(id->ePlate); if(!p1) continue;
232 p->GetAffineXY()->Transform( p1->GetAffineXY() );
233 p->SetZlayer(p1->Z(),p1->Zmin(),p1->Zmax());
234 // ToDo: add also AFFTXTY?
235 }
236}
float Zmin() const
Definition: EdbLayer.h:80
float Zmax() const
Definition: EdbLayer.h:81
EdbID * GetID(Int_t i)
Definition: EdbScanSet.h:46
ss
Definition: energy.C:62

◆ TransformSidesIntoBrickRS()

Int_t EdbScanSet::TransformSidesIntoBrickRS ( )

to put layers 1&2 of each plate into brick reference system - to have absolute transformations for microtracks

206{
208 Int_t npl=eB.Npl();
209 if(npl>0)
210 for(Int_t i=0; i<npl; i++)
212 return npl;
213}
void TransformSidesIntoPlateRS()
Definition: EdbBrick.cxx:31

◆ UpdateBrickWithP2P()

void EdbScanSet::UpdateBrickWithP2P ( EdbLayer la,
int  plate1,
int  plate2 
)

update the brick geometry with aff & dz defined in layer (plate1->plate2) in this version assumed that plate1 and plate2 are consequtive in the brick Note eReferencePlate should be correctly defined

276{
280
281 int pfrom = TMath::Abs(plate1-eReferencePlate)<TMath::Abs(plate2-eReferencePlate)? plate2 : plate1;
282 int step = ((pfrom-eReferencePlate)>0)? 1: -1; //the corrections to be propagated in the opposite direction from ref
283
284 Log(3,"UpdateBrickWithP2P","from,step: %d %d ref: %d Z(%d)=%f",pfrom,step,eReferencePlate,pfrom,GetPlate(pfrom)->Z());
285
286 float dz = (GetDZP2P(plate1,plate2) - la.Zcorr());
287
288 EdbAffine2D *affxy = la.GetAffineXY();
289 if(pfrom==plate2) {
290 dz *= -1.;
291 affxy->Invert();
292 }
293 //Log(1,"UpdateBrickWithP2P","%d->%d dz= %f",plate1,plate2, dz);
294 int ipl=pfrom;
295 for(int i=0; i<eB.Npl(); i++) {
296 EdbPlateP *p = GetPlate(ipl); if(!p) continue;
297 p->SetZlayer( p->Z() + dz, p->Zmin(), p->Zmax() );
298 p->GetAffineXY()->Transform(affxy);
299 ipl+=step;
300 }
301
302 EdbPlateP *plate = GetPlate(plate1);
303 plate->SetShrinkage( GetPlate(plate1)->Shr() * la.Shr() );
304 //plate->GetAffineTXTY()->Print();
305 plate->GetAffineTXTY()->Transform( la.GetAffineTXTY() );
306 //plate->GetAffineTXTY()->Print();
307
308 plate->ApplyCorrectionsLocal(la.Map());
309
310 //Log( 1, "UpdateBrickWithP2P","after: Z(%d)=%f", pfrom, GetPlate(pfrom)->Z() );
311
312}
float Zcorr() const
Definition: EdbLayer.h:90
float Shr() const
Definition: EdbLayer.h:89
EdbAffine2D * GetAffineTXTY()
Definition: EdbLayer.h:120
EdbCorrectionMap & Map()
Definition: EdbLayer.h:72
Float_t GetDZP2P(Int_t p1, Int_t p2)
Definition: EdbScanSet.cxx:263

◆ UpdateGap()

void EdbScanSet::UpdateGap ( float  dz,
int  pid1,
int  pid2 
)

update gap between plates pid1 and pid2. All updates will be done in a way to have the reference plate unchanged

509{
512 EdbLayer la; la.SetZcorr(dz);
513 UpdateBrickWithP2P(la,pid1,pid2);
514}
Definition: EdbLayer.h:39
void SetZcorr(float zcorr)
Definition: EdbLayer.h:106
void UpdateBrickWithP2P(EdbLayer &la, int plate1, int plate2)
Definition: EdbScanSet.cxx:275

◆ UpdateIDS()

void EdbScanSet::UpdateIDS ( int  brick,
int  ma,
int  mi 
)

update all identifiers in dataset with (brik,ma,mi) keeping plateids unchanged

38{
40
41 eID.Set(brick, eID.ePlate, ma, mi);
42 TIter next( &(eIDS) );
43 TObject *obj=0;
44 while((obj = next())) {
45 EdbID *id=(EdbID*)obj;
46 id->Set(brick, id->ePlate, ma, mi);
47 }
48}
BRICK brick
Definition: RecDispMC.C:103
bool Set(const char *id_string)
Definition: EdbID.cxx:19
Int_t ePlate
Definition: EdbID.h:11

◆ ValidSide()

bool EdbScanSet::ValidSide ( int  side) const
inline
52{ if(side>-1&&side<3) return 1; return 0;}

◆ WriteGeom()

void EdbScanSet::WriteGeom ( const char *  fname)
421{
422 FILE *f = 0;
423 if(file) f = fopen(file, "w");
424 if(!f) {Log(1,"EdbScanSet::WriteGeom","Can not open file %s", file); return; }
425 EdbID *id1=0, *id2=0;
426 int n = eIDS.GetEntries();
427 fprintf(f,"%d\n", n);
428 for( Int_t i=0; i<eIDS.GetEntries(); i++ ) {
429 id1 = (EdbID *)eIDS.At(i);
430 if(i==0) fprintf(f,"%d.%d.%d.%d\n", id1->eBrick, id1->ePlate, id1->eMajor, id1->eMinor );
431 if( n>1 && i<n-1) id2 = (EdbID *)eIDS.At(i+1); else continue;
432 fprintf(f,"%f\n", GetDZP2P(id1->ePlate, id2->ePlate));
433 fprintf(f,"%d.%d.%d.%d\n", id2->eBrick, id2->ePlate, id2->eMajor, id2->eMinor );
434 }
435
436 if(file) fclose(f);
437}
Int_t eMinor
Definition: EdbID.h:13
Int_t eMajor
Definition: EdbID.h:12

◆ WriteIDS()

Int_t EdbScanSet::WriteIDS ( const char *  file = 0)
365{
366 FILE *f = 0;
367 if(file) f = fopen(file, "w");
368 EdbID *id;
369 for( Int_t i=0; i<eIDS.GetEntries(); i++ ) {
370 id = (EdbID *)eIDS.At(i);
371 if(file) fprintf(f,"%d.%d.%d.%d\n", id->eBrick, id->ePlate, id->eMajor, id->eMinor );
372 else printf("%d.%d.%d.%d\n", id->eBrick, id->ePlate, id->eMajor, id->eMinor );
373 }
374
375 if(file) fclose(f);
376 return eIDS.GetEntries();
377}

◆ Zlayer()

float EdbScanSet::Zlayer ( int  plate,
int  side 
)
inline
54{ if(GetLayer(plate,side)) return GetLayer(plate,side)->Z() + GetPlate(plate)->Z(); else return 0; }
EdbLayer * GetLayer(Int_t p, Int_t side)
Definition: EdbScanSet.h:53

Member Data Documentation

◆ eB

EdbBrickP EdbScanSet::eB

all layers of the brick defined here

◆ eID

EdbID EdbScanSet::eID

id of the scanset itself

◆ eIDS

TList EdbScanSet::eIDS

list of the identifiers to be processed

◆ ePC

TObjArray EdbScanSet::ePC

Plate Couples. Each couple represented as EdbPlateP where 2 sides in reality corresponds to 2 plates.

◆ ePID

TIndexCell EdbScanSet::ePID

correspondance between index in eB and the plate id

◆ eReferencePlate

Int_t EdbScanSet::eReferencePlate

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