StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
MakeSsdWaferOnGlobal.C
1 class St_db_Maker;
2 St_db_Maker *dbMk = 0;
3 class TTable;
4 TTable *table = 0;
5 //________________________________________________________________________________
6 void Load() {
7  if (gClassTable->GetID("StDbManager") < 0) {
8  // Baseline shared libraries
9  // gSystem->Load("libTable");
10  gSystem->Load("St_base");
11  gSystem->Load("StChain");
12  gSystem->Load("StUtilities");
13  // DB-specific libs
14  // may only need libStDb_Tables.so
15  // but may need St_Tables.so... so I'm using this one
16  // gSystem->Load("libStDb_Tables.so");
17  gSystem->Load("libmysqlclient");
18  gSystem->Load("St_Tables.so");
19  gSystem->Load("StDbLib.so");
20  gSystem->Load("StDbBroker.so");
21  gSystem->Load("St_db_Maker.so");
22  }
23  dbMk = new St_db_Maker("db","$STAR/StarDb","$PWD/StarDb");
24  // dbMk = new St_db_Maker("db","MySQL:StarDb","$STAR/StarDb","$PWD/StarDb");
25  dbMk->SetDebug(1);
26  // dbMk->SetFlavor("ofl+sim");
27  // dbMk->SetFlavor("simu","ssdWafersPosition");
28  dbMk->Init();
29 }
30 //________________________________________________________________________________
31 //void Db(const Char_t *tabNam = "Calibrations/tpc/noiseElim",
32 // void DbS(const Char_t *tabNam =
33 // "Survey/ssd/LadderOnSurvey",Int_t date = 20051101, Int_t time = 0
34 void MakeSsdWaferOnGlobal(Int_t date = 20050316, Int_t time = 195612){
35  if (dbMk == 0) Load();
36  dbMk->SetDebug(2);
37  dbMk->SetDateTime(date,time);
38  // dbMk->SetFlavor("ofl+laserDV","tpcDriftVelocity");
39  // dbMk->SetMaxEntryTime(20040520,0);
40  // to browse 1 database, use this one
41  TDataSet *set = dbMk->GetDataBase("Geometry/ssd");
42  if (! set) return; // Positioning of the SSD:
43  St_Survey *SsdOnGlobal = (St_Survey *) set->Find("SsdOnGlobal");
44  if (! SsdOnGlobal) {cout << "SsdOnGlobal has not been found" << endl; return;}
45  St_Survey *SsdSectorsOnGlobal = (St_Survey *) set->Find("SsdSectorsOnGlobal"); // sectors in the SSD barrel coordinate system
46  if (! SsdSectorsOnGlobal) {cout << "SsdSectorsOnGlobal has not been found" << endl; return;}
47  St_Survey *SsdLaddersOnSectors = (St_Survey *) set->Find("SsdLaddersOnSectors");// ladders in the SSD sector coordinate systems
48  if (! SsdLaddersOnSectors) {cout << "SsdLaddersOnSectors has not been found" << endl; return;}
49  St_Survey *SsdWafersOnLadders = (St_Survey *) set->Find("SsdWafersOnLadders"); // wafers in the SSD ladder coordinate systems
50  // if (! SsdLadderAdjustment) {cout << "SsdLadderAdjustment has not been found" << endl; return;}
51  if (! SsdWafersOnLadders) {cout << "SsdWafersOnLadders has not been found" << endl; return;}
52  St_ssdWafersPosition *ssdWafersPosition = (St_ssdWafersPosition *) set->Find("ssdWafersPosition");
53  if (! ssdWafersPosition) {cout << "ssdWafersPosition has not been found" << endl; return;}
54  Survey_st *OnGlobal = SsdOnGlobal->GetTable(); // SSD and SVT as whole
55  Survey_st *SectorsOnGlobal = SsdSectorsOnGlobal->GetTable(); // sectors in the SSD barrel coordinate system
56  Survey_st *LaddersOnSectors = SsdLaddersOnSectors->GetTable();// ladders in the SSD sector coordinate systems
57  Survey_st *WafersOnLadders = SsdWafersOnLadders->GetTable(); // wafers in the SSD ladder coordinate systems
58  ssdWafersPosition_st *WafersPosition = ssdWafersPosition->GetTable();
59  Int_t NoSectors = SsdSectorsOnGlobal->GetNRows();
60  Int_t NoLadders = SsdLaddersOnSectors->GetNRows();
61  Int_t NoWafers = ssdWafersPosition->GetNRows();
62  St_ssdWafersPosition *ssdwafer = new St_ssdWafersPosition("ssdWafersPosition",NoWafers);
63  TGeoHMatrix GL, WL,LS,SG,LA,WG;
64  GL.SetRotation(&OnGlobal->r00);
65  GL.SetTranslation(&OnGlobal->t0); //cout << "WL\t"; WL.Print();
66 
67  Int_t num = 0;
68  for (Int_t w = 0; w < NoWafers; w++, WafersPosition++) {
69  ssdWafersPosition_st row = *WafersPosition;
70  WafersOnLadders = SsdWafersOnLadders->GetTable();
71  Int_t Id = 0;
72  for (Int_t i = 0; i < NoWafers; i++,WafersOnLadders++) {
73  if (WafersOnLadders->Id != row.id) continue;
74  Id = row.id;
75  break;
76  }
77  if (! Id ) {cout << "Wafer Id\t" << Id << " has not been found" << endl;}
78  Int_t layer = Id/1000;
79  if (layer > 7) layer = 7;
80  Int_t wafer = (Id - 1000*layer)/100;
81  Int_t ladder = Id%100;
82  //cout << "Id\t "<< Id << " " << 100*wafer + ladder + 1000*layer <<endl;
83  WL.SetRotation(&WafersOnLadders->r00);
84  WL.SetTranslation(&WafersOnLadders->t0); //cout << "WL\t"; WL.Print();
85  LaddersOnSectors = SsdLaddersOnSectors->GetTable();
86  Int_t Ladder = 0;
87  Int_t Sector = 0;
88  for (Int_t l = 0; l < NoLadders; l++, LaddersOnSectors++) {
89  //cout << "LaddersOnSectors Id\t" << LaddersOnSectors->Id << endl;
90  Ladder = LaddersOnSectors->Id%100;
91  if (Ladder == ladder) {
92  Sector = LaddersOnSectors->Id/100;
93  LS.SetRotation(&LaddersOnSectors->r00);
94  LS.SetTranslation(&LaddersOnSectors->t0);
95  //cout << "LS\t"; LS.Print();
96  break;
97  }
98  }
99  if (Sector <= 0 || Sector > 4) {cout << "Sector has not been defined" << endl; continue;}
100  SectorsOnGlobal = SsdSectorsOnGlobal->GetTable();
101  Int_t sector = 0;
102  for (Int_t s = 0; s <NoSectors; s++, SectorsOnGlobal++) {
103  //cout << "SectorsOnGlobal Id\t" << SectorsOnGlobal->Id << endl;
104  if (SectorsOnGlobal->Id != Sector) continue;
105  sector = Sector;
106  SG.SetRotation(&SectorsOnGlobal->r00);
107  SG.SetTranslation(&SectorsOnGlobal->t0); //cout << "Sector\t" << Sector << "\tSG\t"; SG.Print();
108  break;
109  }
110  if (! sector) {cout << "Sector\t" << Sector << " has not been found" << endl; continue;}
111  // WG = SG * LS * WL * LA; //cout << "WG\t"; WG.Print();
112  // WG = SG * LS * WL * LA = SG * ( LS * WL * LA * WL**-1 ) *WL
113  WG = GL * SG * LS * WL; //cout << "WG\t"; WG.Print();
114  ssdWafersPosition_st row;
115  row.id = Id;
116  row.id_shape = 2;
117  row.ladder = ladder;
118  row.layer = layer;
119  num++;
120  row.num_chip = (num-1)%16 + 1;
121  // TGeoHMatrix WGInv = WG.Inverse();
122  // Double_t *wgrot = WGInv.GetRotationMatrix();
123  Double_t *r = WG.GetRotationMatrix();
124  row.driftDirection[0] = r[0]; row.normalDirection[0] = r[1]; row.transverseDirection[0] = r[2];
125  row.driftDirection[1] = r[3]; row.normalDirection[1] = r[4]; row.transverseDirection[1] = r[5];
126  row.driftDirection[2] = r[6]; row.normalDirection[2] = r[7]; row.transverseDirection[2] = r[8];
127  Double_t norm;
128  TVector3 d(row.driftDirection); norm = 1/d.Mag(); d *= norm;
129  TVector3 t(row.transverseDirection); norm = 1/t.Mag(); t *= norm;
130  TVector3 n(row.normalDirection);
131  TVector3 c = d.Cross(t);
132  if (c.Dot(n) < 0) c *= -1;
133  d.GetXYZ(row.driftDirection);
134  t.GetXYZ(row.transverseDirection);
135  c.GetXYZ(row.normalDirection);
136 
137  Double_t *wgtr = WG.GetTranslation();
138  // memcpy(row.driftDirection,wgrot, 9*sizeof(Double_t));
139  memcpy(row.centerPosition,wgtr, 3*sizeof(Double_t));
140  ssdwafer->AddAt(&row);
141  }
142  cout << "Create " << Form("ssdWafersPosition.%8i.%06i.C",date,time) << endl;
143  ofstream out;
144  out.open(Form("ssdWafersPosition.%8i.%06i.C",date,time));
145  ssdwafer->SavePrimitive(out,"");
146  out.close();
147 
148 }
Definition: TTable.h:48
virtual TDataSet * Find(const char *path) const
Definition: TDataSet.cxx:362