StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFstFastSimMaker.cxx
1 #include "StFstFastSimMaker.h"
2 
3 #include "St_base/StMessMgr.h"
4 
5 #include "StEvent/StEvent.h"
6 #include "StEvent/StRnDHit.h"
7 #include "StEvent/StFstHit.h"
8 #include "StEvent/StFstHitCollection.h"
9 #include "StEvent/StRnDHitCollection.h"
10 
11 #include "tables/St_g2t_fts_hit_Table.h"
12 #include "tables/St_g2t_track_Table.h"
13 
14 #include "StarClassLibrary/StThreeVectorF.hh"
15 #include "StarGenerator/UTIL/StarRandom.h"
16 
17 #include "TCanvas.h"
18 #include "TCernLib.h"
19 #include "TH2F.h"
20 #include "TLine.h"
21 #include "TRandom3.h"
22 #include "TString.h"
23 #include "TVector2.h"
24 #include "TVector3.h"
25 
26 #include <array>
27 #define _USE_MATH_DEFINES
28 #include <cmath>
29 #include <map>
30 #include <algorithm>
31 
32 
33 // lets not polute the global scope
34 namespace FstGlobal{
35  // converts the position error to a cov mat
36  StMatrixF Hack1to6(const StHit *stHit);
37  // used to convert between uniform and 1 sigma error
38  constexpr float SQRT12 = sqrt(12.0);
39 
40  // Disk segmentation
41  //
42  float RMIN[] = {0.95 * 4.3, 0.95 * 4.3, 0.95 * 4.3, 0.95 * 5.0, 0.95 * 5.0, 0.95 * 5.0};
43  float RMAX[] = {1.05 * 15.0, 1.05 * 25.0, 1.05 * 25.0, 1.05 * 28.0, 1.05 * 28.0, 1.05 * 28.0};
44 
45  //NEXT IS only for disk ARRAY 456 with the radius from 5 to 28.
46  float RSegment[] = {5., 7.875, 10.75, 13.625, 16.5, 19.375, 22.25, 25.125, 28.};
47 
48  // controls some extra output
49  const bool verbose = false;
50 
51  // key type for lookup tables on disk r-phi strip
52  typedef std::tuple<int, int, int> FstKeyTriple;
53 
54 
55 }
56 
57 StFstFastSimMaker::StFstFastSimMaker(const Char_t *name)
58  : StMaker{name},
59  mNumR{8},
60  mNumPHI{128},
61  mNumSEC{12},
62  mInEff{0},
63  mHist{false},
64  mGEANTPassthrough{false},
65  mQAFileName(0),
66  hTrutHitYXDisk(0),
67  hTrutHitRDisk(0),
68  hTrutHitRShower{0},
69  hTrutHitPhiDisk(0),
70  hTrutHitPhiZ(0),
71  hRecoHitYXDisk(0),
72  hRecoHitRDisk(0),
73  hRecoHitPhiDisk(0),
74  hRecoHitPhiZ(0),
75  hGlobalDRDisk(0),
76  hGlobalZ(0),
77  h2GlobalXY(0),
78  h2GlobalSmearedXY(0),
79  h2GlobalDeltaXY(0),
80  h3GlobalDeltaXYDisk(0),
81  h3GlobalDeltaXYR(0) { }
82 
83 int StFstFastSimMaker::Init() {
84 
85  if(mHist){
86  fOut = new TFile(mQAFileName.Data(), "RECREATE");
87  AddHist(hTrutHitYXDisk = new TH3F("hTrutHitYXDisk", "Global hits before segmentation", 151, -75.5, 75.5, 151, -75.5, 75.5, 10, 0, 10));
88  AddHist(hTrutHitRDisk = new TH2F("hTrutHitRDisk", "Global hits before segmentation", 400, 0, 40, 10, 0, 10));
89  AddHist(hTrutHitRShower[0] = new TH2F("hTrutHitRShower_4", "Global hits before segmentation", 400, 0, 40, 20, -10, 10));
90  AddHist(hTrutHitRShower[1] = new TH2F("hTrutHitRShower_5", "Global hits before segmentation", 400, 0, 40, 20, -10, 10));
91  AddHist(hTrutHitRShower[2] = new TH2F("hTrutHitRShower_6", "Global hits before segmentation", 400, 0, 40, 20, -10, 10));
92  AddHist(hTrutHitPhiDisk = new TH2F("hTrutHitPhiDisk", "Global hits before segmentation", 360, 0, 360, 10, 0, 10));
93  AddHist(hTrutHitPhiZ = new TH2F("hTrutHitPhiZ", "Global hits before segmentation", 360, 0, 360, 6000, 0, 600));
94  AddHist(hRecoHitYXDisk = new TH3F("hRecoHitYXDisk", "Global hits after segmentation", 151, -75.5, 75.5, 151, -75.5, 75.5, 10, 0, 10));
95  AddHist(hRecoHitRDisk = new TH2F("hRecoHitRDisk", "Global hits after segmentation", 400, 0, 40, 10, 0, 10));
96  AddHist(hRecoHitPhiDisk = new TH2F("hRecoHitPhiDisk", "Global hits after segmentation", 360, 0, 360, 10, 0, 10));
97  AddHist(hRecoHitPhiZ = new TH2F("hRecoHitPhiZ", "Global hits after segmentation", 360, 0, 360, 6000, 0, 600));
98  AddHist(hGlobalDRDisk = new TH2F("hGlobalDRDisk", "; Reco. r - MC r [cm]; Events;", 1000, -50, 50, 10, 0, 10));
99  AddHist(hGlobalZ = new TH1F("hGlobalZ", "; Z [cm]; Events;", 6000, 0, 600));
100  AddHist(h3GlobalDeltaXYR = new TH3F("h3GlobalDeltaXYR", ";globalDeltaX; globalDeltaY; R", 300, -0.3, 0.3, 300, -3, 3, 100, 0, 30));
101  AddHist(h2GlobalXY = new TH2F("h2GlobalXY", ";globalX; globalY", 1510, -75.5, 75.5, 1510, -75.5, 75.5));
102  AddHist(h2GlobalSmearedXY = new TH2F("h2GlobalSmearedXY", ";globalSmearedX; globalSmearedY", 1510, -75.5, 75.5, 1510, -75.5, 75.5));
103  AddHist(h2GlobalDeltaXY = new TH2F("h2GlobalDeltaXY", ";globalDeltaX; globalDeltaY", 151, -75.5, 75.5, 151, -75.5, 75.5));
104  AddHist(h3GlobalDeltaXYDisk = new TH3F("h3GlobalDeltaXYDisk", ";globalDeltaX; globalDeltaY; Disk", 151, -75.5, 75.5, 151, -75.5, 75.5, 10, 0, 10));
105  }
106  return StMaker::Init();
107 }
108 
109 void StFstFastSimMaker::SetDisk(const int i, const float rmn, const float rmx) {
110  FstGlobal::RMIN[i] = rmn;
111  FstGlobal::RMAX[i] = rmx;
112 }
113 
115  LOG_DEBUG << "StFstFastSimMaker::Make" << endm;
116 
117  // Get the existing StEvent, or add one if it doesn't exist.
118  StEvent *event = static_cast<StEvent *>(GetDataSet("StEvent"));
119  if (!event) {
120  event = new StEvent;
121  AddData(event);
122  LOG_DEBUG << "Creating StEvent" << endm;
123  }
124 
125  if (0 == event->rndHitCollection()) {
126  event->setRnDHitCollection(new StRnDHitCollection());
127  LOG_DEBUG << "Creating StRnDHitCollection for FTS" << endm;
128  }
129 
130  // Get pointer to an existing StFstHitCollection if any
131  StFstHitCollection *fstHitCollection = event->fstHitCollection();
132  // If no fst hit collection, create one
133  if (!fstHitCollection) {
134  fstHitCollection = new StFstHitCollection();
135  event->setFstHitCollection(fstHitCollection);
136  LOG_DEBUG << "Make() - Added new StFstHitCollection to this StEvent" << endm;
137  }
138 
139  // Digitize GEANT FTS hits
140  FillSilicon(event);
141 
142  return kStOk;
143 }
144 
145 /* Fill an event with StFtsHits. */
146 /* This should fill StFtsStrip for realistic simulator and let clustering fill StFtsHit */
147 /* For now skipping StFtsStrip and clustering, and fill StFtsHits directly here*/
148 
149 void StFstFastSimMaker::FillSilicon(StEvent *event) {
150 
151  StRnDHitCollection *fsicollection = event->rndHitCollection();
152 
153  const int MAXR = mNumR;
154  const int MAXPHI = mNumPHI * mNumSEC;
155 
156  if ( mGEANTPassthrough ){
157  LOG_INFO << "FST Hits using GEANT xyz directly (no raster etc.)" << endm;
158  }
159 
160  // maps for hit and energy for each disk's r-phi strip
161  std::map< FstGlobal::FstKeyTriple, StRnDHit* > hitMap;
162  std::map< FstGlobal::FstKeyTriple, double > energySum;
163  std::map< FstGlobal::FstKeyTriple, double > energyMax;
164 
165  // Read the g2t table
166  St_g2t_fts_hit *hitTable = static_cast<St_g2t_fts_hit *>(GetDataSet("g2t_fsi_hit"));
167  if (!hitTable) {
168  LOG_INFO << "g2t_fsi_hit table is empty" << endm;
169  return; // Nothing to do
170  }
171 
172  const Int_t nHits = hitTable->GetNRows();
173  LOG_DEBUG << "g2t_fsi_hit table has " << nHits << " hits" << endm;
174  const g2t_fts_hit_st *hit = hitTable->GetTable();
175 
176  StPtrVecRnDHit hits;
177 
178  // track table
179  St_g2t_track *trkTable = static_cast<St_g2t_track *>(GetDataSet("g2t_track"));
180  if (!trkTable) {
181  LOG_INFO << "g2t_track table is empty" << endm;
182  return; // Nothing to do
183  }
184 
185  const Int_t nTrks = trkTable->GetNRows();
186 
187  LOG_DEBUG << "g2t_track table has " << nTrks << " tracks" << endm;
188 
189  const g2t_track_st *trk = trkTable->GetTable();
190 
191 
192  int count = 0;
193  for (Int_t i = 0; i < nHits; ++i) {
194 
195  hit = (g2t_fts_hit_st *)hitTable->At(i);
196 
197  // skip bad hits
198  if (!hit)
199  continue;
200 
201 
202  int volume_id = hit->volume_id;
203  LOG_DEBUG << "volume_id = " << volume_id << endm;
204  int disk = volume_id / 1000; // disk id
205  int wedge = (volume_id % 1000) / 10; // wedge id
206  int sensor = volume_id % 10; // sensor id
207 
208  // used as an index for various arrays
209  size_t disk_index = disk - 1;
210 
211  LOG_DEBUG << "disk = " << disk << ", wedge = " << wedge << ", sensor = " << sensor << endm;
212 
213  // skip non-FST hits
214  if ( disk > 6 ) continue;
215 
216 
217  double energy = hit->de;
218  int track = hit->track_p;
219 
220  trk = (g2t_track_st *)trkTable->At(track);
221  int isShower = false;
222  if (trk)
223  isShower = trk->is_shower;
224 
225  // This z-offset is used to shift the hits
226  // to the center of the FST where the tracking planes are defined
227  const double z_delta = 1.755;
228  // hit coordinates
229  double x = hit->x[0];
230  double y = hit->x[1];
231  double z = hit->x[2] + z_delta;
232 
233  if (z > 200)
234  continue; // skip large disks
235 
236  double r = sqrt(x * x + y * y);
237  double p = atan2(y, x);
238 
239  // wrap an angle between 0 and 2pi
240  auto wrapAngle = [&]( double angle ) {
241  angle = fmod( angle, 2.0 * M_PI );
242  if ( angle < 0 )
243  angle += 2.0 * M_PI;
244  return angle;
245  };
246 
247  p = wrapAngle( p );
248 
249  LOG_DEBUG << "r = " << r << " p=" << p << endm;
250  LOG_DEBUG << "RMIN = " << FstGlobal::RMIN[disk_index] << " RMAX= " << FstGlobal::RMAX[disk_index] << endm;
251 
252  // Cuts made on the r value to ensure it is within limits
253  if (r < FstGlobal::RMIN[disk_index] || r > FstGlobal::RMAX[disk_index])
254  continue;
255 
256  // Strip numbers on rastered value
257  int r_index = floor(MAXR * (r - FstGlobal::RMIN[disk_index]) / (FstGlobal::RMAX[disk_index] - FstGlobal::RMIN[disk_index]));
258  for (int ii = 0; ii < MAXR; ii++)
259  if (r > FstGlobal::RSegment[ii] && r <= FstGlobal::RSegment[ii + 1])
260  r_index = ii;
261 
262  // Phi number
263  int phi_index = int(MAXPHI * p / 2.0 / M_PI);
264 
265  if (r_index >= 8)
266  continue;
267 
268  if (MAXR)
269  assert(r_index < MAXR);
270  if (MAXPHI)
271  assert(phi_index < MAXPHI);
272 
273  StRnDHit *fsihit = nullptr;
274 
275  // key of this disk's r & phi strip
276  auto threeKey = std::tie( disk_index, r_index, phi_index );
277 
278  if (hitMap.count( threeKey ) == 0) { // New hit
279 
280  if (FstGlobal::verbose){
281  LOG_INFO << Form("NEW d=%1d xyz=%8.4f %8.4f %8.4f r=%8.4f phi=%8.4f iR=%2d iPhi=%4d dE=%8.4f[MeV] truth=%d",
282  disk, x, y, z, r, p, r_index, phi_index, energy * 1000.0, track)
283  << endm;
284  }
285 
286  count++;
287  fsihit = new StRnDHit();
288  fsihit->setDetectorId(kFtsId);
289  fsihit->setLayer(disk);
290  fsihit->setLadder(wedge);
291  fsihit->setWafer(sensor);
292  //
293  // Set position and position error based on radius-constant bins
294  //
295  double p0 = (phi_index + 0.5) * 2.0 * M_PI / double(MAXPHI);
296  double dp = 2.0 * M_PI / double(MAXPHI) / FstGlobal::SQRT12;
297 
298  // ONLY valid for the disk array 456, no difference for each disk
299  double r0 = (FstGlobal::RSegment[r_index] + FstGlobal::RSegment[r_index + 1]) * 0.5;
300  double dr = FstGlobal::RSegment[r_index + 1] - FstGlobal::RSegment[r_index];
301 
302  double x0 = r0 * cos(p0);
303  double y0 = r0 * sin(p0);
304  assert(TMath::Abs(x0) + TMath::Abs(y0) > 0);
305  double dz = 0.03 / FstGlobal::SQRT12;
306  double er = dr / FstGlobal::SQRT12;
307  fsihit->setPosition(StThreeVectorF(x0, y0, z));
308  // pass the GEANT hits through without modification
309  if ( mGEANTPassthrough ){
310  fsihit->setPosition(StThreeVectorF(x, y, z));
311  }
312 
313  fsihit->setPositionError(StThreeVectorF(er, dp, dz));
314  // set covariance matrix
315  fsihit->setErrorMatrix(&FstGlobal::Hack1to6(fsihit)[0][0]);
316 
317  fsihit->setCharge(energy);
318  fsihit->setIdTruth(track, 100);
319  hits.push_back(fsihit);
320 
321  hitMap[ threeKey ] = fsihit;
322  energySum[ threeKey ] = energy;
323  energyMax[ threeKey ] = energy;
324 
325  if (FstGlobal::verbose){
326  LOG_INFO << Form("NEW d=%1d xyz=%8.4f %8.4f %8.4f ", disk, x, y, z) << endm;
327  LOG_INFO << Form("smeared xyz=%8.4f %8.4f %8.4f ", fsihit->position().x(), fsihit->position().y(), fsihit->position().z()) << endm;
328  }
329 
330  if(mHist){
331  TVector2 hitpos_mc(x, y);
332  TVector2 hitpos_rc(fsihit->position().x(), fsihit->position().y());
333 
334  hTrutHitYXDisk->Fill(x, y, disk);
335  hTrutHitRDisk->Fill(hitpos_mc.Mod(), disk);
336 
337  if (disk == 4)
338  hTrutHitRShower[0]->Fill(hitpos_mc.Mod(), isShower);
339  if (disk == 5)
340  hTrutHitRShower[1]->Fill(hitpos_mc.Mod(), isShower);
341  if (disk == 6)
342  hTrutHitRShower[2]->Fill(hitpos_mc.Mod(), isShower);
343 
344  hTrutHitPhiDisk->Fill(hitpos_mc.Phi() * 180.0 / TMath::Pi(), disk);
345  hTrutHitPhiZ->Fill(hitpos_mc.Phi() * 180.0 / TMath::Pi(), z);
346  hRecoHitYXDisk->Fill(fsihit->position().x(), fsihit->position().y(), disk);
347  hRecoHitRDisk->Fill(hitpos_rc.Mod(), disk);
348  hRecoHitPhiDisk->Fill(hitpos_rc.Phi() * 180.0 / TMath::Pi(), disk);
349  hRecoHitPhiZ->Fill(hitpos_rc.Phi() * 180.0 / TMath::Pi(), z);
350  hGlobalDRDisk->Fill(hitpos_rc.Mod() - hitpos_mc.Mod(), disk);
351  hGlobalZ->Fill(fsihit->position().z());
352 
353  // cout << "CHECK : " << fsihit->position().x()-x << " ||| "<< fsihit->position().y()-y << endl;
354  h2GlobalXY->Fill(x, y);
355  h2GlobalSmearedXY->Fill(fsihit->position().x(), fsihit->position().y());
356  h2GlobalDeltaXY->Fill(fsihit->position().x() - x, fsihit->position().y() - y);
357  h3GlobalDeltaXYDisk->Fill(fsihit->position().x() - x, fsihit->position().y() - y, disk);
358 
359  h3GlobalDeltaXYR->Fill(fsihit->position().x() - x, fsihit->position().y() - y, sqrt(pow(fsihit->position().x(), 2) + pow(fsihit->position().y(), 2)));
360  }
361  }
362  else { // Hit on this strip already exists, adding energy to old hit
363  // get hit from the map
364  fsihit = hitMap[ threeKey ] ;
365  fsihit->setCharge(fsihit->charge() + energy);
366 
367  // Add energy to running sum
368  energySum[ threeKey ] = energySum[ threeKey ] + energy;
369 
370  if (energy> energyMax[ threeKey ])
371  energyMax[ threeKey ] = energy;
372 
373  // keep idtruth but dilute it...
374  track = fsihit->idTruth();
375 
376  fsihit->setIdTruth(track, 100 * (energyMax[ threeKey ] / energySum[ threeKey ]));
377  }
378  } // loop on hits
379  int nfsihit = hits.size();
380 
382  LOG_INFO << "FST Fast simulator is using mInEff = " << mInEff << endm;
383  // NOW run back through the hits and add them if they pass an efficiency roll
384  for (int i = 0; i < nfsihit; i++) {
385  double rnd_save = rand.flat();
386  if (rnd_save > mInEff || mGEANTPassthrough){
387  fsicollection->addHit(hits[i]);
388 
389  StFstHit *fHit = new StFstHit(hits[i]->position(), hits[i]->positionError(), 0, hits[i]->charge(), 0);
390  fHit->setIdTruth(hits[i]->idTruth());
391  float r = sqrt(hits[i]->position().x() * hits[i]->position().x() + hits[i]->position().y() * hits[i]->position().y());
392  float phi = atan2(hits[i]->position().y(), hits[i]->position().x());
393  fHit->setLocalPosition( r, phi, hits[i]->position().z() );
394  fHit->setDiskWedgeSensor(hits[i]->layer(), hits[i]->ladder(), hits[i]->wafer());
395  event->fstHitCollection()->addHit(fHit);
396  } else {
397  }
398  }
399  if (FstGlobal::verbose) {
400  LOG_DEBUG << Form("Found %d/%d g2t hits in %d cells, created %d hits with ADC>0 and put %d into StFstHitCollection", count, nHits, nfsihit, fsicollection->numberOfHits(), event->fstHitCollection()->numberOfHits()) << endm;
401  }
402 
403 }
404 //
405 
407  if(mHist){
408  fOut->cd();
409  hTrutHitYXDisk->Write();
410  hTrutHitRDisk->Write();
411  hTrutHitRShower[0]->Write();
412  hTrutHitRShower[1]->Write();
413  hTrutHitRShower[2]->Write();
414  hTrutHitPhiDisk->Write();
415  hTrutHitPhiZ->Write();
416  hRecoHitYXDisk->Write();
417  hRecoHitRDisk->Write();
418  hRecoHitPhiDisk->Write();
419  hRecoHitPhiZ->Write();
420  hGlobalDRDisk->Write();
421  hGlobalZ->Write();
422  h3GlobalDeltaXYR->Write();
423  h2GlobalXY->Write();
424  h2GlobalSmearedXY->Write();
425  h2GlobalDeltaXY->Write();
426  h3GlobalDeltaXYDisk->Write();
427  fOut->Close();
428  }
429  return kStOK;
430 }
431 
432 //_____________________________________________________________________________
433 StMatrixF FstGlobal::Hack1to6(const StHit *stHit) {
434  // X = R*cos(Fi), Y=R*sin(Fi), Z = z
435  // dX/dR = ( cos(Fi) ,sin(Fi),0)
436  // dX/dFi = (-R*sin(Fi), R*cos(Fi),0)
437  // dX/dZ = ( 0, 0,1)
438 
439  auto hiPos = stHit->position();
440  auto hiErr = stHit->positionError();
441  double Rxy = sqrt(hiPos[0] * hiPos[0] + hiPos[1] * hiPos[1]);
442  double cosFi = hiPos[0] / Rxy;
443  double sinFi = hiPos[1] / Rxy;
444  double T[3][3] = {{cosFi, -Rxy * sinFi, 0}, {sinFi, Rxy * cosFi, 0}, {0, 0, 1}};
445  double Ginp[6] = {hiErr[0] * hiErr[0], 0, hiErr[1] * hiErr[1], 0, 0, hiErr[2] * hiErr[2]};
446  double Gout[6];
447 
448  TCL::trasat(T[0], Ginp, Gout, 3, 3);
449  StMatrixF mtxF(3, 3);
450 
451  for (int i = 0, li = 0; i < 3; li += ++i) {
452  for (int j = 0; j <= i; j++) {
453  mtxF[i][j] = Gout[li + j];
454  mtxF[j][i] = mtxF[i][j];
455  }
456  }
457 
458  return mtxF;
459 }
static StarRandom & Instance()
Obtain the single instance of the random number generator.
Definition: StarRandom.cxx:87
Definition: StHit.h:125
virtual void AddData(TDataSet *data, const char *dir=".data")
User methods.
Definition: StMaker.cxx:332
A class for providing random number generation.
Definition: StarRandom.h:30
Double_t flat() const
Return a random number uniformly distributed between 0 and 1.
Definition: StarRandom.cxx:68
static float * trasat(const float *a, const float *s, float *r, int m, int n)
Definition: TCernLib.cxx:540
Definition: Stypes.h:40
void SetDisk(const int i, const float rmn, const float rmx)
Set min/max active radii for each disk.
Definition: Stypes.h:41