StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
MakeSvtWaferOnGlobal.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("St_base");
10  gSystem->Load("StChain");
11  gSystem->Load("StUtilities");
12  // DB-specific libs
13  // may only need libStDb_Tables.so
14  // but may need St_Tables.so... so I'm using this one
15  // gSystem->Load("libStDb_Tables.so");
16  gSystem->Load("St_Tables.so");
17  gSystem->Load("StDbLib.so");
18  gSystem->Load("StDbBroker.so");
19  gSystem->Load("St_db_Maker.so");
20  }
21  dbMk = new St_db_Maker("db","$STAR/StarDb","$PWD/StarDb");
22  // dbMk = new St_db_Maker("db","MySQL:StarDb","$STAR/StarDb","$PWD/StarDb");
23  dbMk->SetDebug(1);
24  // dbMk->SetFlavor("ofl+sim");
25  dbMk->SetFlavor("simu","svtWafersPosition");
26  dbMk->Init();
27 }
28 //________________________________________________________________________________
29 //void Db(const Char_t *tabNam = "Calibrations/tpc/noiseElim",
30 // void DbS(const Char_t *tabNam =
31 // "Survey/svt/LadderOnSurvey",Int_t date = 20051101, Int_t time = 0
32 void MakeSvtWaferOnGlobal(Int_t date = 20050101, Int_t time = 65 ){
33  TGeoHMatrix GL, WL,LSU,LSH,SHG,WG;
34  if (dbMk == 0) Load();
35  dbMk->SetDebug(2);
36  dbMk->SetDateTime(date,time);
37  // dbMk->SetFlavor("ofl+laserDV","tpcDriftVelocity");
38  // dbMk->SetMaxEntryTime(20040520,0);
39  // to browse 1 database, use this one
40  TDataSet *set = dbMk->GetDataBase("Geometry/ssd");
41  if (! set) return; // Positioning of the SSD:
42  St_Survey *SsdOnGlobal = (St_Survey *) set->Find("SsdOnGlobal");
43  if (! SsdOnGlobal) {cout << "SsdOnGlobal has not been found" << endl; return;}
44  Survey_st *OnGlobal = SsdOnGlobal->GetTable(); // SSD and SVT as whole
45  GL.SetRotation(&OnGlobal->r00);
46  GL.SetTranslation(&OnGlobal->t0); //cout << "WL\t"; WL.Print();
47  set = dbMk->GetDataBase("Geometry/svt");
48  St_Survey *WaferOnLadder = (St_Survey *) set->Find("WaferOnLadder");
49  St_Survey *LadderOnSurvey = (St_Survey *) set->Find("LadderOnSurvey");
50  St_Survey *LadderOnShell = (St_Survey *) set->Find("LadderOnShell");
51  St_Survey *ShellOnGlobal = (St_Survey *) set->Find("ShellOnGlobal");
52  Int_t NW = WaferOnLadder->GetNRows();
53  Int_t NL = LadderOnSurvey->GetNRows();
54  Survey_st *waferOnLadder = WaferOnLadder->GetTable();
55  Survey_st *ladderOnSurvey = LadderOnSurvey->GetTable();
56  Survey_st *ladderOnShell = LadderOnShell->GetTable();
57  Survey_st *shellOnGlobal0 = ShellOnGlobal->GetTable(0);
58  Survey_st *shellOnGlobal1 = ShellOnGlobal->GetTable(1);
59  St_svtWafersPosition *svtwafer = new St_svtWafersPosition("svtWafersPosition",216);
60  svtWafersPosition_st row;
61  for (Int_t i = 0; i < NW; i++, waferOnLadder++)
62  {
63  Int_t Idw = waferOnLadder->Id;
64  WL.SetRotation(&waferOnLadder->r00);
65  WL.SetTranslation(&waferOnLadder->t0);
66  // if (i==0) WL.Print();
67  Int_t wshell = 0;
68  Int_t wbarrel = Idw/1000;
69  Int_t wwafer = (Idw - 1000*wbarrel)/100;
70  Int_t wladder = Idw%100;
71  Int_t wlayer = 2*wbarrel + wladder%2 - 1;
72  // cout << waferOnLadder->Id << " "<< Idw<< " " << 100*wwafer + wladder + 1000*wlayer <<endl;
73  for ( Int_t j = 0; j < NL; j++, ladderOnSurvey++, ladderOnShell++)
74  {
75  Int_t Idl = ladderOnSurvey->Id;
76  Int_t lbarrel = Idl/1000;
77  Int_t lladder = Idl%100;
78  if( wladder == lladder )
79  {
80  LSU.SetRotation(&ladderOnSurvey->r00);
81  LSU.SetTranslation(&ladderOnSurvey->t0);
82  LSH.SetRotation(&ladderOnShell->r00);
83  LSH.SetTranslation(&ladderOnShell->t0);
84  if( (wbarrel == 1 && wladder <= 4) || (wbarrel == 2 && wladder <= 6) || (wbarrel == 3 && wladder <= 8) )
85  {
86  SHG.SetRotation(&shellOnGlobal0->r00);
87  SHG.SetTranslation(&shellOnGlobal0->t0);
88  }else
89  {
90  SHG.SetRotation(&shellOnGlobal1->r00);
91  SHG.SetTranslation(&shellOnGlobal1->t0);
92  }
93  // SsdOnGlobal * ShellOnGlobal * LadderOnShell * LadderOnSurvey * WaferOnLadder
94  WG = GL * SHG * LSH * LSU * WL; // WG.Print();
95  // TGeoHMatrix WGInv = WG.Inverse();
96  Double_t *r = WG.GetRotationMatrix();
97  Int_t fail = 0;
98  for (int l = 0; l < 9; l++) {
99  if (TMath::Abs(r[l]) >= 1.000001) fail++;
100  }
101  if (fail) {
102  cout << "===============" << waferOnLadder->Id << " "<< Idw << " " << 100*wwafer + wladder + 1000*wlayer <<endl;
103  cout << "WG\t"; WG.Print();
104  // cout << "SHG\t"; SHG.Print();
105  // cout << "LSH\t"; LSH.Print();
106  // cout << "LSU\t"; LSU.Print();
107  // cout << "WL\t"; WL.Print();
108  }
109  row.driftDirection[0] = r[0]; row.normalDirection[0] = r[1]; row.transverseDirection[0] = r[2];
110  row.driftDirection[1] = r[3]; row.normalDirection[1] = r[4]; row.transverseDirection[1] = r[5];
111  row.driftDirection[2] = r[6]; row.normalDirection[2] = r[7]; row.transverseDirection[2] = r[8];
112  Double_t norm;
113  TVector3 d(row.driftDirection); norm = 1/d.Mag(); d *= norm;
114  TVector3 t(row.transverseDirection); norm = 1/t.Mag(); t *= norm;
115  TVector3 n(row.normalDirection);
116  TVector3 c = d.Cross(t);
117  if (c.Dot(n) < 0) c *= -1;
118  d.GetXYZ(row.driftDirection);
119  t.GetXYZ(row.transverseDirection);
120  c.GetXYZ(row.normalDirection);
121 
122  row.ID = 100*wwafer + wladder + 1000*wlayer;
123  Double_t *wgtr = WG.GetTranslation();
124  memcpy(row.centerPosition,wgtr, 3*sizeof(Double_t));
125  svtwafer->AddAt(&row);
126  break;
127 
128  }
129  }
130  }
131  ofstream out;
132  out.open(Form("svtWafersPosition.%8i.%06i.C",date,time));
133  svtwafer->SavePrimitive(out,"");
134  out.close();
135 
136 }
137 
138 // // table = (TTable *) set->FindByName(gSystem->BaseName(tabNam));
139 // if (table) {
140 // TDatime t[2];
141 // dbMk->GetValidity(table,t);
142 // cout << "==============================================" << endl;
143 // Int_t Nrows = table->GetNRows();
144 // cout << "Found table " << table->GetName() << " with NRows = " << Nrows << " in db" << endl;
145 // cout << "Validity:" << t[0].GetDate() << "/" << t[0].GetTime()
146 // << " ----- " << t[1].GetDate() << "/" << t[1].GetTime() << endl;
147 // if (Nrows > 10) Nrows = 10;
148 // table->Print(0,Nrows);
149 // cout << "==============================================" << endl;
150 // // TString name(gSystem->BaseName(tabNam));
151 // // name += Form(".%06i.%06i.old.root",t[0].GetDate(),t[0].GetTime());
152 // // TFile *f = new TFile(name.Data(),"RECREATE");
153 // // table->Write();
154 // // delete f;
155 // }
156 // else cout << "Table:" << tabNam << " has not been found" << endl;
Definition: TTable.h:48
virtual TDataSet * Find(const char *path) const
Definition: TDataSet.cxx:362