FairRoot/PandaRoot
Public Member Functions | Private Attributes | List of all members
PndDisc Class Reference

#include <PndDisc.h>

Inheritance diagram for PndDisc:

Public Member Functions

 PndDisc ()
 
 PndDisc (const char *name, Bool_t active, Int_t det_id=0)
 
 ~PndDisc ()
 
virtual void Initialize ()
 Initialize will be called after the Geometry is created. More...
 
virtual void ConstructGeometry ()
 
virtual void ConstructOpGeometry ()
 
virtual Bool_t ProcessHits (FairVolume *v=0)
 
virtual void BeginEvent ()
 
virtual void EndOfEvent ()
 Interface implementation - EndofEvent. More...
 
virtual void PreTrack ()
 
virtual void PostTrack ()
 
virtual bool CheckIfSensitive (std::string name)
 Interface implementation - Check if a given volume name belongs to a sensitive volume. More...
 
virtual void Register ()
 Interface implementation - Register the data collections in FairRootManager. More...
 
virtual TClonesArray * GetCollection (Int_t iColl) const
 Interface implementation - Get a data collection by index. More...
 
virtual void Reset ()
 Interface implementation - reset the 'containers'. More...
 
void StorePhotonTracks (Bool_t bval)
 Enable/Disable the storage of photon track information. More...
 
void SetFilterInterval (Double_t const &wl_min_nm_, Double_t const &wl_max_nm_)
 Set the wavelength range of the bandpass filters. More...
 

Private Attributes

Bool_t store_photon_tracks
 
Bool_t track_is_photon
 
TLorentzVector old_momentum
 
TClonesArray * clarr_sensor_hits
 
TClonesArray * clarr_photon_tracks
 hit on a photodetector surface More...
 
TClonesArray * clarr_particle_tracks
 optical photon tracks More...
 
int nextid_clarr_sensor_hits
 tracks of other particles More...
 
int nextid_clarr_photon_tracks
 
int nextid_clarr_particle_tracks
 
FairMCEventHeader * ev_header
 
std::map< int, std::pair< int,
double > > 
last_track_occurence
 
std::set< std::string > names_of_sensitive_volumes
 
std::map< int, double > internal_reflection_angle_of_photons
 
std::map< int, std::pair
< TLorentzVector,
TLorentzVector > > 
photons_entering_optics
 
int design_id
 photon track ids that have been More...
 
Double_t wl_min_nm
 
Double_t wl_max_nm
 

Detailed Description

Definition at line 35 of file PndDisc.h.

Constructor & Destructor Documentation

PndDisc::PndDisc ( )

Definition at line 49 of file PndDisc.cxx.

49  : FairDetector(), clarr_sensor_hits(NULL), clarr_particle_tracks(NULL),
51  ev_header(NULL), design_id(0), wl_min_nm(385.0), wl_max_nm(430.0)
52 {
53 }
FairMCEventHeader * ev_header
Definition: PndDisc.h:91
TClonesArray * clarr_sensor_hits
Definition: PndDisc.h:81
int design_id
photon track ids that have been
Definition: PndDisc.h:97
Double_t wl_max_nm
Definition: PndDisc.h:99
TClonesArray * clarr_particle_tracks
optical photon tracks
Definition: PndDisc.h:83
int nextid_clarr_particle_tracks
Definition: PndDisc.h:89
Double_t wl_min_nm
Definition: PndDisc.h:98
int nextid_clarr_sensor_hits
tracks of other particles
Definition: PndDisc.h:87
PndDisc::PndDisc ( const char *  name,
Bool_t  active,
Int_t  det_id = 0 
)

Definition at line 57 of file PndDisc.cxx.

58  : FairDetector(name, active, det_id), store_photon_tracks(kFALSE),
61  ev_header(NULL), design_id(0), wl_min_nm(385.0), wl_max_nm(430.0)
62 {
63 }
FairMCEventHeader * ev_header
Definition: PndDisc.h:91
TClonesArray * clarr_sensor_hits
Definition: PndDisc.h:81
int design_id
photon track ids that have been
Definition: PndDisc.h:97
Double_t wl_max_nm
Definition: PndDisc.h:99
TClonesArray * clarr_particle_tracks
optical photon tracks
Definition: PndDisc.h:83
int nextid_clarr_particle_tracks
Definition: PndDisc.h:89
Double_t wl_min_nm
Definition: PndDisc.h:98
int nextid_clarr_sensor_hits
tracks of other particles
Definition: PndDisc.h:87
TString name
Bool_t store_photon_tracks
Definition: PndDisc.h:77
PndDisc::~PndDisc ( )

Definition at line 67 of file PndDisc.cxx.

References clarr_particle_tracks, and clarr_sensor_hits.

67  {
68  if(clarr_sensor_hits!=NULL) {
69  clarr_sensor_hits->Delete();
70  delete clarr_sensor_hits;
71  }
72  if(clarr_particle_tracks!=NULL) {
73  clarr_particle_tracks->Delete();
74  delete clarr_particle_tracks;
75  }
76 }
TClonesArray * clarr_sensor_hits
Definition: PndDisc.h:81
TClonesArray * clarr_particle_tracks
optical photon tracks
Definition: PndDisc.h:83

Member Function Documentation

void PndDisc::BeginEvent ( )
virtual

Definition at line 539 of file PndDisc.cxx.

References ev_header, internal_reflection_angle_of_photons, last_track_occurence, and photons_entering_optics.

540 {
541  ev_header = FairMCApplication::Instance()->GetGenerator()->GetEvent();
543  photons_entering_optics.clear();
544  last_track_occurence.clear();
545 }
FairMCEventHeader * ev_header
Definition: PndDisc.h:91
std::map< int, double > internal_reflection_angle_of_photons
Definition: PndDisc.h:94
std::map< int, std::pair< int, double > > last_track_occurence
Definition: PndDisc.h:92
std::map< int, std::pair< TLorentzVector, TLorentzVector > > photons_entering_optics
Definition: PndDisc.h:95
bool PndDisc::CheckIfSensitive ( std::string  name)
virtual

Interface implementation - Check if a given volume name belongs to a sensitive volume.

FAIR geometry construction (FairModule::ExpandNode) querys sensitivity for every single volume name. Volume name format: Name_CopyNr (e.g. Radiator_0). For simplicity, copy numbers are neglected and all instances ('copys') of a specific volume will be flagged sensitive/insensitive.

Definition at line 613 of file PndDisc.cxx.

References names_of_sensitive_volumes.

614 {
615  return (names_of_sensitive_volumes.count(name)!=0);
616 }
std::set< std::string > names_of_sensitive_volumes
Definition: PndDisc.h:93
TString name
void PndDisc::ConstructGeometry ( )
virtual

Definition at line 102 of file PndDisc.cxx.

References names_of_sensitive_volumes, and TString.

103 {
104 
105  TString fname = GetGeometryFileName();
106  if(!fname.EndsWith(".root"))
107  {
108  LOG(ERROR) << "Geometry format not supported.";
109  return;
110  }
111 
112  if(!fname.CompareTo("DIRC_GEO_LIF1.root"))
113  {
114  LOG(WARNING) << "Wrong filename: (%s). DiscDIRC_Detector is expecting geo file DIRC_GEO_LRD.root." << fname.Data();
115  }
116 
117  names_of_sensitive_volumes.insert(std::string("DiscDIRC_Radiator"));
118  names_of_sensitive_volumes.insert(std::string("DiscDIRC_prism_a"));
119  names_of_sensitive_volumes.insert(std::string("DiscDIRC_prism_b"));
120  names_of_sensitive_volumes.insert(std::string("DiscDIRC_prism_c"));
121  names_of_sensitive_volumes.insert(std::string("DiscDIRC_focel_a"));
122  names_of_sensitive_volumes.insert(std::string("DiscDIRC_focel_b"));
123  names_of_sensitive_volumes.insert(std::string("DiscDIRC_focel_c"));
124  //names_of_sensitive_volumes.insert(std::string("DiscDIRC_Filter"));
125  names_of_sensitive_volumes.insert(std::string("DiscDIRC_sensor_pmt"));
126 
127  // for debugging:
128  names_of_sensitive_volumes.insert(std::string("container_volume"));
129 
130  ConstructRootGeometry();
131 }
std::set< std::string > names_of_sensitive_volumes
Definition: PndDisc.h:93
void PndDisc::ConstructOpGeometry ( )
virtual

Definition at line 135 of file PndDisc.cxx.

References Double_t, i, wl_max_nm, and wl_min_nm.

136 {
137  int i=0;
138  char str_name[100];
139  LOG(INFO) << "DiscDIRC_Detector::ConstructOpGeometry()";
140 
141  //-------------------------------------------------------
142  // define the surface properties
143  //-------------------------------------------------------
144 
145  gMC->DefineOpSurface("DiscDIRC_focel_mirror_surf", kGlisur, kDielectric_metal, kPolished, 0.0);
146  {
147  Double_t e_photon[] = {1.0e-09, 7.0e-09};
148  Double_t r_coeff[] = {0.9, 0.9};
149  Double_t efficiency[] = {1., 1.};
150  gMC->SetMaterialProperty("DiscDIRC_focel_mirror_surf", "REFLECTIVITY", 2, e_photon, r_coeff);
151  gMC->SetMaterialProperty("DiscDIRC_focel_mirror_surf", "EFFICIENCY", 2, e_photon, efficiency);
152  }
153 
154  gMC->DefineOpSurface("DiscDIRC_filter_a_surf", kGlisur, kDielectric_dielectric, kPolished, 0.0);
155  {
156  Double_t band_a_min = 1.239841939*1.e-06/wl_max_nm;
157  Double_t band_a_max = 1.239841939*1.e-06/wl_min_nm;
158  Double_t e_photon_a[] = {1.0e-09, band_a_min, band_a_min+1E-11, band_a_max, band_a_max+1E-11, 7.0e-09};
159  Double_t transmittance[] = {1E-10, 1E-10, 1., 1., 1E-10, 1E-10};
160  gMC->SetMaterialProperty("DiscDIRC_filter_a_surf", "TRANSMITTANCE", 6, e_photon_a, transmittance);
161  }
162 
163  //-------------------------------------------------------
164  // attach surfaces to the geometry setup
165  //-------------------------------------------------------
166 
167  gMC->SetBorderSurface("rad_interface_0_1","DiscDIRC_Radiator",0,"DiscDIRC_Radiator",1,"DiscDIRC_focel_mirror_surf");
168  gMC->SetBorderSurface("rad_interface_1_0","DiscDIRC_Radiator",1,"DiscDIRC_Radiator",0,"DiscDIRC_focel_mirror_surf");
169 
170  gMC->SetBorderSurface("rad_interface_1_2","DiscDIRC_Radiator",1,"DiscDIRC_Radiator",2,"DiscDIRC_focel_mirror_surf");
171  gMC->SetBorderSurface("rad_interface_2_1","DiscDIRC_Radiator",2,"DiscDIRC_Radiator",1,"DiscDIRC_focel_mirror_surf");
172 
173  gMC->SetBorderSurface("rad_interface_2_3","DiscDIRC_Radiator",2,"DiscDIRC_Radiator",3,"DiscDIRC_focel_mirror_surf");
174  gMC->SetBorderSurface("rad_interface_3_2","DiscDIRC_Radiator",3,"DiscDIRC_Radiator",2,"DiscDIRC_focel_mirror_surf");
175 
176  gMC->SetBorderSurface("rad_interface_3_0","DiscDIRC_Radiator",3,"DiscDIRC_Radiator",0,"DiscDIRC_focel_mirror_surf");
177  gMC->SetBorderSurface("rad_interface_0_3","DiscDIRC_Radiator",0,"DiscDIRC_Radiator",3,"DiscDIRC_focel_mirror_surf");
178 
179  gMC->SetBorderSurface("rad_hole_interface_0","DiscDIRC_Radiator",0,"DiscDIRC_hole",0,"DiscDIRC_focel_mirror_surf");
180  gMC->SetBorderSurface("rad_hole_interface_1","DiscDIRC_Radiator",1,"DiscDIRC_hole",0,"DiscDIRC_focel_mirror_surf");
181  gMC->SetBorderSurface("rad_hole_interface_2","DiscDIRC_Radiator",2,"DiscDIRC_hole",0,"DiscDIRC_focel_mirror_surf");
182  gMC->SetBorderSurface("rad_hole_interface_3","DiscDIRC_Radiator",3,"DiscDIRC_hole",0,"DiscDIRC_focel_mirror_surf");
183 
184  int n_fel = 27;
185  for(i=0; i<4*n_fel; i++)
186  {
187  // filters:
188  sprintf(str_name, "DiscDIRC_focel_a_window_a_%d", i);
189  gMC->SetBorderSurface(str_name, "DiscDIRC_sensor_window", i, "DiscDIRC_focel_a", i, "DiscDIRC_filter_a_surf");
190  sprintf(str_name, "DiscDIRC_focel_a_window_b_%d", i);
191  gMC->SetBorderSurface(str_name, "DiscDIRC_focel_a", i, "DiscDIRC_sensor_window", i, "DiscDIRC_filter_a_surf");
192 
193  sprintf(str_name, "DiscDIRC_focel_b_window_a_%d", i);
194  gMC->SetBorderSurface(str_name, "DiscDIRC_sensor_window", i, "DiscDIRC_focel_b", i, "DiscDIRC_filter_a_surf");
195  sprintf(str_name, "DiscDIRC_focel_b_window_b_%d", i);
196  gMC->SetBorderSurface(str_name, "DiscDIRC_focel_b", i, "DiscDIRC_sensor_window", i, "DiscDIRC_filter_a_surf");
197 
198  sprintf(str_name, "DiscDIRC_focel_c_window_a_%d", i);
199  gMC->SetBorderSurface(str_name, "DiscDIRC_sensor_window", i, "DiscDIRC_focel_c", i, "DiscDIRC_filter_a_surf");
200  sprintf(str_name, "DiscDIRC_focel_c_window_b_%d", i);
201  gMC->SetBorderSurface(str_name, "DiscDIRC_focel_c", i, "DiscDIRC_sensor_window", i, "DiscDIRC_filter_a_surf");
202 
203  // mirror coating
204  sprintf(str_name, "DiscDIRC_focel_a_mirror_a_%d", i);
205  gMC->SetBorderSurface(str_name, "DiscDIRC_focel_mirror_a", i, "DiscDIRC_focel_a", i, "DiscDIRC_focel_mirror_surf");
206  sprintf(str_name, "DiscDIRC_focel_a_mirror_b_%d", i);
207  gMC->SetBorderSurface(str_name, "DiscDIRC_focel_a", i, "DiscDIRC_focel_mirror_a", i, "DiscDIRC_focel_mirror_surf");
208 
209  sprintf(str_name, "DiscDIRC_focel_b_mirror_a_%d", i);
210  gMC->SetBorderSurface(str_name, "DiscDIRC_focel_mirror_b", i, "DiscDIRC_focel_b", i, "DiscDIRC_focel_mirror_surf");
211  sprintf(str_name, "DiscDIRC_focel_b_mirror_b_%d", i);
212  gMC->SetBorderSurface(str_name, "DiscDIRC_focel_b", i, "DiscDIRC_focel_mirror_b", i, "DiscDIRC_focel_mirror_surf");
213 
214  sprintf(str_name, "DiscDIRC_focel_c_mirror_a_%d", i);
215  gMC->SetBorderSurface(str_name, "DiscDIRC_focel_mirror_c", i, "DiscDIRC_focel_c", i, "DiscDIRC_focel_mirror_surf");
216  sprintf(str_name, "DiscDIRC_focel_c_mirror_b_%d", i);
217  gMC->SetBorderSurface(str_name, "DiscDIRC_focel_c", i, "DiscDIRC_focel_mirror_c", i, "DiscDIRC_focel_mirror_surf");
218  }
219 
220 
221  {
222  Double_t sellmeier_coeff[6] = {0.473115591, 0.631038719, 0.906404498,
223  0.012995717, 0.0041280992, 98.7685322 };
224 
225  const int n_entries = 800 - 200;
226  Double_t photon_momenta[n_entries];
227  Double_t rayleigh_length[n_entries];
228  Double_t lambda_um;
229  for(i=0, lambda_um = 0.8; i<n_entries; lambda_um-=0.001, i++)
230  {
231  photon_momenta[i] = 1.239841939*1E-6/lambda_um;
232  Double_t rindex = n_phase(lambda_um, sellmeier_coeff);
233  rayleigh_length[i] = 1E6/(45962092085.3903 * pow(rindex, 8.)/pow(lambda_um*1000.0, 4.)) *10.0;
234  }
235 
236  gMC->SetMaterialProperty( gMC->MediumId("FusedSilica"), "RAYLEIGH", n_entries, photon_momenta, rayleigh_length );
237  gMC->SetMaterialProperty( gMC->MediumId("FusedSilicaB"), "RAYLEIGH", n_entries, photon_momenta, rayleigh_length );
238  }
239 
240 }
Int_t i
Definition: run_full.C:25
Double_t wl_max_nm
Definition: PndDisc.h:99
Double_t
Double_t wl_min_nm
Definition: PndDisc.h:98
void PndDisc::EndOfEvent ( )
virtual

Interface implementation - EndofEvent.

Definition at line 551 of file PndDisc.cxx.

References clarr_particle_tracks, i, Print(), and Reset().

552 {
553  int n_entries = clarr_particle_tracks->GetEntries();
554  for(int i=0; i<n_entries; i++)
555  {
557  }
558  Reset();
559 }
MechFsc Print()
Int_t i
Definition: run_full.C:25
TClonesArray * clarr_particle_tracks
optical photon tracks
Definition: PndDisc.h:83
virtual void Reset()
Interface implementation - reset the 'containers'.
Definition: PndDisc.cxx:596
TClonesArray * PndDisc::GetCollection ( Int_t  iColl) const
virtual

Interface implementation - Get a data collection by index.

Parameters
iColl
Returns
pointer to TClonesArray or NULL if iColl out of range.

Definition at line 578 of file PndDisc.cxx.

References clarr_particle_tracks, and clarr_sensor_hits.

579 {
580  switch(iColl)
581  {
582  case 0: return clarr_sensor_hits;
583  case 1: return clarr_particle_tracks;
584  default: return NULL;
585  }
586 }
TClonesArray * clarr_sensor_hits
Definition: PndDisc.h:81
TClonesArray * clarr_particle_tracks
optical photon tracks
Definition: PndDisc.h:83
void PndDisc::Initialize ( )
virtual

Initialize will be called after the Geometry is created.

Definition at line 86 of file PndDisc.cxx.

References clarr_particle_tracks, clarr_sensor_hits, and Initialize().

87 {
88  // create collections (Register is called after Initialize)
89  clarr_sensor_hits = new TClonesArray("PndDiscSensorMCPoint");
90  clarr_particle_tracks = new TClonesArray("PndDiscParticleMCPoint");
91  clarr_particle_tracks->BypassStreamer(kTRUE);
92 
93  FairDetector::Initialize(); // call base class implementation
94 }
TClonesArray * clarr_sensor_hits
Definition: PndDisc.h:81
TClonesArray * clarr_particle_tracks
optical photon tracks
Definition: PndDisc.h:83
Mvd Initialize()
void PndDisc::PostTrack ( )
virtual

Definition at line 258 of file PndDisc.cxx.

259 {
260 #ifdef DISC_DIRC_VERBOSE
261  fLogger->GetLogger()->Info(MESSAGE_ORIGIN,
262  "tracknumber %d", gMC->GetStack()->GetCurrentTrackNumber() );
263 #endif
264 }
void PndDisc::PreTrack ( )
virtual

Definition at line 248 of file PndDisc.cxx.

249 {
250 #ifdef DISC_DIRC_VERBOSE
251  fLogger->GetLogger()->Info(MESSAGE_ORIGIN,
252  "tracknumber %d", gMC->GetStack()->GetCurrentTrackNumber() );
253 #endif
254 }
Bool_t PndDisc::ProcessHits ( FairVolume *  v = 0)
virtual

Definition at line 268 of file PndDisc.cxx.

References Bool_t, clarr_particle_tracks, Double_t, ev_header, gGeoManager, internal_reflection_angle_of_photons, PndDiscParticleMCPoint::is_primary, last_track_occurence, mom, PndDiscParticleMCPoint::mom_out, nextid_clarr_particle_tracks, nextid_clarr_sensor_hits, p1, p2, PndDiscSensorMCPoint::photon_entering_momentum, PndDiscSensorMCPoint::photon_entering_pos, photons_entering_optics, Pi, pos, PndDiscParticleMCPoint::pos_out, pt(), track, TString, PndDiscParticleMCPoint::volume_id, wl_max_nm, and wl_min_nm.

269 {
270  static TVector3 pos_in, mom_in;
271  static Double_t integrated_energy_deposit;
272 
273  const Int_t pid_optical_photon = 50000050;
274 
275  static const Int_t sensor_volume_id = gMC->VolId("DiscDIRC_sensor_pmt"); // find a good place for one time initialization
276  static const Int_t radiator_volume_id = gMC->VolId("DiscDIRC_Radiator"); // find a good place for one time initialization
277  static const Int_t prism_a_volume_id = gMC->VolId("DiscDIRC_prism_a");
278  static const Int_t prism_b_volume_id = gMC->VolId("DiscDIRC_prism_b");
279  static const Int_t prism_c_volume_id = gMC->VolId("DiscDIRC_prism_c");
280 
281  //static const Int_t sensor_volume_uid = // [R.K. 09/2018] unused variable
282  gGeoManager->GetUID("DiscDIRC_sensor_pmt"); // find a good place for one time initialization
283 
284  static Bool_t entered_inside = kFALSE;
285  static Int_t entered_track_id = -1;
286  static Int_t continued_track = -1;
287 
288  TParticle * track = gMC->GetStack()->GetCurrentTrack();
289 
290  // Some Logic Test:
291  static Int_t current_track;
292  if(gMC->IsTrackEntering())
293  {
294  current_track = gMC->GetStack()->GetCurrentTrackNumber();
295  }
296  else
297  {
298  if(current_track!=gMC->GetStack()->GetCurrentTrackNumber())
299  {
300  LOG(WARNING) << "LOGIC ERROR: ProcessHits called with track (%d) that has not entered the volume!" << gMC->GetStack()->GetCurrentTrackNumber();
301  }
302  }
303 
304  // get volume information:
305  Int_t volume_id, copy_no;
306  volume_id = gMC->CurrentVolID(copy_no);
307 
308  if(gMC->TrackPid() == pid_optical_photon)
309  {
310  //---------------------------------------------------
311  // handle optical photons
312  //---------------------------------------------------
313  TLorentzVector pos, mom;
314  gMC->TrackPosition(pos);
315  gMC->TrackMomentum(mom);
316 
317  Double_t mom_mag = mom.Vect().Mag();
318 
319  if(volume_id == radiator_volume_id)
320  {
321  std::map<int, double>::const_iterator it = internal_reflection_angle_of_photons.find(current_track);
323  {
324  double internal_reflection_angle = TMath::Pi()/2.0 - mom.Vect().Theta();
325  internal_reflection_angle_of_photons.insert(std::pair<int,double>(current_track, internal_reflection_angle));
326  }
327  }
328  else if(volume_id == prism_a_volume_id ||
329  volume_id == prism_b_volume_id ||
330  volume_id == prism_c_volume_id)
331  {
332  std::map<int, std::pair<TLorentzVector,TLorentzVector> >::iterator it = photons_entering_optics.find(current_track);
333  if(it == photons_entering_optics.end()) {
334  TLorentzVector pos_local, mom_local;
335  TString vol_path(gMC->CurrentVolPath());
336  gGeoManager->cd(vol_path.Data());
337  TGeoHMatrix * transform = gGeoManager->GetCurrentMatrix();
338  Double_t p1[3];
339  Double_t p2[3];
340  pos.Vect().GetXYZ(p1);
341  transform->MasterToLocal(p1, p2);
342  pos_local.SetXYZT(p2[0], p2[1], p2[2], 0.0);
343 
344  mom.Vect().GetXYZ(p1);
345  transform->MasterToLocalVect(p1, p2);
346  mom_local.SetXYZT(p2[0], p2[1], p2[2], 0.0);
347  photons_entering_optics.insert(std::pair< int, std::pair<TLorentzVector,TLorentzVector> >(current_track, std::pair<TLorentzVector,TLorentzVector>(pos_local, mom_local)));
348  }
349  }
350 
351  if( (mom_mag<1.239841939*1.e-06/wl_max_nm) || (mom_mag>1.239841939*1.e-06/wl_min_nm) )
352  gMC->StopTrack();
353 
354  //Photon is in focusing element
355  if(gMC->IsTrackEntering() && volume_id==sensor_volume_id)
356  {
357  TString vol_path(gMC->CurrentVolPath());
358  gGeoManager->cd(vol_path.Data());
359  TGeoHMatrix * transform = gGeoManager->GetCurrentMatrix();
360 
361  Double_t p1[3];
362  Double_t p2[3];
363  pos.Vect().GetXYZ(p1);
364  transform->MasterToLocal(p1, p2);
365  pos.SetXYZT(p2[0], p2[1], p2[2], 0.0);
366 
367  mom.Vect().GetXYZ(p1);
368  transform->MasterToLocalVect(p1, p2);
369  mom.SetXYZT(p2[0], p2[1], p2[2], 0.0);
370 
371  // Verbose output
372  if(fVerboseLevel > 0)
373  LOG(INFO) <<
374  "track->T(): "<<track->T()<<"\tgMC->TrackTime(): "<<gMC->TrackTime()<<"\tev_header->GetT(): "
375  <<ev_header->GetT()<<"\n";
376 
377  //double internal_reflection_angle = 0.0; //[R.K. 01/2017] unused variable
378  std::map<int, double>::iterator it = internal_reflection_angle_of_photons.find(current_track);
379 
381  {
382  LOG(WARNING) << "No registered total internal reflection angle for photon track id "<<current_track<<"\n";
383  } else {
384  //internal_reflection_angle = it->second; //[R.K. 01/2017] unused variable
386  }
387 
388  // store the information
389  Int_t volume_uid = gGeoManager->GetUID(gMC->VolName(volume_id));
390 
391  std::cout << "Sensor number: " << copy_no << ", volume id: " << volume_uid << ", photon position: x = " << pos[0] << " y = " << pos[1] << " z = " << pos[2] << std::endl;
392 
394  new((*clarr_sensor_hits)[nextid_clarr_sensor_hits++])
396  gMC->GetStack()->GetCurrentTrackNumber(),
397  copy_no,
398  volume_uid,
399  it->second,
400  pos.Vect(),
401  mom.Vect(),
402  gMC->TrackTime()*1E9,
403  gMC->TrackLength(),
404  gMC->Edep(),
405  track->T()
406  );
407 
408 
409  std::map<int, std::pair<TLorentzVector,TLorentzVector> >::iterator it_optics = photons_entering_optics.find(current_track);
410  if(it_optics == photons_entering_optics.end()) {
411  LOG(WARNING) << "No registered optics entrance data for photon track id "<<current_track<<"\n";
412  } else {
413  pt->photon_entering_pos = it_optics->second.first.Vect();
414  pt->photon_entering_momentum = it_optics->second.second.Vect();
415  }
416 
417  gMC->StopTrack();
418  }
419  }
420  else
421  {
422  //-------------------------------------------------------
423  // Handle non-optical tracks
424  //-------------------------------------------------------
425  TLorentzVector pos, mom;
426  gMC->TrackPosition(pos);
427  gMC->TrackMomentum(mom);
428 
429  if(gMC->IsTrackEntering() || (!gMC->IsTrackEntering() && entered_track_id<0))
430  {
431  if(entered_track_id >= 0)
432  LOG(FATAL) << "Track entering but processing of preceding track was not finished !!!";
433 
434  entered_track_id = gMC->GetStack()->GetCurrentTrackNumber();
435  if(entered_track_id<0)
436  LOG(FATAL) << "Track entering has neg id !!!";
437 
438  std::map<int, std::pair<int, double> >::iterator it = last_track_occurence.find(entered_track_id);
439 
440  continued_track = -1;
441  if(it != last_track_occurence.end()) {
443  if( volume_id == pt->volume_id || copy_no == pt->GetDetectorID() || gMC->TrackTime() == it->second.second ) {
444  continued_track = it->second.first;
445  integrated_energy_deposit = ((PndDiscParticleMCPoint*)clarr_particle_tracks->At(continued_track))->GetEnergyLoss();
446 
447  if(fVerboseLevel > 0) LOG(INFO) << Form("Continued: Track id %d, (vol %d, cpn %d, z=%f, T=%g, pdg %d)", entered_track_id, volume_id, copy_no, pos.Z(), gMC->TrackTime(), gMC->TrackPid());
448  }
449  }
450 
451  if(continued_track < 0)
452  {
453  pos_in = pos.Vect();
454  mom_in = mom.Vect();
455  integrated_energy_deposit = 0.0;
456  entered_inside = gMC->IsNewTrack();
457  if(fVerboseLevel > 0) LOG(INFO) << Form("Entered: Track id %d, (vol %d, cpn %d, z=%f, T=%g, pdg %d)", entered_track_id, volume_id, copy_no, pos.Z(), gMC->TrackTime() ,gMC->TrackPid());
458  }
459  }
460 
461  if(entered_track_id >= 0) {
462  if(entered_track_id != gMC->GetStack()->GetCurrentTrackNumber() )
463  LOG(FATAL) << "Track id changed !!! Need a stack !!!";
464  integrated_energy_deposit += gMC->Edep();
465  if(fVerboseLevel > 0) LOG(INFO) << "Add eloss for: Track id " << entered_track_id;
466  }
467 
468 
469 
470  if(gMC->IsTrackExiting() || gMC->IsTrackStop() || gMC->IsTrackDisappeared() || !gMC->IsTrackAlive())
471  {
472  // sanity check:
473  if(entered_track_id != gMC->GetStack()->GetCurrentTrackNumber() )
474  LOG(FATAL) <<"Track exiting without having entered (entered id = "<< entered_track_id<<")";
475 
476  if(fVerboseLevel > 0) LOG(INFO) << Form("Exiting: Track id %d (vol %d, cpn %d, Inside = %d, stop = %d, alive = %d, z=%f, T=%g)", entered_track_id, volume_id, copy_no, gMC->IsTrackInside(), gMC->IsTrackStop(), gMC->IsTrackAlive(), pos.Z(), gMC->TrackTime());
477 
478  if(gMC->IsTrackStop() || !gMC->IsTrackInside() || gMC->IsTrackDisappeared() || !gMC->IsTrackAlive())
479  {
480 
481  std::map<int, std::pair<int, double> >::iterator it = last_track_occurence.find(entered_track_id);
482  if(it == last_track_occurence.end())
483  last_track_occurence.insert(std::pair<int, std::pair<int, double> >(
484  entered_track_id, std::pair<int, double>(nextid_clarr_particle_tracks, gMC->TrackTime())
485  )
486  );
487  else
488  it->second.second = gMC->TrackTime();
489 
490  if(continued_track >= 0)
491  {
493  pt->pos_out = pos.Vect();
494  pt->mom_out = mom.Vect();
495  pt->SetTime(gMC->TrackTime()*1E9);
496  pt->SetLength(gMC->TrackLength());
497  pt->SetEnergyLoss(integrated_energy_deposit);
498  pt->is_primary = gMC->GetStack()->GetCurrentTrack()->IsPrimary();
499  }
500  else
501  {
502  Int_t volume_uid = gGeoManager->GetUID(gMC->VolName(volume_id));
503  new((*clarr_particle_tracks)[nextid_clarr_particle_tracks++])
505  entered_track_id,
506  copy_no,
507  volume_uid,
508  pos_in,
509  mom_in,
510  pos.Vect(),
511  mom.Vect(),
512  gMC->TrackTime()*1E9,
513  gMC->TrackLength(),
514  integrated_energy_deposit,
515  gMC->TrackCharge(),
516  gMC->TrackMass(),
517  gMC->TrackPid(),
518  entered_inside,
519  gMC->GetStack()->GetCurrentTrack()->IsPrimary()
520  );
521 
522 
523 
524  }
525  entered_track_id = -1;
526  }
527  else
528  {
529  LOG(INFO) << Form("Veto: Track id %d (Inside = %d)", entered_track_id, gMC->IsTrackInside());
530  }
531  }
532  }
533 
534  return kTRUE;
535 }
FairMCEventHeader * ev_header
Definition: PndDisc.h:91
TVector3 pos
std::map< int, double > internal_reflection_angle_of_photons
Definition: PndDisc.h:94
Double_t wl_max_nm
Definition: PndDisc.h:99
PndRiemannTrack track
Definition: RiemannTest.C:33
TClonesArray * clarr_particle_tracks
optical photon tracks
Definition: PndDisc.h:83
Double_t mom
Definition: plot_dirc.C:14
TGeoManager * gGeoManager
TString pt(TString pts, TString exts="px py pz")
Definition: invexp.C:133
std::map< int, std::pair< int, double > > last_track_occurence
Definition: PndDisc.h:92
Double_t
std::map< int, std::pair< TLorentzVector, TLorentzVector > > photons_entering_optics
Definition: PndDisc.h:95
int nextid_clarr_particle_tracks
Definition: PndDisc.h:89
Double_t wl_min_nm
Definition: PndDisc.h:98
int nextid_clarr_sensor_hits
tracks of other particles
Definition: PndDisc.h:87
TPad * p2
Definition: hist-t7.C:117
Int_t volume_id
FairMCPoint forces the implementation.
TPad * p1
Definition: hist-t7.C:116
Double_t Pi
void PndDisc::Register ( )
virtual

Interface implementation - Register the data collections in FairRootManager.

Definition at line 564 of file PndDisc.cxx.

References Bool_t, clarr_particle_tracks, and clarr_sensor_hits.

565 {
566  Bool_t serialized=kTRUE; // write arrays to file?
567  FairRootManager::Instance()->Register("DiscSensorMCPoint","DiscPoint", clarr_sensor_hits, serialized);
568  FairRootManager::Instance()->Register("DiscParticleMCPoint","DiscPoint", clarr_particle_tracks, serialized);
569 }
TClonesArray * clarr_sensor_hits
Definition: PndDisc.h:81
TClonesArray * clarr_particle_tracks
optical photon tracks
Definition: PndDisc.h:83
void PndDisc::Reset ( )
virtual

Interface implementation - reset the 'containers'.

This function is abstract in FairDetector, but is not called by the underlying framework as one would expect. The user has to call it explicitly.

Definition at line 596 of file PndDisc.cxx.

References clarr_particle_tracks, clarr_sensor_hits, nextid_clarr_particle_tracks, and nextid_clarr_sensor_hits.

Referenced by EndOfEvent().

597 {
598  if(clarr_sensor_hits != NULL) clarr_sensor_hits->Clear();
599  if(clarr_particle_tracks != NULL) clarr_particle_tracks->Clear();
602 }
TClonesArray * clarr_sensor_hits
Definition: PndDisc.h:81
TClonesArray * clarr_particle_tracks
optical photon tracks
Definition: PndDisc.h:83
int nextid_clarr_particle_tracks
Definition: PndDisc.h:89
int nextid_clarr_sensor_hits
tracks of other particles
Definition: PndDisc.h:87
void PndDisc::SetFilterInterval ( Double_t const &  wl_min_nm_,
Double_t const &  wl_max_nm_ 
)

Set the wavelength range of the bandpass filters.

Definition at line 634 of file PndDisc.cxx.

References wl_max_nm, and wl_min_nm.

Referenced by sim(), and sim_complete().

635 {
636  wl_min_nm = wl_min_nm_;
637  wl_max_nm = wl_max_nm_;
638 }
Double_t wl_max_nm
Definition: PndDisc.h:99
Double_t wl_min_nm
Definition: PndDisc.h:98
void PndDisc::StorePhotonTracks ( Bool_t  bval)

Enable/Disable the storage of photon track information.

The container for photon track information is created anyway, but it will not be filled if store_photon_tracks == false.

Definition at line 626 of file PndDisc.cxx.

References store_photon_tracks.

627 {
628  store_photon_tracks=bval;
629 }
Bool_t store_photon_tracks
Definition: PndDisc.h:77

Member Data Documentation

TClonesArray* PndDisc::clarr_particle_tracks
private

optical photon tracks

Definition at line 83 of file PndDisc.h.

Referenced by EndOfEvent(), GetCollection(), Initialize(), ProcessHits(), Register(), Reset(), and ~PndDisc().

TClonesArray* PndDisc::clarr_photon_tracks
private

hit on a photodetector surface

Definition at line 82 of file PndDisc.h.

TClonesArray* PndDisc::clarr_sensor_hits
private

Definition at line 81 of file PndDisc.h.

Referenced by GetCollection(), Initialize(), Register(), Reset(), and ~PndDisc().

int PndDisc::design_id
private

photon track ids that have been

Definition at line 97 of file PndDisc.h.

FairMCEventHeader* PndDisc::ev_header
private

Definition at line 91 of file PndDisc.h.

Referenced by BeginEvent(), and ProcessHits().

std::map<int, double> PndDisc::internal_reflection_angle_of_photons
private

Definition at line 94 of file PndDisc.h.

Referenced by BeginEvent(), and ProcessHits().

std::map<int, std::pair<int, double> > PndDisc::last_track_occurence
private

Definition at line 92 of file PndDisc.h.

Referenced by BeginEvent(), and ProcessHits().

std::set<std::string> PndDisc::names_of_sensitive_volumes
private

Definition at line 93 of file PndDisc.h.

Referenced by CheckIfSensitive(), and ConstructGeometry().

int PndDisc::nextid_clarr_particle_tracks
private

Definition at line 89 of file PndDisc.h.

Referenced by ProcessHits(), and Reset().

int PndDisc::nextid_clarr_photon_tracks
private

Definition at line 88 of file PndDisc.h.

int PndDisc::nextid_clarr_sensor_hits
private

tracks of other particles

Definition at line 87 of file PndDisc.h.

Referenced by ProcessHits(), and Reset().

TLorentzVector PndDisc::old_momentum
private

Definition at line 79 of file PndDisc.h.

std::map<int, std::pair<TLorentzVector,TLorentzVector> > PndDisc::photons_entering_optics
private

Definition at line 95 of file PndDisc.h.

Referenced by BeginEvent(), and ProcessHits().

Bool_t PndDisc::store_photon_tracks
private

Definition at line 77 of file PndDisc.h.

Referenced by StorePhotonTracks().

Bool_t PndDisc::track_is_photon
private

Definition at line 78 of file PndDisc.h.

Double_t PndDisc::wl_max_nm
private

Definition at line 99 of file PndDisc.h.

Referenced by ConstructOpGeometry(), ProcessHits(), and SetFilterInterval().

Double_t PndDisc::wl_min_nm
private

Definition at line 98 of file PndDisc.h.

Referenced by ConstructOpGeometry(), ProcessHits(), and SetFilterInterval().


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