FairRoot/PandaRoot
Public Member Functions | Protected Attributes | List of all members
PndDiscTaskReconstruction Class Reference

#include <PndDiscTaskReconstruction.h>

Inheritance diagram for PndDiscTaskReconstruction:
PndPersistencyTask

Public Member Functions

 PndDiscTaskReconstruction ()
 
 PndDiscTaskReconstruction (const char *name)
 
virtual ~PndDiscTaskReconstruction ()
 
virtual InitStatus Init ()
 
virtual InitStatus ReInit ()
 
virtual void Exec (Option_t *opt)
 
virtual void FinishEvent ()
 
virtual void FinishTask ()
 
void RunTimeBased (Bool_t time_based)
 
void SetFlag_ExportPatterns (Bool_t flag)
 
void SetAvgWavelength (Double_t const &val)
 
Double_t const & GetAvgWavelength () const
 
void SetMinWavelength (Double_t const &val)
 
Double_t const & GetMinWavelength () const
 
double gauss (int x, double mean, double rms)
 
double mean (std::vector< double > values)
 
double deviation (double mean, std::vector< double > values)
 
void SetPersistency (Bool_t val=kTRUE)
 
Bool_t GetPersistency ()
 

Protected Attributes

TString branch_name_digits
 Branch name where digitized hits are stored. More...
 
TString folder_name_digits
 Folder name for digits. More...
 
TClonesArray * tclarr_digits
 
TClonesArray * tclarr_digits_out
 
TClonesArray * tclarr_particles_out
 
TClonesArray * tclarr_particles
 to cache the pointer to particle MC TClonesArray returned by IO manager. More...
 
TClonesArray * tclarr_tracks
 to cache the pointer to fitted tracks TClonesArray returned by IO manager. More...
 
TClonesArray * tclarr_recon_results
 results of reconstruction (pattern matching) More...
 
Bool_t is_time_based
 Time based buffering on/off. More...
 
Bool_t flag_export_patterns
 Write the pattern hypothesis and measured patterns to file. More...
 
StopTime start_functor
 
StopTime stop_functor
 
Double_t average_wl
 
Double_t minimum_wl
 

Detailed Description

Definition at line 23 of file PndDiscTaskReconstruction.h.

Constructor & Destructor Documentation

PndDiscTaskReconstruction::PndDiscTaskReconstruction ( )

Definition at line 410 of file PndDiscTaskReconstruction.cxx.

References branch_name_digits, and folder_name_digits.

410  : PndPersistencyTask("DiscDircTaskReconstruction"), tclarr_digits(NULL), tclarr_particles(NULL), tclarr_tracks(NULL), is_time_based(kTRUE), flag_export_patterns(kFALSE), average_wl(420.), minimum_wl(385.)
411 {
412  branch_name_digits = "DiscDigit";
413  folder_name_digits = "DiscDIRC";
414 }
TString branch_name_digits
Branch name where digitized hits are stored.
TString folder_name_digits
Folder name for digits.
TClonesArray * tclarr_particles
to cache the pointer to particle MC TClonesArray returned by IO manager.
TClonesArray * tclarr_tracks
to cache the pointer to fitted tracks TClonesArray returned by IO manager.
Bool_t flag_export_patterns
Write the pattern hypothesis and measured patterns to file.
Bool_t is_time_based
Time based buffering on/off.
PndDiscTaskReconstruction::PndDiscTaskReconstruction ( const char *  name)
PndDiscTaskReconstruction::~PndDiscTaskReconstruction ( )
virtual

Definition at line 418 of file PndDiscTaskReconstruction.cxx.

419 {
420 }

Member Function Documentation

double PndDiscTaskReconstruction::deviation ( double  mean,
std::vector< double >  values 
)

Definition at line 393 of file PndDiscTaskReconstruction.cxx.

References i, and sqrt().

394 {
395  double average = 0;
396 
397  for(size_t i = 0; i < values.size(); i++)
398  {
399  average += (values[i]-mean_value)*(values[i]-mean_value);
400  }
401 
402  average = sqrt(1/(double(values.size())-1)*average);
403 
404  return average;
405 }
Int_t i
Definition: run_full.C:25
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
void PndDiscTaskReconstruction::Exec ( Option_t *  opt)
virtual

Use particle track array as driving quantity and retrieve pattern digits by using time-based-simulation.

Definition at line 486 of file PndDiscTaskReconstruction.cxx.

References acos(), alpha, angle, angle_sensor, cos(), dx, dy, E0, PndDiscDigitizedHit::GetDetectorID(), PndDiscDigitizedHit::GetPixelNumber(), PndDiscDigitizedHit::GetReadoutID(), PndDiscDigitizedHit::GetTdcTime(), PndDiscReconResult::hypothesis, i, lambda(), n, p, PndDiscReconResult::particle, phi, Pi, PndDiscReconResult::pixel, pos_fel_x, pos_fel_y, pz, PndDiscReconResult::sensor, sin(), sqrt(), t, tclarr_digits, tclarr_digits_out, tclarr_particles, tclarr_particles_out, tclarr_recon_results, theta, PndDiscReconResult::time, x0, and y0.

487 {
488  //Reading out entries
489  int particles = tclarr_particles->GetEntries();
490  int entries = tclarr_digits->GetEntries();
491 
492  // time window around track time:
493  //const double time_window_t1 = -5.0; //[R.K. 01/2017] unused variable?
494  //const double time_window_t2 = 40.0; //[R.K. 01/2017] unused variable?
495 
496 
497  //Information of primary particle
498  std::vector<double> x0,y0,px,py,pz,p;
499  std::vector<double> theta,phi;
500  std::vector<double> t;
501 
502  std::vector<std::vector<double > > mean_time;
503  std::vector<std::vector<double > > deviation_time;
504 
505 
506 
507  //double average; //[R.K. 01/2017] unused variable?
508 
509  //-----------------------------------------------------------------------
510  // Reading out primary particle information
511  //-----------------------------------------------------------------------
512 
513  for(int i = 0; i < particles; i++)
514  {
516 
517  x0.push_back(particle_mc_point->GetX());
518  y0.push_back(particle_mc_point->GetY());
519 
520  px.push_back(particle_mc_point->GetPx());
521  py.push_back(particle_mc_point->GetPy());
522  pz.push_back(particle_mc_point->GetPz());
523 
524  p.push_back(sqrt(px[i]*px[i]+py[i]*py[i]+pz[i]*pz[i]));
525 
526  t.push_back(particle_mc_point->GetTime());
527 
528  theta.push_back(acos(pz[i]/sqrt(px[i]*px[i]+py[i]*py[i]+pz[i]*pz[i])));
529 
530  std::cout << "Primary particles:" << std::endl;
531  std::cout << "Particle pos: x = "<< x0[i] << std::endl;
532  std::cout << "Particle pos: y = "<< y0[i] << std::endl;
533  std::cout << "Particle mom: px = "<< px[i] << std::endl;
534  std::cout << "Particle mom: py = "<< py[i] << std::endl;
535  std::cout << "Particle mom: pz = "<< pz[i] << std::endl;
536  std::cout << "Particle abs mom: p = "<< p[i] << std::endl;
537  std::cout << "Particle time: t = "<< t[i] << std::endl;
538  std::cout << "Particle angle: theta[deg] = " << theta[i]*180/TMath::Pi() << std::endl;
539  std::cout << std::endl;
540  }
541 
542  //-------------------------------------------------------------------------
543  // Caculating theoretical hit pattern for pion, kaon, proton
544  //-------------------------------------------------------------------------
545 
546  PndDiscReconResult *recon_result = new PndDiscReconResult();
547 
548  for(int i = 0; i < particles; i++)
549  {
550  for(int j = 0; j < 3; j++)
551  {
552  std::cout << "Particle: " << j << std::endl;
553  for(int k = 0; k < 108; k++)
554  {
555  double E = sqrt(p[i]*1000*p[i]*1000 + E0[j]*E0[j]);
556  double gamma = E/E0[j];
557  double beta = sqrt(1-1/(gamma*gamma));
558  double n = 1.47;
559  double cherenkov = acos(1/(n*beta));
560 
561  double dx = pos_fel_x[k] - x0[i];
562  double dy = pos_fel_y[k] - y0[i];
563 
564  double phirel = acos((dx*px[i]+dy*py[i])/(sqrt(dx*dx+dy*dy)*sqrt(px[i]*px[i]+py[i]*py[i])));
565 
566  double angle;
567 
568  if(div(k,27).quot == 0 || div(k,27).quot == 3)
569  {
570  angle = angle_sensor[k]*TMath::Pi()/180;
571  }
572  if(div(k,27).quot == 1 || div(k,27).quot == 2)
573  {
574  angle = 2*TMath::Pi() - angle_sensor[k]*TMath::Pi()/180;
575  }
576 
577  double alpha = acos(dx/(sqrt(dx*dx+dy*dy))) - angle;
578 
579  double delta = sin(theta[i])*cos(phirel);
580  double lambda = delta*delta+cos(theta[i])*cos(theta[i]);
581  double varphi = acos(delta*cos(cherenkov)/lambda+sqrt((cos(theta[i])*cos(theta[i])-cos(cherenkov)*cos(cherenkov))/lambda+(delta*cos(cherenkov)/lambda)*(delta*cos(cherenkov)/lambda)));
582  double varphidash = atan(tan(varphi)/cos(alpha));
583 
584  int pixel = (TMath::Pi()/2 - varphidash-0.8530931655)/0.0035321867;
585 
586  double time = (sqrt(dx*dx+dy*dy))/100/3e8*n/cos(varphidash)*1e9;
587 
588  //Results are only valid for pixels between 0 and 100
589  if(pixel >= 0 && pixel <= 100)
590  {
591 
592  std::cout << "Sensor: " << k << " Pixel: " << pixel << " Time: " << time << std::endl;
593 
594  recon_result->hypothesis = j;
595  recon_result->sensor = k;
596  recon_result->pixel = pixel;
597  recon_result->time = time;
598  recon_result->particle = i;
599 
600  new ((*tclarr_recon_results)[tclarr_recon_results->GetEntriesFast()]) PndDiscReconResult(*recon_result);
601  }
602  else
603  {
604 
605 
606  recon_result->hypothesis = j;
607  recon_result->sensor = k;
608  recon_result->pixel = -1;
609  recon_result->time = time;
610  recon_result->particle = i;
611 
612  new ((*tclarr_recon_results)[tclarr_recon_results->GetEntriesFast()]) PndDiscReconResult(*recon_result);
613  }
614  }
615  std::cout << std::endl;
616  }
617  }
618 
619 
620 
621  //-------------------------------------------------------------------------
622  // Copying digitized hits and particle information
623  //-------------------------------------------------------------------------
624 
625  for(int i = 0; i < entries; i++)
626  {
628  digits = new((*tclarr_digits_out)[tclarr_digits_out->GetEntriesFast()]) PndDiscDigitizedHit(*digits);
629  }
630 
631  for(int i = 0; i < particles; i++)
632  {
634  particles_out = new((*tclarr_particles_out)[tclarr_particles_out->GetEntriesFast()]) PndDiscParticleMCPoint(*particles_out);
635  }
636 
637 
638 
639  //Reading out hit information and calculating Cherenkov angle
640  for(int i = 0; i < entries; i++)
641  {
642  PndDiscDigitizedHit* sensor_mc_point = (PndDiscDigitizedHit*)tclarr_digits->At(i);
643 
644  int detector_id = sensor_mc_point->GetDetectorID();
645  int sensor_id = 27*detector_id + sensor_mc_point->GetReadoutID(); //Sensor ID of hit
646  int pixel = sensor_mc_point->GetPixelNumber(); // Pixel number of hit
647 
648  double time = sensor_mc_point->GetTdcTime();
649 
650  //Difference vector between particle and sensor element
651  //double dx = pos_fel_x[sensor_id] - x0[0];
652  //double dy = pos_fel_y[sensor_id] - y0[0];
653 
654  //double angle = 0; //[R.K. 01/2017] unused variable
655 
656  //if(detector_id == 0 || detector_id == 3) //[R.K. 01/2017] unused variable
657  //{
658  //angle = angle_sensor[sensor_id]*TMath::Pi()/180; //[R.K. 01/2017] unused variable
659  //}
660  //if(detector_id == 1 || detector_id == 2) //[R.K. 01/2017] unused variable
661  //{
662  //angle = 2*TMath::Pi() - angle_sensor[sensor_id]*TMath::Pi()/180; //[R.K. 01/2017] unused variable
663  //}
664 
665  //double phirel = acos((dx*px[0]+dy*py[0])/(sqrt(dx*dx+dy*dy)*sqrt(px[0]*px[0]+py[0]*py[0])));
666  //double alpha = acos(dx/(sqrt(dx*dx+dy*dy))) - angle;
667  //double varphidash = TMath::Pi()/2-(0.00347196*pixel + 0.852934);
668  //double varphi = atan(tan(varphidash)*cos(alpha));
669 
670  //double cherenkov = acos(sin(theta[0])*cos(phirel)*cos(varphi)+cos(theta[0])*sin(varphi)); //Reconstructed Cherenkov angle
671 
672  std::cout << "Hit informations:" << std::endl;
673  std::cout << "TDC Time [ns]: " << time << std::endl;
674  std::cout << "Detector ID: " << detector_id << std::endl;
675  std::cout << "Sensor ID: " << sensor_id << std::endl;
676  std::cout << "Sensor coor: x = " << pos_fel_x[sensor_id] << ", y = " << pos_fel_y[sensor_id] << std::endl;
677  std::cout << "Pixel no: " << pixel << std::endl;
678  std::cout << std::endl;
679  }
680 
681 
682  std::cout << "DiscDIRC_TaskReconstruction::Exec\n";
683 }
friend F32vec4 acos(const F32vec4 &a)
Definition: P4_F32vec4.h:113
Double_t x0
Definition: checkhelixhit.C:70
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
double dy
Int_t i
Definition: run_full.C:25
const Int_t & GetReadoutID() const
Double_t lambda(Double_t x, Double_t y, Double_t z)
Definition: drawdal.C:48
double angle_sensor[108]
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:111
int n
const Double_t & GetTdcTime() const
Double_t p
Definition: anasim.C:58
Double_t y0
Definition: checkhelixhit.C:71
TClonesArray * tclarr_particles
to cache the pointer to particle MC TClonesArray returned by IO manager.
const Int_t & GetDetectorID() const
double pos_fel_y[108]
const Int_t & GetPixelNumber() const
double dx
double E0[3]
TClonesArray * tclarr_recon_results
results of reconstruction (pattern matching)
TTree * t
Definition: bump_analys.C:13
double alpha
Definition: f_Init.h:9
Double_t angle
double pos_fel_x[108]
Double_t Pi
double pz[39]
Definition: pipisigmas.h:14
void PndDiscTaskReconstruction::FinishEvent ( )
virtual

Definition at line 686 of file PndDiscTaskReconstruction.cxx.

References tclarr_digits_out, tclarr_particles_out, and tclarr_recon_results.

687 {
688  // called after all Tasks did their Exec() and the data is copied to the file
689  if(tclarr_recon_results != NULL) tclarr_recon_results->Clear();
690  if(tclarr_digits_out != NULL) tclarr_digits_out->Clear();
691  if(tclarr_particles_out != NULL) tclarr_particles_out->Clear();
692 
693  FinishEvents();
694 }
TClonesArray * tclarr_recon_results
results of reconstruction (pattern matching)
void PndDiscTaskReconstruction::FinishTask ( )
virtual

Definition at line 697 of file PndDiscTaskReconstruction.cxx.

698 {
699 }
double PndDiscTaskReconstruction::gauss ( int  x,
double  mean,
double  rms 
)

Definition at line 371 of file PndDiscTaskReconstruction.cxx.

References exp(), Pi, and sqrt().

372 {
373  return 1/sqrt(2*TMath::Pi()*rms*rms)*exp(-0.5*((x-mean_value)/rms)*((x-mean_value)/rms));
374 }
friend F32vec4 exp(const F32vec4 &a)
Definition: P4_F32vec4.h:109
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
Double_t x
Double_t Pi
Double_t const& PndDiscTaskReconstruction::GetAvgWavelength ( ) const
inline

Definition at line 41 of file PndDiscTaskReconstruction.h.

References average_wl.

Double_t const& PndDiscTaskReconstruction::GetMinWavelength ( ) const
inline

Definition at line 44 of file PndDiscTaskReconstruction.h.

References minimum_wl.

Bool_t PndPersistencyTask::GetPersistency ( )
inlineinherited

Definition at line 32 of file PndPersistencyTask.h.

References PndPersistencyTask::fPersistency.

Referenced by PndLmdPixelHitProducerFast::GetPersistance(), PndMdtDigitization::Init(), PndMdtHitProducerIdeal::Init(), PndMdtClusterTask::Init(), PndFtsHitProducerRealFast::Init(), PndRichHitProducer::Init(), PndSttHitProducerRealFast::Init(), Init(), PndSttHelixHitProducer::Init(), PndDiscTaskPID::Init(), PndIdealTrackFinder::Init(), PndSttMvdGemTracking::Init(), PndMdtTrkProducer::Init(), PndFtsHitProducerRealFull::Init(), PndLmdPixelClusterTask::Init(), PndSttHitProducerRealFull::Init(), PndLmdStripClusterTask::Init(), PndEmcApdHitProducer::Init(), PndMissingPzCleanerTask::Init(), PndEmcMakeRecoHit::Init(), PndEmcMakeClusterOnline::Init(), PndTrackSmearTask::Init(), PndEmcFWEndcapTimebasedWaveforms::Init(), PndSttHitProducerIdeal::Init(), PndEmcFWEndcapDigi::Init(), PndFtsHitProducerIdeal::Init(), PndEmcMakeCluster::Init(), PndMdtPointsToWaveform::Init(), PndDiscTaskDigitization::Init(), PndEmcMakeDigi::Init(), PndSdsTimeWalkCorrTask::Init(), PndLmdPixelHitProducerFast::Init(), PndDrcHitFinder::Init(), PndRichHitFinder::Init(), PndEmcMakeCorr::Init(), PndFtofHitProducerIdeal::Init(), PndEmcHitsToWaveform::Init(), PndSciTDigiTask::Init(), PndDrcHitProducerIdeal::Init(), PndSdsHitProducerIdeal::Init(), PndSciTHitProducerIdeal::Init(), PndRecoMultiKalmanTask2::Init(), PndEmcHitProducer::Init(), PndDrcHitProducerReal::Init(), PndDskFLGHitProducerIdeal::Init(), PndEmcTmpWaveformToDigi::Init(), PndDrcDigiTask::Init(), PndEmcWaveformToDigi::Init(), PndSttMatchTracks::Init(), PndEmcWaveformToCalibratedDigi::Init(), PndTrkTracking2::Init(), PndSttFindTracks::Init(), PndEmcMultiWaveformToCalibratedDigi::Init(), PndRecoKalmanTask2::Init(), PndDrcTimeDigiTask::Init(), PndEmcExpClusterSplitter::Init(), PndFtsHoughTrackerTask::Init(), PndSdsNoiseProducer::Init(), PndEmcPhiBumpSplitter::Init(), PndSdsIdealRecoTask::Init(), PndSdsHybridHitProducer::Init(), PndRecoMultiKalmanTask::Init(), PndSdsIdealClusterTask::Init(), PndRecoKalmanTask::Init(), PndSdsStripHitProducerDif::Init(), PndGemDigitize::Init(), PndSdsStripHitProducer::Init(), PndGemFindHits::Init(), PndSdsPixelClusterTask::Init(), PndSdsStripClusterTask::Init(), PndMvdGemTrackFinderOnHits::Init(), PndBarrelTrackFinder::Init(), PndEmcFullDigiTask::PndEmcFullDigiTask(), PndEmcMakeBump::PndEmcMakeBump(), PndUnassignedHitsTask::RegisterBranches(), PndMvdClusterTask::SetPersistance(), PndMvdDigiTask::SetPersistance(), PndEmcMakeBump::SetStorageOfData(), and PndEmcFullDigiTask::StoreDigi().

32 { return fPersistency; }
InitStatus PndDiscTaskReconstruction::Init ( )
virtual

Definition at line 435 of file PndDiscTaskReconstruction.cxx.

References branch_name_digits, PndPersistencyTask::GetPersistency(), tclarr_digits, tclarr_digits_out, tclarr_particles, tclarr_particles_out, and tclarr_recon_results.

436 {
437  LOG(INFO) << "PndDiscTaskReconstruction::Init()";
438 
439  // Get IO manager instance:
440  FairRootManager* io_manager = FairRootManager::Instance();
441  if(!io_manager)
442  {
443  LOG(FATAL) << "FairRootManager instance is NULL !!!";
444  return kFATAL;
445  }
446 
447  // Get branch with digitized hits:
448  tclarr_digits = (TClonesArray*) io_manager->GetObject(branch_name_digits);
449  if(!tclarr_digits)
450  {
451  LOG(ERROR) << "Branch " << branch_name_digits.Data() << " is not accessible through FairRootManager.";
452  return kERROR;
453  }
454 
455  tclarr_particles = (TClonesArray*) io_manager->GetObject("DiscMCTruthTracks");
456  if(!tclarr_particles)
457  {
458  LOG(ERROR) << "No DiscDIRC_ParticleMCPoint collection registered with FairRootManager";
459  return kFATAL;
460  }
461 
462 
463  tclarr_recon_results = new TClonesArray("PndDiscReconResult");
464  tclarr_digits_out = new TClonesArray("PndDiscDigitizedHit");
465  tclarr_particles_out = new TClonesArray("PndDiscParticleMCPoint");
466 
467  io_manager->Register("DiscPatternPrediction", "DiscDircDetector", tclarr_recon_results, GetPersistency());
468  io_manager->Register("DiscDigitizedHit","DiscDircDetector", tclarr_digits_out, GetPersistency());
469  io_manager->Register("DiscRealTracks","DiscDircDetector", tclarr_particles_out, GetPersistency());
470 
471 
472 
473  static bool recon_initialized = false;
474  if(recon_initialized) return kSUCCESS;
475 
476  return kSUCCESS;
477 }
TString branch_name_digits
Branch name where digitized hits are stored.
TClonesArray * tclarr_particles
to cache the pointer to particle MC TClonesArray returned by IO manager.
TClonesArray * tclarr_recon_results
results of reconstruction (pattern matching)
double PndDiscTaskReconstruction::mean ( std::vector< double >  values)

Definition at line 378 of file PndDiscTaskReconstruction.cxx.

References i.

379 {
380  double average = 0;
381 
382  for(size_t i = 0; i < values.size(); i++)
383  {
384  average += values[i];
385  }
386 
387  return average/double(values.size());
388 }
Int_t i
Definition: run_full.C:25
InitStatus PndDiscTaskReconstruction::ReInit ( )
virtual

Definition at line 426 of file PndDiscTaskReconstruction.cxx.

427 {
428  LOG(INFO) << "DiscDircTaskReconstruction::ReInit()";
429  return kSUCCESS;
430 }
void PndDiscTaskReconstruction::RunTimeBased ( Bool_t  time_based)
inline

Definition at line 37 of file PndDiscTaskReconstruction.h.

References is_time_based.

37 { is_time_based = time_based; }
Bool_t is_time_based
Time based buffering on/off.
void PndDiscTaskReconstruction::SetAvgWavelength ( Double_t const &  val)
inline

Definition at line 40 of file PndDiscTaskReconstruction.h.

References average_wl, and val.

40 {average_wl = val; }
Double_t val[nBoxes][nFEBox]
Definition: createCalib.C:11
void PndDiscTaskReconstruction::SetFlag_ExportPatterns ( Bool_t  flag)
inline

Definition at line 38 of file PndDiscTaskReconstruction.h.

References flag_export_patterns.

38 {flag_export_patterns = flag; }
Bool_t flag_export_patterns
Write the pattern hypothesis and measured patterns to file.
void PndDiscTaskReconstruction::SetMinWavelength ( Double_t const &  val)
inline

Definition at line 43 of file PndDiscTaskReconstruction.h.

References minimum_wl, and val.

43 {minimum_wl = val; }
Double_t val[nBoxes][nFEBox]
Definition: createCalib.C:11
void PndPersistencyTask::SetPersistency ( Bool_t  val = kTRUE)
inlineinherited

Definition at line 31 of file PndPersistencyTask.h.

References PndPersistencyTask::fPersistency, and val.

Referenced by barrelTrackFinder(), digi_complete(), digi_complete_newSTT(), digiOnly_complete(), PndBarrelTrackFinder::PndBarrelTrackFinder(), PndCATracking::PndCATracking(), PndDrcHitFinder::PndDrcHitFinder(), PndEmc2DLocMaxFinder::PndEmc2DLocMaxFinder(), PndEmcExpClusterSplitter::PndEmcExpClusterSplitter(), PndEmcFullDigiTask::PndEmcFullDigiTask(), PndEmcFWEndcapDigi::PndEmcFWEndcapDigi(), PndEmcFWEndcapTimebasedWaveforms::PndEmcFWEndcapTimebasedWaveforms(), PndEmcHitProducer::PndEmcHitProducer(), PndEmcHitsToWaveform::PndEmcHitsToWaveform(), PndEmcMakeBump::PndEmcMakeBump(), PndEmcMakeCluster::PndEmcMakeCluster(), PndEmcMakeClusterOnline::PndEmcMakeClusterOnline(), PndEmcMakeDigi::PndEmcMakeDigi(), PndEmcMakeRecoHit::PndEmcMakeRecoHit(), PndEmcMultiWaveformToCalibratedDigi::PndEmcMultiWaveformToCalibratedDigi(), PndEmcPhiBumpSplitter::PndEmcPhiBumpSplitter(), PndEmcTmpWaveformToDigi::PndEmcTmpWaveformToDigi(), PndEmcWaveformToCalibratedDigi::PndEmcWaveformToCalibratedDigi(), PndEmcWaveformToDigi::PndEmcWaveformToDigi(), PndFtofHitProducerIdeal::PndFtofHitProducerIdeal(), PndFtsCATracking::PndFtsCATracking(), PndFtsHitProducerIdeal::PndFtsHitProducerIdeal(), PndFtsHitProducerRealFast::PndFtsHitProducerRealFast(), PndFtsHitProducerRealFull::PndFtsHitProducerRealFull(), PndFtsHoughTrackerTask::PndFtsHoughTrackerTask(), PndGemDigitize::PndGemDigitize(), PndGemFindHits::PndGemFindHits(), PndIdealTrackFinder::PndIdealTrackFinder(), PndLmdPixelClusterTask::PndLmdPixelClusterTask(), PndLmdPixelHitProducerFast::PndLmdPixelHitProducerFast(), PndMdtClusterTask::PndMdtClusterTask(), PndMdtDigitization::PndMdtDigitization(), PndMdtHitProducerIdeal::PndMdtHitProducerIdeal(), PndMdtPointsToWaveform::PndMdtPointsToWaveform(), PndMdtTrkProducer::PndMdtTrkProducer(), PndMissingPzCleanerTask::PndMissingPzCleanerTask(), PndMvdGemTrackFinderOnHits::PndMvdGemTrackFinderOnHits(), PndMvdHitProducerIdeal::PndMvdHitProducerIdeal(), PndMvdPixelClusterTask::PndMvdPixelClusterTask(), PndMvdTimeWalkCorrTask::PndMvdTimeWalkCorrTask(), PndMvdToPix4ClusterTask::PndMvdToPix4ClusterTask(), PndRecoKalmanTask::PndRecoKalmanTask(), PndRecoKalmanTask2::PndRecoKalmanTask2(), PndRecoMultiKalmanTask::PndRecoMultiKalmanTask(), PndRecoMultiKalmanTask2::PndRecoMultiKalmanTask2(), PndRichHitFinder::PndRichHitFinder(), PndRichHitProducer::PndRichHitProducer(), PndSciTDigiTask::PndSciTDigiTask(), PndSciTHitProducerIdeal::PndSciTHitProducerIdeal(), PndSdsHitProducerIdeal::PndSdsHitProducerIdeal(), PndSdsHybridHitProducer::PndSdsHybridHitProducer(), PndSdsIdealClusterTask::PndSdsIdealClusterTask(), PndSdsIdealRecoTask::PndSdsIdealRecoTask(), PndSdsNoiseProducer::PndSdsNoiseProducer(), PndSdsPixelClusterTask::PndSdsPixelClusterTask(), PndSdsStripClusterTask::PndSdsStripClusterTask(), PndSdsStripHitProducer::PndSdsStripHitProducer(), PndSdsTimeWalkCorrTask::PndSdsTimeWalkCorrTask(), PndSttFindTracks::PndSttFindTracks(), PndSttHelixHitProducer::PndSttHelixHitProducer(), PndSttHitProducerIdeal::PndSttHitProducerIdeal(), PndSttHitProducerRealFast::PndSttHitProducerRealFast(), PndSttHitProducerRealFull::PndSttHitProducerRealFull(), PndSttMatchTracks::PndSttMatchTracks(), PndSttMvdGemTracking::PndSttMvdGemTracking(), PndTrackSmearTask::PndTrackSmearTask(), PndTrkTracking2::PndTrkTracking2(), reco(), reco_complete(), reco_complete_gf2(), reco_complete_newSTT(), reco_complete_sec(), recoideal_complete(), PndMvdClusterTask::SetPersistance(), PndMvdDigiTask::SetPersistance(), PndLmdPixelHitProducerFast::SetPersistance(), PndSdsHitProducerIdeal::SetPersistance(), PndSttMvdGemTracking::SetPersistenc(), PndMdtClusterTask::SetPersistence(), PndSttHelixHitProducer::SetPersistence(), PndMissingPzCleanerTask::SetPersistence(), PndFtsHitProducerRealFast::SetPersistence(), PndFtsHitProducerRealFull::SetPersistence(), PndSttHitProducerRealFull::SetPersistence(), PndSttHitProducerIdeal::SetPersistence(), PndSttHitProducerRealFast::SetPersistence(), PndFtsHitProducerIdeal::SetPersistence(), PndTrackSmearTask::SetPersistence(), PndSciTHitProducerIdeal::SetPersistence(), PndIdealTrackFinder::SetPersistence(), PndSttMatchTracks::SetPersistence(), PndSttFindTracks::SetPersistence(), PndFtsHoughTrackerTask::SetPersistence(), PndTrkTracking2::SetPersistence(), PndEmcMakeRecoHit::SetStorageOfData(), PndEmcFWEndcapDigi::SetStorageOfData(), PndEmcMakeClusterOnline::SetStorageOfData(), PndEmcFWEndcapTimebasedWaveforms::SetStorageOfData(), PndEmcMakeDigi::SetStorageOfData(), PndMdtPointsToWaveform::SetStorageOfData(), PndEmc2DLocMaxFinder::SetStorageOfData(), PndEmcMakeCluster::SetStorageOfData(), PndEmcHitsToWaveform::SetStorageOfData(), PndEmcMakeBump::SetStorageOfData(), PndEmcTmpWaveformToDigi::SetStorageOfData(), PndEmcWaveformToDigi::SetStorageOfData(), PndEmcWaveformToCalibratedDigi::SetStorageOfData(), PndEmcMultiWaveformToCalibratedDigi::SetStorageOfData(), PndEmcExpClusterSplitter::SetStorageOfData(), PndEmcPhiBumpSplitter::SetStorageOfData(), standard_tracking(), and PndEmcFullDigiTask::StoreDigi().

31 { fPersistency = val; }
Double_t val[nBoxes][nFEBox]
Definition: createCalib.C:11

Member Data Documentation

Double_t PndDiscTaskReconstruction::average_wl
protected

Definition at line 69 of file PndDiscTaskReconstruction.h.

Referenced by GetAvgWavelength(), and SetAvgWavelength().

TString PndDiscTaskReconstruction::branch_name_digits
protected

Branch name where digitized hits are stored.

Definition at line 51 of file PndDiscTaskReconstruction.h.

Referenced by Init(), and PndDiscTaskReconstruction().

Bool_t PndDiscTaskReconstruction::flag_export_patterns
protected

Write the pattern hypothesis and measured patterns to file.

Definition at line 64 of file PndDiscTaskReconstruction.h.

Referenced by SetFlag_ExportPatterns().

TString PndDiscTaskReconstruction::folder_name_digits
protected

Folder name for digits.

Definition at line 52 of file PndDiscTaskReconstruction.h.

Referenced by PndDiscTaskReconstruction().

Bool_t PndDiscTaskReconstruction::is_time_based
protected

Time based buffering on/off.

Definition at line 62 of file PndDiscTaskReconstruction.h.

Referenced by RunTimeBased().

Double_t PndDiscTaskReconstruction::minimum_wl
protected

Definition at line 70 of file PndDiscTaskReconstruction.h.

Referenced by GetMinWavelength(), and SetMinWavelength().

StopTime PndDiscTaskReconstruction::start_functor
protected

Definition at line 66 of file PndDiscTaskReconstruction.h.

StopTime PndDiscTaskReconstruction::stop_functor
protected

Definition at line 67 of file PndDiscTaskReconstruction.h.

TClonesArray* PndDiscTaskReconstruction::tclarr_digits
protected

Definition at line 54 of file PndDiscTaskReconstruction.h.

Referenced by Exec(), and Init().

TClonesArray* PndDiscTaskReconstruction::tclarr_digits_out
protected

Definition at line 55 of file PndDiscTaskReconstruction.h.

Referenced by Exec(), FinishEvent(), and Init().

TClonesArray* PndDiscTaskReconstruction::tclarr_particles
protected

to cache the pointer to particle MC TClonesArray returned by IO manager.

Definition at line 57 of file PndDiscTaskReconstruction.h.

Referenced by Exec(), and Init().

TClonesArray* PndDiscTaskReconstruction::tclarr_particles_out
protected

Definition at line 56 of file PndDiscTaskReconstruction.h.

Referenced by Exec(), FinishEvent(), and Init().

TClonesArray* PndDiscTaskReconstruction::tclarr_recon_results
protected

results of reconstruction (pattern matching)

Definition at line 59 of file PndDiscTaskReconstruction.h.

Referenced by Exec(), FinishEvent(), and Init().

TClonesArray* PndDiscTaskReconstruction::tclarr_tracks
protected

to cache the pointer to fitted tracks TClonesArray returned by IO manager.

Definition at line 58 of file PndDiscTaskReconstruction.h.


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