CbmRoot
Loading...
Searching...
No Matches
CbmMuchDigitizeGem.cxx
Go to the documentation of this file.
1/* Copyright (C) 2009-2021 St. Petersburg Polytechnic University, St. Petersburg
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Apar Agarwal, Vikas Singhal, Ekata Nandy, Volker Friese, Evgeny Kryshen [committer] */
4
33
34// Includes from MUCH
35#include "CbmMuchDigitizeGem.h"
36
37#include "CbmMuchDigi.h"
38#include "CbmMuchModuleGem.h"
41#include "CbmMuchPad.h"
42#include "CbmMuchPadRadial.h"
44#include "CbmMuchPoint.h"
46#include "CbmMuchRecoDefs.h"
47#include "CbmMuchSector.h"
48#include "CbmMuchSectorRadial.h"
50#include "CbmMuchStation.h"
51
52// Includes from base
53#include "FairEventHeader.h"
54#include "FairMCEventHeader.h"
55#include "FairMCPoint.h"
56#include "FairRootManager.h"
57#include "FairRunAna.h"
58#include "FairRunSim.h"
59
60#include <Logger.h>
61
62// Includes from Cbm
63#include "CbmMCTrack.h"
64#include "TCanvas.h"
65#include "TChain.h"
66#include "TDatabasePDG.h"
67#include "TFile.h"
68#include "TH1D.h"
69#include "TObjArray.h"
70#include "TRandom.h"
71
72#include <TGeoManager.h>
73
74#include <cassert>
75#include <cstring>
76#include <iomanip>
77#include <iostream>
78#include <sstream>
79#include <vector>
80
81using std::cout;
82using std::endl;
83using std::fixed;
84using std::map;
85using std::right;
86using std::setprecision;
87using std::setw;
88using std::string;
89
90// ----- Default constructor -------------------------------------------
92 : CbmDigitize<CbmMuchDigi>("MuchDigitizeGem")
93 ,
94 // fgDeltaResponse(),
95 fAlgorithm(1)
96 , fGeoScheme(CbmMuchGeoScheme::Instance())
97 , fDigiFile()
98 , fPoints(NULL)
99 , fMCTracks(NULL)
100 ,
101 // fDigis(NULL),
102 // fDigiMatches(NULL),
103 fNFailed(0)
104 , fNOutside(0)
105 , fNMulti(0)
106 , fFlag(0)
107 , fNADCChannels(-1)
108 , fQMax(-1)
109 , fQThreshold(-1)
110 , fMeanNoise(1500)
111 , fSpotRadius(-1)
112 , fMeanGasGain(-1)
113 , fDTime(-1)
114 , fDeadPadsFrac(-1)
115 , fTimer()
116 , fMcChain(NULL)
117 , fDeadTime(400)
118 , fDriftVelocity(-1)
119 ,
120 // fPeakingTime(20),
121 // fRemainderTime(40),
123 , fNTimeBins(200)
124 , fTOT(0)
125 , fTotalDriftTime(-1)
126 ,
127 // fTotalDriftTime(0.3/fDriftVelocity*10000), // 40 ns
128 fSigma()
129 , fMPV()
130 , fIsLight(1)
131 ,
132 // fIsLight = 1 (default) Store Light CbmMuchDigiMatch in
133 // output branch, fIsLight = 0 Create Heavy CbmMuchDigiMatch
134 // with fSignalShape info.
136 , fTimeDigiFirst(-1)
137 , fTimeDigiLast(-1)
138 , fNofPoints(0)
139 , fNofSignals(0)
140 , fNofDigis(0)
141 , fNofEvents(0)
142 , fNofPointsTot(0.)
143 , fNofSignalsTot(0.)
144 , fNofDigisTot(0.)
145 , fTimeTot()
146 ,
147 // hPriElAfterDriftpath(NULL),
148 fPerPadNoiseRate(10e-9)
150 , fNoiseCharge(nullptr)
152{
153 fSigma[0] = new TF1("sigma_e", "pol6", -5, 10);
154 fSigma[0]->SetParameters(sigma_e);
155
156 fSigma[1] = new TF1("sigma_mu", "pol6", -5, 10);
157 fSigma[1]->SetParameters(sigma_mu);
158
159 fSigma[2] = new TF1("sigma_p", "pol6", -5, 10);
160 fSigma[2]->SetParameters(sigma_p);
161
162 fMPV[0] = new TF1("mpv_e", "pol6", -5, 10);
163 fMPV[0]->SetParameters(mpv_e);
164
165 fMPV[1] = new TF1("mpv_mu", "pol6", -5, 10);
166 fMPV[1]->SetParameters(mpv_mu);
167
168 fMPV[2] = new TF1("mpv_p", "pol6", -5, 10);
169 fMPV[2]->SetParameters(mpv_p);
170 Reset();
171 fNoiseCharge = new TF1("Noise Charge", "TMath::Gaus(x, [0], [1])", fQThreshold,
172 fQMax / 10); // noise function to calculate charge for noise hit.
173 // mean=fQThreashold(10000),fQMax=500000
174}
175// -------------------------------------------------------------------------
176
177// -------------------------------------------------------------------------
178CbmMuchDigitizeGem::CbmMuchDigitizeGem(const char* digiFileName, Int_t flag)
179 : CbmDigitize<CbmMuchDigi>("MuchDigitizeGem")
180 ,
181 // fgDeltaResponse(),
182 fAlgorithm(1)
183 , fGeoScheme(CbmMuchGeoScheme::Instance())
184 , fDigiFile(digiFileName)
185 , fPoints(NULL)
186 , fMCTracks(NULL)
187 ,
188 // fDigis(NULL),
189 // fDigiMatches(NULL),
190 fNFailed(0)
191 , fNOutside(0)
192 , fNMulti(0)
193 , fFlag(flag)
194 , fNADCChannels(-1)
195 , fQMax(-1)
196 , fQThreshold(-1)
197 , fMeanNoise(1500)
198 , fSpotRadius(-1)
199 , fMeanGasGain(-1)
200 , fDTime(-1)
201 , fDeadPadsFrac(-1)
202 , fTimer()
203 , fMcChain(NULL)
204 , fDeadTime(400)
205 , fDriftVelocity(-1)
206 ,
207 // fPeakingTime(20),
208 // fRemainderTime(40),
210 , fNTimeBins(200)
211 , fTOT(0)
212 , fTotalDriftTime(-1)
213 ,
214 // fTotalDriftTime(0.3/fDriftVelocity*10000), // 40 ns
215 fSigma()
216 , fMPV()
217 , fIsLight(1)
218 ,
219 // fIsLight = 1 (default) Store Light CbmMuchDigiMatch in
220 // output branch, fIsLight = 0 Create Heavy CbmMuchDigiMatch
221 // with fSignalShape info.
223 , fTimeDigiFirst(-1)
224 , fTimeDigiLast(-1)
225 , fNofPoints(0)
226 , fNofSignals(0)
227 , fNofDigis(0)
228 , fNofEvents(0)
229 , fNofPointsTot(0.)
230 , fNofSignalsTot(0.)
231 , fNofDigisTot(0.)
232 ,
233 // hPriElAfterDriftpath(NULL),
234 fTimeTot()
235 , fPerPadNoiseRate(10e-9)
237 , fNoiseCharge(nullptr)
239{
240 fSigma[0] = new TF1("sigma_e", "pol6", -5, 10);
241 fSigma[0]->SetParameters(sigma_e);
242
243 fSigma[1] = new TF1("sigma_mu", "pol6", -5, 10);
244 fSigma[1]->SetParameters(sigma_mu);
245
246 fSigma[2] = new TF1("sigma_p", "pol6", -5, 10);
247 fSigma[2]->SetParameters(sigma_p);
248
249 fMPV[0] = new TF1("mpv_e", "pol6", -5, 10);
250 fMPV[0]->SetParameters(mpv_e);
251
252 fMPV[1] = new TF1("mpv_mu", "pol6", -5, 10);
253 fMPV[1]->SetParameters(mpv_mu);
254
255 fMPV[2] = new TF1("mpv_p", "pol6", -5, 10);
256 fMPV[2]->SetParameters(mpv_p);
257 Reset();
258 fNoiseCharge = new TF1("Noise Charge", "TMath::Gaus(x, [0], [1])", fQThreshold,
259 fQMax / 10); // noise function to calculate charge for noise hit.
260 // mean=fQThreashold(10000),fQMax=500000
261}
262// -------------------------------------------------------------------------
263
264// ----- Destructor ----------------------------------------------------
266{
267
268 delete fSigma[0];
269 delete fSigma[1];
270 delete fSigma[2];
271 delete fMPV[0];
272 delete fMPV[1];
273 delete fMPV[2];
274}
275// -------------------------------------------------------------------------
276
277// ----- Private method Reset -------------------------------------------
283// -------------------------------------------------------------------------
284
285// ----- Setting Detector related parameters
286// -------------------------------------------
288{
289 Int_t DetectorType = module->GetDetectorType();
290 // cout<<" dec type "<<DetectorType<<endl;
291 switch (DetectorType) {
292 case 3: // For GEM
293 {
294 fNADCChannels = 32; // Total ADC Channels
295 fQMax = 500000; // Maximum charge = 80 fC for GEM
296 fQThreshold = 12500; // Minimum charge threshold = 2fC for GEM
297 fMeanNoise = 1500; // mean noise
298 fSpotRadius = 0.05; // Spot radius = 500 um for GEM
299 fMeanGasGain = 5000; // Mean gas gain = 5000 for GEM
300 fDTime = 3;
301 fDeadPadsFrac = 0;
302 fDriftVelocity = 100; // Drift Velocity = 100 um/ns
303 auto DriftVolumeWidth = module->GetSize().Z(); // Drift gap
304 fTotalDriftTime = DriftVolumeWidth / fDriftVelocity * 10000; // Drift time (ns)
305 } break;
306 case 4: { // For RPC
307 fNADCChannels = 32; // Total ADC Channels
308 fQMax = 812500; // Maximum charge = 130 fC for RPC
309 fQThreshold = 187500; // Minimum charge threshold = 30 fC for RPC
310 fMeanNoise = 1500; // mean noise
311 fSpotRadius = 0.2; // Spot radius = 2 mm for RPC
312 fMeanGasGain = 3e4; // Mean gas gain = 3e4 for RPC
313 fDTime = 3;
314 fDeadPadsFrac = 0;
315 fDriftVelocity = 120; // Drift velocity = 120 um/ns
316 auto DriftVolumeWidth = module->GetSize().Z(); // Drift gap
317 fTotalDriftTime = DriftVolumeWidth / fDriftVelocity * 10000; // Drift time (ns)
318
319 } break;
320 }
321 // set parameters for GEM
322}
323// -------------------------------------------------------------------------
324
325// ----- Private method Init -------------------------------------------
327{
328
329 // Screen output
330 std::cout << std::endl;
331 LOG(info) << "==========================================================";
332 LOG(info) << GetName() << ": Initialisation\n";
333
334 FairRootManager* ioman = FairRootManager::Instance();
335 if (!ioman) Fatal("Init", "No FairRootManager");
336
337 if (fEventMode) {
338 LOG(info) << fName << ": Using event-by-event mode";
339 }
340
341 // hPriElAfterDriftpathgem = new TH1D("hPriElAfterDriftpathgem"," GEM:primary
342 // electron ",300,0,300); hPriElAfterDriftpathrpc = new
343 // TH1D("hPriElAfterDriftpathrpc"," RPC:primary electron ",300,0,300);
344 // hadcGEM = new TH1F("hadcGEM","GEM:ADC",70,-0.5,69.5);
345 // hadcRPC = new TH1F("hadcRPC","RPC:ADC",70,-0.5,69.5);
346
347 // hdrifttimegem = new TH1D("hdrifttimegem"," GEM ",100,0,50);
348 // hdrifttimerpc = new TH1D("hdrifttimerpc"," RPC ",100,0,50);
349
350 // Get geometry version tag
351 gGeoManager->CdTop();
352 TGeoNode* cave = gGeoManager->GetCurrentNode(); // cave
353 TString geoTag;
354 for (Int_t iNode = 0; iNode < cave->GetNdaughters(); iNode++) {
355 TString name = cave->GetDaughter(iNode)->GetVolume()->GetName();
356 if (name.Contains("much", TString::kIgnoreCase)) {
357 if (name.Contains("mcbm", TString::kIgnoreCase)) {
358 geoTag = TString(name(5, name.Length()));
359 }
360 else {
361 geoTag = TString(name(5, name.Length() - 5));
362 }
363 // geoTag = "v17b"; // modified ekata
364 // cout << " geo tag " << geoTag << endl;
365 LOG(info) << fName << ": MUCH geometry tag is " << geoTag;
366 break;
367 } //? node is MUCH
368 } //# top level nodes
369
370 // Set the parameter file and the flag, if not done in constructor
371 if (fDigiFile.IsNull()) {
372 if (geoTag.IsNull())
373 LOG(fatal) << fName
374 << ": no parameter file specified and no MUCH node found in "
375 "geometry!";
376 fDigiFile = gSystem->Getenv("VMCWORKDIR");
377 // TODO: (VF) A better naming convention for the geometry tag and the
378 // corresponding parameter file is surely desirable.
379 if (geoTag.Contains("mcbm", TString::kIgnoreCase)) {
380 fDigiFile += "/parameters/much/much_" + geoTag + "_digi_sector.root";
381 }
382 else {
383 fDigiFile += "/parameters/much/much_" + geoTag(0, 4) + "_digi_sector.root";
384 }
385 // fDigiFile
386 // ="/home/ekata/CbmRoot/AUG18/parameters/much/much_v17b_digi_sector.root";
387
388 LOG(info) << fName << ": Using parameter file " << fDigiFile;
389
390 fFlag = (geoTag.Contains("mcbm", TString::kIgnoreCase) ? 1 : 0);
391 LOG(info) << fName << ": Using flag " << fFlag << (fFlag ? " (mcbm) " : "(standard)");
392 }
393
394 // Initialize GeoScheme
396 TFile* oldFile = gFile;
397 TDirectory* oldDir = gDirectory;
398 TFile* file = new TFile(fDigiFile);
399 LOG_IF(fatal, !file->IsOpen()) << fName << ": parameter file " << fDigiFile << " does not exist!";
400 TObjArray* stations = file->Get<TObjArray>("stations");
401 LOG_IF(fatal, !stations) << "No TObjArray stations found in file " << fDigiFile;
402 file->Close();
403 file->Delete();
405 gFile = oldFile;
406 gDirectory = oldDir;
407 fGeoScheme->Init(stations, fFlag);
408
409 // For X, Y position correction according to different geometry file and its
410 // tag.
411 if (fDigiFile.Contains("mcbm")) {
412 if (fDigiFile.Contains("v19a")) {
413 fGemTX = 18.5;
414 fGemTY = 80.5;
415 }
416 else if (fDigiFile.Contains("v20a")) {
417 fGemTX = 18.5;
418 fGemTY = 80.0;
419 }
420 else if (fDigiFile.Contains("v22j")) { // for high intensity runs during
421 // June 2022 mMuCh in acceptance
422 fGemTX = 8.0;
423 fGemTY = 81.5;
424 fRpcTX = 66.25; // RPC introduced during 2022 data taking
425 fRpcTY = -70.5;
426 }
427 else if (fDigiFile.Contains("v22k")) { // during benchmark runs. mMuCh is
428 // out of acceptance
429 fGemTX = -44.5;
430 fGemTY = 81.5;
431 fRpcTX = 48.0;
432 fRpcTY = -70.0;
433 }
434 else if (fDigiFile.Contains("v24a")) { // Similar GEM orientation as CBM
435 fGemTX = 33.0;
436 fGemTY = 0.0;
437 }
438 else if (fDigiFile.Contains("v25a")) { // Similar GEM orientation as CBM
439 fGemTX = 36.5;
440 fGemTY = 0.0;
441 }
442 }
443
444 // Determine drift volume width
445 // Double_t driftVolumeWidth = 0.4; // cm - default and will over written by
446 // Detector Parameters function.
447
448 // Double_t DriftVolumeWidth;
449 for (Int_t i = 0; i < fGeoScheme->GetNStations(); i++) {
450 CbmMuchStation* station = fGeoScheme->GetStation(i);
451 // CbmMuchStation* station = fGeoScheme->GetStation(3);
452 if (station->GetNLayers() <= 0) continue;
453 CbmMuchLayerSide* side = station->GetLayer(0)->GetSide(0);
454 if (side->GetNModules() <= 0) continue;
455 CbmMuchModule* module = side->GetModule(0);
456 // if (module->GetDetectorType()!=4)continue;
457 if (module->GetDetectorType() != 3 && module->GetDetectorType() != 4) continue;
458 DetectorParameters(module);
459 // cout<<" dec type "<<module->GetDetectorType()<<" drift width
460 // "<<DriftVolumeWidth<<endl;
461 }
462
463 // cout<<" dec type "<<module->GetDetectorType()<<" drift width
464 // "<<driftVolumeWidth<<endl;
465 /*
466 driftVolumeWidth = module->GetSize().Z();
467 cout<<"driftVolumeWidth "<<driftVolumeWidth<<endl;
468 break;
469
470 }
471 fTotalDriftTime = driftVolumeWidth/fDriftVelocity*10000; // [ns];
472*/
473
474 // Reading MC point as Event by event for time based digi generation also.
475 // Get input array of MuchPoints
476 fPoints = (TClonesArray*) ioman->GetObject("MuchPoint");
477 assert(fPoints);
478
479 // Get input array of MC tracks
480 fMCTracks = (TClonesArray*) ioman->GetObject("MCTrack");
481 assert(fMCTracks);
482
483 // Register output arrays in vector
484 /*
485 fDigis = new std::vector<CbmMuchDigi>();
486 ioman->RegisterAny("MuchDigi",fDigis, IsOutputBranchPersistent("MuchDigi"));
487 if ( fCreateMatches ) {
488 fDigiMatches = new std::vector<CbmMatch>();
489 ioman->RegisterAny("MuchDigiMatch",fDigiMatches,
490 IsOutputBranchPersistent("MuchDigiMatch"));
491 }
492 */
494
495 // Register output arrays in TClonesArray which is replaced with vector
496 /*fDigis = new TClonesArray("CbmMuchDigi", 1000);
497 ioman->Register("MuchDigi", "Digital response in MUCH", fDigis, kTRUE);
498 if ( fCreateMatches ) {
499 fDigiMatches = new TClonesArray("CbmMatch", 1000);
500 ioman->Register("MuchDigiMatch", "Digi Match in MUCH", fDigiMatches, kTRUE);
501 }*/
502 // fgDeltaResponse is used in the CbmMuchSignal for analysing the Signal
503 // Shape, it is generated once in the digitizer and then be used by each
504 // CbmMuchSignal. For reducing the time therefore generated once in the
505 // CbmMuchDigitize Gem and not generated in the CbmMuchSignal
506 // Set response on delta function
507
508 // Not using fSignalShape as consuming memory. Commening all the field
509 // related.
510
511 // Int_t nShapeTimeBins=Int_t(gkResponsePeriod/gkResponseBin);
512 // fgDeltaResponse.Set(nShapeTimeBins);
513 // for (Int_t i=0;i<fgDeltaResponse.GetSize();i++){
514 // Double_t time = i*gkResponseBin;
515 // if (time<=fPeakingTime) fgDeltaResponse[i]=time/fPeakingTime;
516 // else fgDeltaResponse[i] = exp(-(time-fPeakingTime)/fRemainderTime);
517 //}
518
519 // --- Read list of inactive channels
520 if (!fInactiveChannelFileName.IsNull()) {
521 LOG(info) << GetName() << ": Reading inactive channels from " << fInactiveChannelFileName;
522 auto result = ReadInactiveChannels();
523 if (!result.second) {
524 LOG(error) << GetName() << ": Error in reading from file! Task will be inactive.";
525 return kFATAL;
526 }
527 LOG(info) << GetName() << ": " << std::get<0>(result) << " lines read from file, " << fInactiveChannels.size()
528 << " channels set inactive";
529 }
530
531 // --- Enable histogram if want to enalyze Noise spectrum.
532 // noise = new TH1D("noise", "Noise Generated per Event NoiseRate 10e-8", 100
533 // , 0 , 200);
534 LOG(info) << GetName() << ": Initialisation successful";
535 LOG(info) << "==========================================================";
536 std::cout << std::endl;
537
538 return kSUCCESS;
539}
540// -------------------------------------------------------------------------
541
542// ----- Public method Exec --------------------------------------------
544{
545
546 // --- Start timer and reset counters
547 fTimer.Start();
548 Reset();
549
550 // --- Event number and time
551 GetEventInfo();
552 LOG(debug) << GetName() << ": Processing event " << fCurrentEvent << " from input " << fCurrentInput
553 << " at t = " << fCurrentEventTime << " ns with " << fPoints->GetEntriesFast() << " MuchPoints ";
554
555 // ReadAndRegister(fCurrentEventTime);
556
557 for (Int_t iPoint = 0; iPoint < fPoints->GetEntriesFast(); iPoint++) {
558 const CbmMuchPoint* point = (const CbmMuchPoint*) fPoints->At(iPoint);
559 LOG(debug1) << GetName() << ": Processing MCPoint " << iPoint;
560 assert(point);
561 ExecPoint(point, iPoint);
562 fNofPoints++;
563 } // MuchPoint loop
564
565 // --------------------NOISE Generated before
566 // ReadAndRegister---------------------- //
570 fNofNoiseTot += nNoise;
571 // noise->Fill(nNoise);
572 LOG(info) << "+ " << setw(20) << GetName() << ": Generated " << nNoise
573 << " noise signals from t = " << fPreviousEventTime << " ns to " << fCurrentEventTime << " ns and "
574 << fNofNoiseTot << " total noise generated till now.";
575 LOG(debug3) << "+ " << setw(20) << GetName() << ": Generated " << fNofNoiseSignals
576 << " noise signals for this time slice from t = " << fPreviousEventTime << " ns to "
577 << fCurrentEventTime << "ns";
579 }
580
581 if (fEventMode)
582 ReadAndRegister(-1);
583 else
585
586 // --- Event log
587 LOG(info) << "+ " << setw(15) << GetName() << ": Event " << setw(6) << right << fCurrentEvent << " at " << fixed
588 << setprecision(3) << fCurrentEventTime << " ns, points: " << fNofPoints << ", signals: " << fNofSignals
589 << ", digis: " << fNofDigis << ". Exec time " << setprecision(6) << fTimer.RealTime() << " s.";
590
591 fTimer.Stop();
592 fNofEvents++;
596 fTimeTot += fTimer.RealTime();
597
598} // -----------------------------------------------------------------------------------------
599
600//================================Generate
601// Noise==============================================//
603{
604 LOG(debug) << "+ " << setw(20) << GetName() << ": Previous event time " << t1 << " Current event time " << t2
605 << " ns.";
606 if (t1 > t2) {
607 LOG(debug) << "+ "
608 << ": Previous event time " << t1 << " is greater than Current event time " << t2
609 << ". No electronics noise signal generated for " << fCurrentEvent << " event.";
610 return 0;
611 };
612 Int_t numberofstations = fGeoScheme->GetNStations();
613 auto StationNoise = 0;
614 for (Int_t i = 0; i < numberofstations; i++) {
615 CbmMuchStation* station = fGeoScheme->GetStation(i);
616 auto numberoflayers = station->GetNLayers();
617 if (numberoflayers <= 0) continue;
618
619 auto LayerNoise = 0;
620 for (Int_t j = 0; j < numberoflayers; j++) {
621 CbmMuchLayerSide* side = station->GetLayer(j)->GetSide(0);
622 auto numberofmodules = side->GetNModules();
623 if (numberofmodules <= 0) continue;
624
625 auto FrontModuleNoise = 0;
626 for (auto k = 0; k < numberofmodules; k++) {
627 CbmMuchModuleGem* module = (CbmMuchModuleGem*) (side->GetModule(k));
628 if (module->GetDetectorType() != 4 && module->GetDetectorType() != 3) continue;
629 FrontModuleNoise += GenerateNoisePerModule(module, t1, t2);
630 }
631 side = station->GetLayer(j)->GetSide(1);
632 numberofmodules = side->GetNModules();
633 if (numberofmodules <= 0) continue;
634 auto BackModuleNoise = 0;
635 for (auto k = 0; k < numberofmodules; k++) {
636 CbmMuchModuleGem* module = (CbmMuchModuleGem*) side->GetModule(k);
637 if (module->GetDetectorType() != 4 && module->GetDetectorType() != 3) continue;
638 BackModuleNoise += GenerateNoisePerModule(module, t1, t2);
639 }
640 LayerNoise += FrontModuleNoise + BackModuleNoise;
641 }
642 LOG(debug1) << "+ " << setw(20) << GetName() << ": Generated " << LayerNoise << " noise signals in station " << i
643 << " from t = " << fPreviousEventTime << " ns to " << fCurrentEventTime << " ns";
644 StationNoise += LayerNoise;
645 }
646 return StationNoise;
647}
648
649//================================Generate
650// Noise==============================================//
651
653{
654 auto NumberOfPad = module->GetNPads();
655 Double_t nNoiseMean = fPerPadNoiseRate * NumberOfPad * (t2 - t1);
656 Int_t nNoise = gRandom->Poisson(nNoiseMean);
657 LOG(debug) << "+ " << setw(20) << GetName() << ": Number of noise signals : " << nNoise << " in one module. ";
658 for (Int_t iNoise = 0; iNoise <= nNoise; iNoise++) {
659 Int_t padnumber = Int_t(gRandom->Uniform(Double_t(NumberOfPad)));
660 CbmMuchPad* pad = module->GetPad(padnumber);
661 Double_t NoiseTime = gRandom->Uniform(t1, t2);
662 Double_t charge = fNoiseCharge->GetRandom();
663 while (charge < 0)
664 charge = fNoiseCharge->GetRandom();
665 AddNoiseSignal(pad, NoiseTime, charge);
666 }
667 // noise->Fill(nNoise);
668 return nNoise;
669}
670
671//=================Add a signal to the buffer=====================//
672
673void CbmMuchDigitizeGem::AddNoiseSignal(CbmMuchPad* pad, Double_t time, Double_t charge)
674{
675 assert(pad);
676 LOG(debug3) << GetName() << ": Receiving signal " << charge << " in channel " << pad->GetAddress() << " at time "
677 << time << "ns";
678 // LOG(debug) << GetName() << ": discarding signal in dead channel "
679 // << channel;
680 // return;
681 CbmMuchSignal* signal = new CbmMuchSignal(pad->GetAddress(), time);
682 // signal->SetTimeStart(time);
683 // signal->SetTimeStop(time+fDeadTime);
684 signal->SetCharge((UInt_t) charge);
685 UInt_t address = pad->GetAddress();
686 CbmMuchReadoutBuffer::Instance()->Fill(address, signal);
688 LOG(debug3) << "+ " << setw(20) << GetName()
689 << ": Registered a Noise CbmMuchSignal into the "
690 "CbmMuchReadoutBuffer. Number of Noise Signal generated "
692}
693//====================End of Noise part=================//
694
695// -------------------------------------------------------------------------
696// Read all the Signal from CbmMuchReadoutBuffer, convert the analog signal into
697// the digital response and register Output according to event by event mode
698// and Time based mode.
700{
701 std::vector<CbmMuchSignal*> SignalList;
702 // Event Time should be passed with the Call
703 /*Double_t eventTime = -1.;
704 if(fDaq){
705 eventTime = FairRun::Instance()->GetEventHeader()->GetEventTime();
706 }*/
707
708 Int_t ReadOutSignal = CbmMuchReadoutBuffer::Instance()->ReadOutData(eventTime, SignalList);
709 LOG(debug) << GetName() << ": Number of Signals read out from Buffer " << ReadOutSignal << " and SignalList contain "
710 << SignalList.size() << " entries.";
711
712 for (std::vector<CbmMuchSignal*>::iterator LoopOver = SignalList.begin(); LoopOver != SignalList.end(); LoopOver++) {
713 CbmMuchDigi* digi = ConvertSignalToDigi(*LoopOver);
714 // assert(digi);
715 if (!digi) {
716 LOG(debug2) << GetName() << ": Digi not created as signal is below threshold.";
717 }
718 else {
719 LOG(debug2) << GetName() << ": New digi: sector = " << CbmMuchAddress::GetSectorIndex(digi->GetAddress())
720 << " channel= " << CbmMuchAddress::GetChannelIndex(digi->GetAddress());
721
722 if (fCreateMatches) {
723 CbmMatch* digiMatch = new CbmMatch(*(*LoopOver)->GetMatch()); // must be copied from signal
724 SendData(digi->GetTime(), digi, digiMatch);
725 }
726 else {
727 SendData(digi->GetTime(), digi);
728 }
729 fNofDigis++;
730 }
731 }
732
733 LOG(debug) << GetName() << ": " << fNofDigis << (fNofDigis == 1 ? " digi " : " digis ") << "created and sent to DAQ ";
734 if (fNofDigis)
735 LOG(debug) << "( from " << fixed << setprecision(3) << fTimeDigiFirst << " ns to " << fTimeDigiLast << " ns )";
736 LOG(debug);
737
738 // After digis are created from signals the signals have to be removed
739 // Otherwise there is a huge memeory leak
740 for (auto signal : SignalList) {
741 delete (signal);
742 }
743
744} //----ReadAndRegister -------
745
746// Convert Signal into the Digi with appropriate methods.
747
749{
750
751 // Setting Parameters for RPC or GEM (10 lines)
752 auto address = signal->GetAddress();
753 auto StationIndex = CbmMuchAddress::GetStationIndex(address);
754 CbmMuchStation* station = fGeoScheme->GetStation(StationIndex);
755 auto LayerIndex = CbmMuchAddress::GetLayerIndex(address);
756 auto SideIndex = CbmMuchAddress::GetLayerSideIndex(address);
757 CbmMuchLayerSide* side = station->GetLayer(LayerIndex)->GetSide(SideIndex);
758 auto ModuleIndex = CbmMuchAddress::GetLayerSideIndex(address);
759 CbmMuchModule* module = side->GetModule(ModuleIndex);
760 assert(module);
761 DetectorParameters(module);
762
763 // cout<<"Det Type: "<<module->GetDetectorType()<<" fQThreshold
764 // "<<fQThreshold<<" fqmax "<<fQMax<<" drift time "<<fTotalDriftTime<<" drift
765 // vel "<<fDriftVelocity<<" spot radius "<<fSpotRadius<<endl;
766 // cout<<"fQThreshold "<<fQThreshold<<" fqmax "<<fQMax<<endl;
767 // signal below threshold should be discarded.
768 if (signal->GetCharge() < 0)
769 return (NULL); // Before type casting signed int into unsigned int need to
770 // check for -ve value otherwise after casing -ve int value
771 // will become a large +ve value.
772 else if ((unsigned int) signal->GetCharge() < fQThreshold)
773 return (NULL);
774
775 Long64_t TimeStamp = signal->GetTimeStamp();
776 // Int_t TimeStamp = signal->GetTimeStamp(fQThreshold);
777 // if (TimeStamp < 0) return (NULL);//entire signal is below threshold, no
778 // digi generation.
779
780 CbmMuchDigi* digi = new CbmMuchDigi();
781 digi->SetAddress(signal->GetAddress());
782 // Charge in number of electrons, need to be converted in ADC value
783 // digi->SetAdc((signal->GetCharge())*fNADCChannels/fQMax);//Charge should
784 // be computed as per Electronics Response.
785 Float_t adcValue;
786 adcValue = ((signal->GetCharge() - fQThreshold) * fNADCChannels) / (fQMax - fQThreshold);
787 digi->SetAdc(adcValue + 1); // Charge should be computed as per Electronics
788 // Response. *** modified by Ekata Nandy***
789 // ADC channels number is 32.GEM & RPC has different charge threshold value
790 // and dynamic range, so SetAdc has been changed accordingly. ADC value starts
791 // from 1 to 32. ADC 0 has been excluded as it gives wrong x,y,t in Hit
792 // finder.@modified by Ekata Nandy
793
794 // cout<<" dynamic range "<<(fQMax -fQThreshold)<<" adc "<<digi->GetAdc()<<"
795 // channels "<<fNADCChannels<< "charge "<<signal->GetCharge()<<endl;
796 /* if(module->GetDetectorType()==3)
797 {hadcGEM->Fill(digi->GetAdc());}
798 if(module->GetDetectorType()==4)
799 {hadcRPC->Fill(digi->GetAdc());}
800 */
801 digi->SetTime(TimeStamp);
802
803 // Update times of first and last digi
804 fTimeDigiFirst = fNofDigis ? TMath::Min(fTimeDigiFirst, Double_t(TimeStamp)) : TimeStamp;
805 fTimeDigiLast = TMath::Max(fTimeDigiLast, Double_t(TimeStamp));
806
807 // digi->SetPileUp();
808 // digi->SetDiffEvent();
809 return (digi);
810}
811
812// -------------------------------------------------------------------------
814{
815
816 // --- In event-by-event mode, the analogue buffer should be empty.
817 if (fEventMode) {
818 if (CbmMuchReadoutBuffer::Instance()->GetNData()) {
819 LOG(info) << fName << ": " << CbmMuchReadoutBuffer::Instance()->GetNData() << " signals in readout buffer";
820 LOG(fatal) << fName << ": Readout buffer is not empty at end of run "
821 << "in event-by-event mode!";
822 } //? non-empty buffer
823 } //? event-by-event mode
824
825 else { // time-based mode
826 fTimer.Start();
827 std::cout << std::endl;
828 LOG(info) << GetName() << ": Finish run";
829 Reset();
830 LOG(info) << fName << ": " << CbmMuchReadoutBuffer::Instance()->GetNData() << " signals in readout buffer";
831 ReadAndRegister(-1.); // -1 means process all data
832 LOG(info) << setw(15) << GetName() << ": Finish, points " << fNofPoints << ", signals: " << fNofSignals
833 << ", digis: " << fNofDigis << ". Exec time " << setprecision(6) << fTimer.RealTime() << " s.";
834 LOG(info) << fName << ": " << CbmMuchReadoutBuffer::Instance()->GetNData() << " signals in readout buffer";
835 fTimer.Stop();
839 fTimeTot += fTimer.RealTime();
840 } //? time-based mode
841 // noise->Draw();
842 std::cout << std::endl;
843 LOG(info) << "=====================================";
844 LOG(info) << GetName() << ": Run summary";
845 LOG(info) << "Events processed : " << fNofEvents;
846 LOG(info) << "MuchPoint / event : " << setprecision(1) << fNofPointsTot / Double_t(fNofEvents);
847 LOG(info) << "MuchSignal / event : " << fNofSignalsTot / Double_t(fNofEvents);
848 //<< " / " << fNofSignalsBTot / Double_t(fNofEvents)
849 LOG(info) << "MuchDigi / event : " << fNofDigisTot / Double_t(fNofEvents);
850 LOG(info) << "Digis per point : " << setprecision(6) << fNofDigisTot / fNofPointsTot;
851 LOG(info) << "Digis per signal : " << fNofDigisTot / fNofSignalsTot;
852 LOG(info) << "Noise digis / event : " << fNofNoiseTot / Double_t(fNofEvents);
853 LOG(info) << "Noise fraction : " << fNofNoiseTot / fNofDigisTot;
854
855 LOG(info) << "Real time per event : " << fTimeTot / Double_t(fNofEvents) << " s";
856 LOG(info) << "=====================================";
857
858 /*
860 TFile* oldFile = gFile;
861 TDirectory* oldDir = gDirectory;
862 TFile *f1 =new TFile ("pri_el_info.root","RECREATE");
863 hPriElAfterDriftpathgem->Write();
864 hPriElAfterDriftpathrpc->Write();
865 hadcGEM->Write();
866 hadcRPC->Write();
867 f1->Close();
869 gFile = oldFile;
870 gDirectory = oldDir;
871 */
872 // if (fDaq) ReadAndRegister(-1.);
873}
874// -------------------------------------------------------------------------
875
876// ------- Private method ExecAdvanced -------------------------------------
878{
879
880 // std::cout<<" start execution "<<iPoint<<std::endl;
881 TVector3 v1, v2, dv;
882 point->PositionIn(v1);
883 point->PositionOut(v2);
884 dv = v2 - v1;
885
886 Bool_t Status;
887 Int_t detectorId = point->GetDetectorID();
888 CbmMuchModule* module = fGeoScheme->GetModuleByDetId(detectorId);
889 DetectorParameters(module);
890
891 int detectortype = module->GetDetectorType();
892 // cout<<" dec type === "<<detectortype<<endl;
893 if (fAlgorithm == 0) { // Simple digitization
894 TVector3 v = 0.5 * (v1 + v2);
895 CbmMuchPad* pad = 0;
896 if (module->GetDetectorType() == 1) {
898 pad = module1->GetPad(v[0], v[1]);
899 if (pad) printf("x0=%f,y0=%f\n", pad->GetX(), pad->GetY());
900 }
901 else if (module->GetDetectorType() == 3 || module->GetDetectorType() == 4) {
903 pad = module3->GetPad(v[0], v[1]);
904 }
905 if (!pad) return kFALSE;
906 AddCharge(pad, fQMax, iPoint, point->GetTime(), 0);
907 return kTRUE;
908 }
909
910 // Start of advanced digitization
911 Int_t nElectrons = Int_t(GetNPrimaryElectronsPerCm(point, detectortype) * dv.Mag());
912 if (nElectrons < 0) return kFALSE;
913 // cout<<" nElectrons "<<nElectrons<<endl;
914 /*
915if(detectortype==3)
916{//cout<<" nElectrons "<<nElectrons<<endl;
917hPriElAfterDriftpathgem->Fill(nElectrons);}
918
919if(detectortype==4)
920{//cout<<" nElectrons "<<nElectrons<<endl;
921hPriElAfterDriftpathrpc->Fill(nElectrons);}
922*/
923
924 Double_t time = point->GetTime();
925
926 if (module->GetDetectorType() == 1) {
928 map<CbmMuchSector*, Int_t> firedSectors;
929 for (Int_t i = 0; i < nElectrons; i++) {
930 Double_t aL = gRandom->Rndm();
931 Double_t driftTime = (1 - aL) * fTotalDriftTime;
932 TVector3 ve = v1 + dv * aL;
933 UInt_t ne = GasGain();
934 Double_t x = ve.X();
935 Double_t y = ve.Y();
936 Double_t x1 = x - fSpotRadius;
937 Double_t x2 = x + fSpotRadius;
938 Double_t y1 = y - fSpotRadius;
939 Double_t y2 = y + fSpotRadius;
940 Double_t s = 4 * fSpotRadius * fSpotRadius;
941 firedSectors[module1->GetSector(x1, y1)] = 0;
942 firedSectors[module1->GetSector(x1, y2)] = 0;
943 firedSectors[module1->GetSector(x2, y1)] = 0;
944 firedSectors[module1->GetSector(x2, y2)] = 0;
945 for (map<CbmMuchSector*, Int_t>::iterator it = firedSectors.begin(); it != firedSectors.end(); it++) {
946 CbmMuchSector* sector = (*it).first;
947 if (!sector) continue;
948 for (Int_t iPad = 0; iPad < sector->GetNChannels(); iPad++) {
949 CbmMuchPad* pad = sector->GetPadByChannelIndex(iPad);
950 Double_t xp0 = pad->GetX();
951 Double_t xpd = pad->GetDx() / 2.;
952 Double_t xp1 = xp0 - xpd;
953 Double_t xp2 = xp0 + xpd;
954 if (x1 > xp2 || x2 < xp1) continue;
955 Double_t yp0 = pad->GetY();
956 Double_t ypd = pad->GetDy() / 2.;
957 Double_t yp1 = yp0 - ypd;
958 Double_t yp2 = yp0 + ypd;
959 if (y1 > yp2 || y2 < yp1) continue;
960 Double_t lx = x1 > xp1 ? (x2 < xp2 ? x2 - x1 : xp2 - x1) : x2 - xp1;
961 Double_t ly = y1 > yp1 ? (y2 < yp2 ? y2 - y1 : yp2 - y1) : y2 - yp1;
962 AddCharge(pad, UInt_t(ne * lx * ly / s), iPoint, time, driftTime);
963 }
964 } // loop fired sectors
965 firedSectors.clear();
966 }
967 }
968
969 if (module->GetDetectorType() == 3 || module->GetDetectorType() == 4) {
970 fAddressCharge.clear();
972 if (!module3) {
973 LOG(debug) << GetName() << ": Not Processing MCPoint " << iPoint << " because it is not on any GEM module.";
974 return 1;
975 }
976 CbmMuchSectorRadial* sFirst = (CbmMuchSectorRadial*) module3->GetSectorByIndex(0); // First sector
977 if (!sFirst) {
978 LOG(debug) << GetName() << ": Not Processing MCPoint " << iPoint << " because it is on the module " << module3
979 << " but not the first sector. " << sFirst;
980 return 1;
981 }
982 CbmMuchSectorRadial* sLast =
983 (CbmMuchSectorRadial*) module3->GetSectorByIndex(module3->GetNSectors() - 1); // Last sector
984
985 if (!sLast) {
986 LOG(debug) << GetName() << ": Not Processing MCPoint " << iPoint
987 << " because it is not the last sector of module." << module3;
988 return 1;
989 }
990 Double_t rMin = sFirst->GetR1(); // Mimimum radius of the Sector//5
991 Double_t rMax = sLast->GetR2(); // Maximum radius of the Sector//35
992
993 // cout<<rMin<<" Yeah "<<rMax<<endl;
994 // std::cout<<"min Rad "<<rMin<<" max Rad "<<rMax<<std::endl;
995 // Calculating drifttime once for one track or one MCPoint, not for all the
996 // Primary Electrons generated during DriftGap.
997
998 Double_t driftTime = -1;
999 while (driftTime < 0) {
1000
1001 Double_t aL = gRandom->Gaus(0.5, 0.133); // Generting random number for calculating Drift Time.
1002 // cout<<"Det Type "<<detectortype<<"fTotalDriftTime
1003 // "<<fTotalDriftTime<<endl;
1004 driftTime = (1 - aL) * fTotalDriftTime;
1005 }
1006
1007 for (Int_t i = 0; i < nElectrons; i++) { // Looping over all the primary electrons
1008 Double_t RandomNumberForPrimaryElectronPosition = gRandom->Rndm();
1009 TVector3 ve = v1 + dv * RandomNumberForPrimaryElectronPosition;
1010
1011 //------------------------Added by O. Singh 11.12.2017 for
1012 // mCbm-------------------------
1013 Double_t r = 0.0, phi = 0.0;
1014 if (fFlag == 1) { // mCbm
1015 if (module->GetDetectorType() == 3) // GEM
1016 {
1017 ve.SetX(ve.X() - fGemTX);
1018 ve.SetY(ve.Y() - fGemTY);
1019 if (fDigiFile.Contains("v25a")) // GEM
1020 {
1021 if (CbmMuchAddress::GetStationIndex(module3->GetDetectorId()) == 1) {
1022 // During May 2025, GEM-2 at mCBM will be Real Size Station-2 Module
1023 ve.SetX(ve.X() - fGemTX + 24.0); //Nilay
1024 ve.SetY(ve.Y() - fGemTY);
1025 }
1026 }
1027 }
1028 else if (module->GetDetectorType() == 4) // RPC
1029 {
1030 ve.SetX(ve.X() - fRpcTX);
1031 ve.SetY(ve.Y() - fRpcTY);
1032 }
1033 else {
1034 LOG(error) << "Unknown detector type";
1035 }
1036 }
1037 r = ve.Perp();
1038 phi = ve.Phi();
1039 //--------------------------------------------------------------------------
1040 UInt_t ne = GasGain(); // Number of secondary electrons
1041 Double_t r1 = r - fSpotRadius;
1042 Double_t r2 = r + fSpotRadius;
1043 Double_t phi1 = phi - fSpotRadius / r;
1044 Double_t phi2 = phi + fSpotRadius / r;
1045 // cout<<" fSpotRadius "<<fSpotRadius<<endl;
1046 if (r1 < rMin && r2 > rMin) { // Adding charge to the pad which is on Lower Boundary
1047 Status = AddCharge(sFirst, UInt_t(ne * (r2 - rMin) / (r2 - r1)), iPoint, time, driftTime, phi1, phi2);
1048 if (!Status)
1049 LOG(debug) << GetName() << ": Processing MCPoint " << iPoint << " in which Primary Electron : " << i
1050 << " not contributed charge. ";
1051 continue;
1052 }
1053 if (r1 < rMax && r2 > rMax) { // Adding charge to the pad which is on Upper Boundary
1054 Status = AddCharge(sLast, UInt_t(ne * (rMax - r1) / (r2 - r1)), iPoint, time, driftTime, phi1, phi2);
1055 if (!Status)
1056 LOG(debug) << GetName() << ": Processing MCPoint " << iPoint << " in which Primary Electron : " << i
1057 << " not contributed charge. ";
1058 continue;
1059 }
1060 if (r1 < rMin && r2 < rMin) continue;
1061 if (r1 > rMax && r2 > rMax) continue;
1062
1063 CbmMuchSectorRadial* s1 = module3->GetSectorByRadius(r1);
1064 CbmMuchSectorRadial* s2 = module3->GetSectorByRadius(r2);
1065
1066 if (s1 == s2) {
1067 Status = AddCharge(s1, ne, iPoint, time, driftTime, phi1, phi2);
1068 if (!Status)
1069 LOG(debug3) << GetName() << ": Processing MCPoint " << iPoint << " in which Primary Electron : " << i
1070 << " not contributed charge. ";
1071 }
1072 else { // Adding praportionate charge to both the pad
1073 Status = AddCharge(s1, UInt_t(ne * (s1->GetR2() - r1) / (r2 - r1)), iPoint, time, driftTime, phi1, phi2);
1074 if (!Status)
1075 LOG(debug3) << GetName() << ": Processing MCPoint " << iPoint << " in which Primary Electron : " << i
1076 << " not contributed charge. ";
1077 Status = AddCharge(s2, UInt_t(ne * (r2 - s2->GetR1()) / (r2 - r1)), iPoint, time, driftTime, phi1, phi2);
1078 if (!Status)
1079 LOG(debug3) << GetName() << ": Processing MCPoint " << iPoint << " in which Primary Electron : " << i
1080 << " not contributed charge. ";
1081 }
1082 }
1083
1084 // Generate CbmMuchSignal for each entry of fAddressCharge and store in the
1085 // CbmMuchReadoutBuffer
1086 if (!BufferSignals(iPoint, time, driftTime))
1087 LOG(debug3) << GetName() << ": Processing MCPoint " << iPoint << " nothing is buffered. ";
1088 fAddressCharge.clear();
1089 LOG(debug1) << GetName() << ": fAddressCharge size is " << fAddressCharge.size() << " Cleared fAddressCharge. ";
1090 }
1091 // std::cout<<" Execution completed for point # "<<iPoint<<std::endl;
1092 return kTRUE;
1093}
1094// -------------------------------------------------------------------------
1095
1096// -------------------------------------------------------------------------
1098{
1099 // cout<<" mean gain "<<fMeanGasGain<<endl;
1100 Double_t gasGain = -fMeanGasGain * TMath::Log(1 - gRandom->Rndm());
1101 if (gasGain < 0.) gasGain = 1e6;
1102 return (Int_t) gasGain;
1103}
1104// -------------------------------------------------------------------------
1105
1106// -------------------------------------------------------------------------
1107Double_t CbmMuchDigitizeGem::Sigma_n_e(Double_t Tkin, Double_t mass)
1108{
1109 Double_t logT;
1110 if (mass < 0.1) {
1111 logT = log(Tkin * 0.511 / mass);
1112 if (logT > 9.21034) logT = 9.21034;
1113 if (logT < min_logT_e) logT = min_logT_e;
1114 return fSigma[0]->Eval(logT);
1115 }
1116 else if (mass >= 0.1 && mass < 0.2) {
1117 logT = log(Tkin * 105.658 / mass);
1118 if (logT > 9.21034) logT = 9.21034;
1119 if (logT < min_logT_mu) logT = min_logT_mu;
1120 return fSigma[1]->Eval(logT);
1121 }
1122 else {
1123 logT = log(Tkin * 938.272 / mass);
1124 if (logT > 9.21034) logT = 9.21034;
1125 if (logT < min_logT_p) logT = min_logT_p;
1126 return fSigma[2]->Eval(logT);
1127 }
1128}
1129// -------------------------------------------------------------------------
1130
1131// -------------------------------------------------------------------------
1132Double_t CbmMuchDigitizeGem::MPV_n_e(Double_t Tkin, Double_t mass)
1133{
1134 Double_t logT;
1135 if (mass < 0.1) {
1136 logT = log(Tkin * 0.511 / mass);
1137 if (logT > 9.21034) logT = 9.21034;
1138 if (logT < min_logT_e) logT = min_logT_e;
1139 return fMPV[0]->Eval(logT);
1140 }
1141 else if (mass >= 0.1 && mass < 0.2) {
1142 logT = log(Tkin * 105.658 / mass);
1143 if (logT > 9.21034) logT = 9.21034;
1144 if (logT < min_logT_mu) logT = min_logT_mu;
1145 return fMPV[1]->Eval(logT);
1146 }
1147 else {
1148 logT = log(Tkin * 938.272 / mass);
1149 if (logT > 9.21034) logT = 9.21034;
1150 if (logT < min_logT_p) logT = min_logT_p;
1151 return fMPV[2]->Eval(logT);
1152 }
1153}
1154// -------------------------------------------------------------------------
1155
1156// -------------------------------------------------------------------------
1157Double_t CbmMuchDigitizeGem::GetNPrimaryElectronsPerCm(const CbmMuchPoint* point, int detectortype)
1158{
1159 Int_t trackId = point->GetTrackID();
1160 // Int_t eventId = point->GetEventID();
1161 if (trackId < 0) return -1;
1162 /* Commented out on request of A. Senger from 22.01.2014
1163 if (fDaq &&
1164 eventId!=FairRootManager::Instance()->GetInTree()->GetBranch("MCTrack")->GetReadEntry())
1165 FairRootManager::Instance()->GetInTree()->GetBranch("MCTrack")->GetEntry(eventId);
1166 */
1167 CbmMCTrack* mcTrack = (CbmMCTrack*) fMCTracks->At(trackId);
1168
1169 if (!mcTrack) return -1;
1170 Int_t pdgCode = mcTrack->GetPdgCode();
1171
1172 TParticlePDG* particle = TDatabasePDG::Instance()->GetParticle(pdgCode);
1173 // Assign proton hypothesis for unknown particles
1174 if (!particle) particle = TDatabasePDG::Instance()->GetParticle(2212);
1175 if (TMath::Abs(particle->Charge()) < 0.1) return -1;
1176
1177 Double_t m = particle->Mass();
1178 TLorentzVector p;
1179 p.SetXYZM(point->GetPx(), point->GetPy(), point->GetPz(), m);
1180 Double_t Tkin = p.E() - m; // kinetic energy of the particle
1181 Double_t sigma = CbmMuchDigitizeGem::Sigma_n_e(Tkin, m); // sigma for Landau distribution
1182 Double_t mpv = CbmMuchDigitizeGem::MPV_n_e(Tkin, m); // most probable value for Landau distr.
1183 Double_t n;
1184 Double_t mpvRpc = 50.0; // For RPC MPV and sigma is taken to produce final
1185 // landau in accordance with experimental value
1186 Double_t sigmaRpc = 12.0;
1187 // cout<<" dec type function "<<detectortype<<endl;
1188 if (detectortype == 3)
1189 {
1190 n = gRandom->Landau(mpv, sigma);
1191 while (n > 5e4)
1192 n = gRandom->Landau(mpv, sigma); // restrict Landau tail to increase performance
1193 return m < 0.1 ? n / l_e : n / l_not_e;
1194 }
1195 if (detectortype == 4)
1196 {
1197 n = gRandom->Landau(mpvRpc, sigmaRpc);
1198 while (n > 5e4)
1199 n = gRandom->Landau(mpvRpc, sigmaRpc); // restrict Landau tail to increase performance
1200 return n;
1201 }
1202 // for all other cases will return -1.
1203 return -1;
1204}
1205// -------------------------------------------------------------------------
1206// -------------------------------------------------------------------------
1207Bool_t CbmMuchDigitizeGem::AddCharge(CbmMuchSectorRadial* s, UInt_t ne, Int_t /*iPoint*/, Double_t /*time*/,
1208 Double_t /*driftTime*/, Double_t phi1, Double_t phi2)
1209{
1210 CbmMuchPadRadial* pad1 = s->GetPadByPhi(phi1);
1211 if (!pad1)
1212 pad1 = s->GetPadByPhi((phi1 + phi2) / 2.0); // This condition helps us deal with boundary pads
1213 else // Special case if pad size is smaller than spot radius
1214 {
1215 for (Double_t phi = phi1; phi < (phi1 + phi2) / 2.0; phi += 0.1) // This may potentially slow down the code
1216 {
1217 pad1 = s->GetPadByPhi(phi);
1218 if (pad1) break;
1219 }
1220 }
1221 if (!pad1) return kFALSE;
1222 // assert(pad1); has to check if any pad address is NULL
1223 CbmMuchPadRadial* pad2 = s->GetPadByPhi(phi2);
1224 if (!pad2)
1225 pad2 = s->GetPadByPhi((phi1 + phi2) / 2.0); // This condition helps us deal with boundary pads
1226 else // Special case if pad size is smaller than spot radius
1227 {
1228 for (Double_t phi = phi2; phi > (phi1 + phi2) / 2.0; phi -= 0.1) // This may potentially slow down the code
1229 {
1230 pad2 = s->GetPadByPhi(phi);
1231 if (pad2) break;
1232 }
1233 }
1234
1235 if (!pad2) return kFALSE;
1236 // assert(pad2); has to check if any pad address is NULL
1237 if (pad1 == pad2) {
1238 UInt_t address = pad1->GetAddress();
1239 // Finding that if for the same address if already charge stored then add
1240 // the charge.
1241 std::map<UInt_t, UInt_t>::iterator it = fAddressCharge.find(address);
1242 if (it != fAddressCharge.end())
1243 it->second = it->second + ne;
1244 else
1245 fAddressCharge.insert(std::pair<UInt_t, UInt_t>(address, ne));
1246 // AddChargePerMC(pad1,ne,iPoint,time,driftTime);
1247 }
1248 else {
1249 Double_t phi = pad1 ? pad1->GetPhi2() : pad2 ? pad2->GetPhi1() : 0;
1250 UInt_t pad1_ne = UInt_t(ne * (phi - phi1) / (phi2 - phi1));
1251
1252 UInt_t address = pad1->GetAddress();
1253 // Finding that if for the same address if already charge stored then add
1254 // the charge.
1255 std::map<UInt_t, UInt_t>::iterator it = fAddressCharge.find(address);
1256 if (it != fAddressCharge.end())
1257 it->second = it->second + pad1_ne;
1258 else
1259 fAddressCharge.insert(std::pair<UInt_t, UInt_t>(address, pad1_ne));
1260 // AddChargePerMC(pad1,pad1_ne ,iPoint,time,driftTime);
1261
1262 // Getting some segmentation fault a
1263 address = pad2->GetAddress();
1264 // Finding that if for the same address if already charge stored then add
1265 // the charge.
1266 it = fAddressCharge.find(address);
1267 if (it != fAddressCharge.end())
1268 it->second = it->second + ne - pad1_ne;
1269 else
1270 fAddressCharge.insert(std::pair<UInt_t, UInt_t>(address, ne - pad1_ne));
1271 // AddChargePerMC(pad2,ne-pad1_ne,iPoint,time,driftTime);
1272 }
1273 return kTRUE;
1274}
1275// -------------------------------------------------------------------------
1276
1277// -------------------------------------------------------------------------
1278// Will remove this AddCharge, only used for simple and Rectangular Geometry.
1279void CbmMuchDigitizeGem::AddCharge(CbmMuchPad* pad, UInt_t charge, Int_t iPoint, Double_t time, Double_t driftTime)
1280{
1281 if (!pad) return;
1282
1283 Long_t AbsTime = fCurrentEventTime + time + driftTime;
1284
1285 // Creating a new Signal, it will be deleted by CbmReadoutBuffer()
1286 CbmMuchSignal* signal = new CbmMuchSignal(pad->GetAddress(), AbsTime);
1287 // signal->SetTimeStart(AbsTime);
1288 // signal->SetTimeStop(AbsTime+fDeadTime);
1289 signal->SetCharge(charge);
1290 // signal->MakeSignalShape(charge,fgDeltaResponse);
1291 signal->AddNoise(fMeanNoise);
1292 UInt_t address = pad->GetAddress();
1293 // match->AddCharge(iPoint,charge,time+driftTime,fgDeltaResponse,time,eventNr,inputNr);
1294 CbmLink link(charge, iPoint, fCurrentMCEntry, fCurrentInput);
1295 // std::cout<<"Before AddLink"<< endl;
1296 (signal->GetMatch())->AddLink(link);
1297 // std::cout<<"After AddLink"<< endl;
1298 // Adding all these temporary signal into the CbmMuchReadoutBuffer
1299 CbmMuchReadoutBuffer::Instance()->Fill(address, signal);
1300 // Increasing number of signal by one.
1301 fNofSignals++;
1302 LOG(debug4) << " Registered the CbmMuchSignal into the CbmMuchReadoutBuffer ";
1303
1304} // end of AddCharge
1305
1306//----------------------------------------------------------
1307Bool_t CbmMuchDigitizeGem::BufferSignals(Int_t iPoint, Double_t time, Double_t driftTime)
1308{
1309
1310 if (!fAddressCharge.size()) {
1311 LOG(debug2) << "Buffering MC Point " << iPoint << " but fAddressCharge size is " << fAddressCharge.size()
1312 << "so nothing to Buffer for this MCPoint.";
1313 return kFALSE;
1314 }
1315 UInt_t AbsTime = fCurrentEventTime + time + driftTime;
1316 LOG(debug2) << GetName() << ": Processing event " << fCurrentEvent << " from input " << fCurrentInput
1317 << " at t = " << fCurrentEventTime << " ns with " << fPoints->GetEntriesFast() << " MuchPoints "
1318 << " and Number of pad hit is " << fAddressCharge.size() << ".";
1319 // Loop on the fAddressCharge to store all the Signals into the
1320 // CbmReadoutBuffer() Generate one by one CbmMuchSignal from the
1321 // fAddressCharge and store them into the CbmMuchReadoutBuffer.
1322 for (auto it = fAddressCharge.begin(); it != fAddressCharge.end(); ++it) {
1323 UInt_t address = it->first;
1324 // Creating a new Signal, it will be deleted by CbmReadoutBuffer()
1325 CbmMuchSignal* signal = new CbmMuchSignal(address, AbsTime);
1326 // signal->SetTimeStart(AbsTime);
1327 // signal->SetTimeStop(AbsTime+fDeadTime);
1328 // signal->MakeSignalShape(it->second,fgDeltaResponse);
1329 signal->SetCharge(it->second);
1330 signal->AddNoise(fMeanNoise);
1331 CbmLink link(signal->GetCharge(), iPoint, fCurrentMCEntry, fCurrentInput);
1332 (signal->GetMatch())->AddLink(link);
1333 // Adding all these temporary signal into the CbmMuchReadoutBuffer
1334 CbmMuchReadoutBuffer::Instance()->Fill(address, signal);
1335 // Increasing number of signal by one.
1336 fNofSignals++;
1337 LOG(debug3) << " Registered the CbmMuchSignal into the CbmMuchReadoutBuffer ";
1338 }
1339
1340 LOG(debug2) << GetName() << ": For MC Point " << iPoint << " buffered " << fAddressCharge.size()
1341 << " CbmMuchSignal into the CbmReadoutBuffer.";
1342 return kTRUE;
1343} // end of BufferSignals
1344// -------------------------------------------------------------------------
1345
ClassImp(CbmConverterManager)
constexpr double min_logT_e
constexpr double sigma_mu[]
constexpr double l_e
constexpr double sigma_p[]
constexpr double l_not_e
constexpr double mpv_e[]
constexpr double mpv_p[]
constexpr double min_logT_p
constexpr double sigma_e[]
constexpr double mpv_mu[]
constexpr double min_logT_mu
fscal v[fmask::Size]
Definition KfSimdPseudo.h:4
friend fvec log(const fvec &a)
float Float_t
int Int_t
bool Bool_t
virtual std::pair< size_t, bool > ReadInactiveChannels()
std::set< uint32_t > fInactiveChannels
void SendData(Double_t time, CbmMuchDigi *digi, CbmMatch *match=nullptr)
int32_t GetPdgCode() const
Definition CbmMCTrack.h:67
static int32_t GetSectorIndex(int32_t address)
static int32_t GetChannelIndex(int32_t address)
static int32_t GetLayerIndex(int32_t address)
static int32_t GetLayerSideIndex(int32_t address)
static int32_t GetStationIndex(int32_t address)
void SetAddress(int32_t address)
void SetTime(uint64_t time)
int32_t GetAddress() const
Definition CbmMuchDigi.h:93
void SetAdc(int32_t adc)
double GetTime() const
Definition CbmMuchDigi.h:94
Int_t fNofDigis
Number of created digis in Exec.
Double_t MPV_n_e(Double_t Tkin, Double_t mass)
Int_t fNofEvents
Total number of events processed.
Int_t fNofSignals
Number of signals.
void DetectorParameters(CbmMuchModule *module)
CbmMuchGeoScheme * fGeoScheme
Double_t fTimeDigiLast
Time of last digi sent to DAQ.
Int_t GenerateNoisePerModule(CbmMuchModuleGem *, Double_t, Double_t)
Double_t fTimeTot
Total execution time.
void AddNoiseSignal(CbmMuchPad *, Double_t, Double_t)
Bool_t ExecPoint(const CbmMuchPoint *point, Int_t)
Double_t fTimeDigiFirst
Time of first digi sent to DAQ.
virtual void Exec(Option_t *opt)
Bool_t BufferSignals(Int_t, Double_t, Double_t)
Bool_t AddCharge(CbmMuchSectorRadial *s, UInt_t ne, Int_t iPoint, Double_t time, Double_t driftTime, Double_t phi1, Double_t phi2)
Double_t GetNPrimaryElectronsPerCm(const CbmMuchPoint *point, int detectortype)
Int_t fNofNoiseTot
Set this boolean variable to True if want to generated Elecronics Noise.
Int_t fNofPoints
Number of points processed in Exec.
Double_t Sigma_n_e(Double_t Tkin, Double_t mass)
Bool_t fGenerateElectronicsNoise
Previous Event Time.
std::map< UInt_t, UInt_t > fAddressCharge
Function to sample the noise charge.
Double_t fNofSignalsTot
Total Number of signals.
Double_t fNofDigisTot
Total number of digis created.
Double_t fNofPointsTot
Total number of points processed.
Int_t GenerateNoise(Double_t, Double_t)
virtual InitStatus Init()
CbmMuchDigi * ConvertSignalToDigi(CbmMuchSignal *)
Int_t GetNModules() const
CbmMuchModule * GetModule(Int_t iModule) const
CbmMuchLayerSide * GetSide(Bool_t side)
CbmMuchPadRadial * GetPad(Double_t x, Double_t y)
CbmMuchSectorRadial * GetSectorByRadius(Double_t r)
CbmMuchSectorRectangular * GetSector(Double_t x, Double_t y)
CbmMuchPadRectangular * GetPad(Double_t x, Double_t y)
CbmMuchSector * GetSectorByIndex(Int_t iSector)
Int_t GetNSectors() const
Int_t GetDetectorId() const
Int_t GetDetectorType() const
Double_t GetPhi2() const
Double_t GetPhi1() const
Double_t GetX() const
Definition CbmMuchPad.h:32
Double_t GetDy() const
Definition CbmMuchPad.h:35
Int_t GetAddress() const
Definition CbmMuchPad.h:31
Double_t GetDx() const
Definition CbmMuchPad.h:34
Double_t GetY() const
Definition CbmMuchPad.h:33
void PositionIn(TVector3 &pos) const
void PositionOut(TVector3 &pos) const
static CbmMuchReadoutBuffer * Instance()
CbmMuchPadRadial * GetPadByPhi(Double_t phi)
Int_t GetNChannels() const
CbmMuchPad * GetPadByChannelIndex(Int_t iChannel) const
Data class for an analog signal in the MUCH Simple data class used in the digitisation process of the...
void AddNoise(UInt_t)
UInt_t GetCharge() const
Long64_t GetTimeStamp()
CbmMatch * GetMatch() const
Int_t GetAddress() const
void SetCharge(UInt_t charge)
CbmMuchLayer * GetLayer(Int_t iLayer) const
Int_t GetNLayers() const
Int_t ReadOutData(Double_t time, std::vector< Data * > &dataList)
void Fill(UInt_t address, Data *data)
virtual Int_t GetNData()