All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Attributes | List of all members
sim::MCShowerRecoAlg Class Reference

#include <MCShowerRecoAlg.h>

Public Member Functions

 MCShowerRecoAlg (fhicl::ParameterSet const &pset)
 Default constructor with fhicl parameters. More...
 
std::unique_ptr< std::vector
< sim::MCShower > > 
Reconstruct (MCRecoPart &part_v, MCRecoEdep &edep_v)
 

Protected Attributes

MCShowerRecoPart fPartAlg
 
bool fDebugMode
 
double fMinShowerEnergy
 
unsigned int fMinNumDaughters
 

Detailed Description

Definition at line 23 of file MCShowerRecoAlg.h.

Constructor & Destructor Documentation

sim::MCShowerRecoAlg::MCShowerRecoAlg ( fhicl::ParameterSet const &  pset)
explicit

Default constructor with fhicl parameters.

Definition at line 35 of file MCShowerRecoAlg.cxx.

36  : fPartAlg(pset.get<fhicl::ParameterSet>("MCShowerRecoPart")),
37  fDebugMode(pset.get<bool>("DebugMode")),
38  fMinShowerEnergy(pset.get<double>("MinShowerEnergy")),
39  fMinNumDaughters(pset.get<unsigned int>("MinNumDaughters"))
40  //##################################################################
41  {
42  }
MCShowerRecoPart fPartAlg
unsigned int fMinNumDaughters

Member Function Documentation

std::unique_ptr< std::vector< sim::MCShower > > sim::MCShowerRecoAlg::Reconstruct ( MCRecoPart part_v,
MCRecoEdep edep_v 
)

Looping through the MCShower daughter's energy depositions

Looping through MCShower daughters

Looping through MCShowers

Definition at line 45 of file MCShowerRecoAlg.cxx.

47  {
48 
49  art::ServiceHandle<geo::Geometry const> geo;
50 
51  auto pindex = details::createPlaneIndexMap();
52 
53  fPartAlg.ConstructShower(part_v);
54  auto result = std::make_unique<std::vector<sim::MCShower>>();
55  auto& mcshower = *result;
56  //std::vector<sim::MCShower> mcshower;
57  // Get shower info from grouped particles
58  const std::vector<unsigned int> shower_index_v = fPartAlg.ShowerMothers();
59  mcshower.reserve(shower_index_v.size());
60  std::vector<size_t> mcs_to_spart_v;
61  mcs_to_spart_v.reserve(shower_index_v.size());
62 
63  bool daughter_stored=false;
64  for(size_t shower_index = 0; shower_index < shower_index_v.size(); ++shower_index) {
65 
66  unsigned int shower_candidate = shower_index_v.at(shower_index);
67  auto const& shower_part = part_v.at(shower_candidate);
68 
69  unsigned int mother_track = part_v.MotherTrackID(shower_candidate);
70  unsigned int ancestor_track = part_v.AncestorTrackID(shower_candidate);
71 
72  if(mother_track == kINVALID_UINT || ancestor_track == kINVALID_UINT)
73 
74  throw cet::exception(__FUNCTION__) << "LOGIC ERROR: mother/ancestor track ID is invalid!";
75 
76  MCMiniPart mother_part;
77  MCMiniPart ancestor_part;
78 
79  unsigned int mother_index = part_v.TrackToParticleIndex(mother_track);
80  unsigned int ancestor_index = part_v.TrackToParticleIndex(ancestor_track);
81 
82  if(mother_index != kINVALID_UINT) mother_part = part_v[mother_index];
83  else mother_part._track_id = mother_track;
84 
85  if(ancestor_index != kINVALID_UINT) ancestor_part = part_v[ancestor_index];
86  else ancestor_part._track_id = ancestor_track;
87 
88  double shower_g4_energy = shower_part._start_mom[3];
89 
90  if(fDebugMode)
91 
92  std::cout << "Found MCShower with mother energy: " << shower_g4_energy << " MeV";
93 
94  // Skip if mother energy is less than the enery threshold
95  if(shower_g4_energy < fMinShowerEnergy) {
96  if(fDebugMode)
97  std::cout << " ... below energy threshold: skipping!"<<std::endl;
98 
99  continue;
100  }else if(shower_part._daughters.size() < fMinNumDaughters) {
101  if(fDebugMode)
102  std::cout << " ... below # daughter particle count threshold: skipping!"<<std::endl;
103 
104  continue;
105  }else if(fDebugMode) {
106  std::cout << " ... condition matched. Storing this MCShower..."<<std::endl;
107  }
108 
109  // Record this MCShower
110  mcs_to_spart_v.push_back(shower_index);
111 
112  if(fDebugMode)
113 
114  std::cout << " Storage index " << mcshower.size() << " => Shower index " << shower_index
115  << std::endl;
116 
117  ::sim::MCShower shower_prof;
118 
119  shower_prof.Origin ( shower_part._origin );
120  shower_prof.PdgCode ( shower_part._pdgcode );
121  shower_prof.TrackID ( shower_part._track_id );
122  shower_prof.Process ( shower_part._process );
123 
124  shower_prof.MotherPdgCode ( mother_part._pdgcode );
125  shower_prof.MotherTrackID ( mother_part._track_id );
126  shower_prof.MotherProcess ( mother_part._process );
127 
128  shower_prof.AncestorPdgCode ( ancestor_part._pdgcode );
129  shower_prof.AncestorTrackID ( ancestor_part._track_id );
130  shower_prof.AncestorProcess ( ancestor_part._process );
131 
132  shower_prof.Start ( MCStep ( shower_part._start_vtx, shower_part._start_mom ) );
133  shower_prof.End ( MCStep ( shower_part._end_vtx, shower_part._end_mom ) );
134  shower_prof.MotherStart ( MCStep ( mother_part._start_vtx, mother_part._start_mom ) );
135  shower_prof.MotherEnd ( MCStep ( mother_part._end_vtx, mother_part._end_mom ) );
136  shower_prof.AncestorStart ( MCStep ( ancestor_part._start_vtx, ancestor_part._start_mom ) );
137  shower_prof.AncestorEnd ( MCStep ( ancestor_part._end_vtx, ancestor_part._end_mom ) );
138 
139  // Daughter list
140  std::vector<unsigned int> daughter_track_id;
141  daughter_track_id.reserve( fPartAlg.ShowerDaughters(shower_index).size() );
142 
143  for(auto const& index : fPartAlg.ShowerDaughters(shower_index))
144 
145  daughter_track_id.push_back( part_v.at(index)._track_id );
146 
147  shower_prof.DaughterTrackID(daughter_track_id);
148 
149  if(!daughter_stored && daughter_track_id.size()>1) daughter_stored=true;
150 
151  mcshower.push_back(shower_prof);
152  }
153 
154  if(fDebugMode)
155  std::cout << " Found " << mcshower.size() << " MCShowers. Now computing DetProfile position..." << std::endl;
156 
157 
158  //
159  // Daughter vtx
160  //
161  std::vector<TLorentzVector> mcs_daughter_vtx_v(mcshower.size(),TLorentzVector(sim::kINVALID_DOUBLE,
165  std::vector<TLorentzVector> mcs_daughter_mom_v ( mcshower.size(), TLorentzVector() );
166 
167  std::vector< std::vector<double> > plane_charge_v ( mcshower.size(), std::vector<double>(3,0) );
168  std::vector< std::vector<double> > plane_dqdx_v ( mcshower.size(), std::vector<double>(3,0) );
169 
170  //For dEdx Calculation
171  std::vector<double> mcs_daughter_dedx_v ( mcshower.size(), 0 );
172  std::vector<double> mcs_daughter_dedxRAD_v ( mcshower.size(), 0 );
173  std::vector<TVector3> mcs_daughter_dir_v ( mcshower.size(), TVector3() );
174 
175  for(size_t mcs_index=0; mcs_index<mcshower.size(); ++mcs_index) {
176 
177  auto& mcs_daughter_vtx = mcs_daughter_vtx_v[mcs_index];
178  auto& mcs_daughter_mom = mcs_daughter_mom_v[mcs_index];
179  auto& mcs_daughter_dedx = mcs_daughter_dedx_v[mcs_index];
180  auto& mcs_daughter_dedxRAD = mcs_daughter_dedxRAD_v[mcs_index];
181  auto& mcs_daughter_dir = mcs_daughter_dir_v[mcs_index];
182  auto& plane_charge = plane_charge_v[mcs_index];
183  auto& plane_dqdx = plane_dqdx_v[mcs_index];
184 
185  for(auto const& daughter_trk_id : mcshower[mcs_index].DaughterTrackID()) {
186 
187  auto const daughter_part_index = part_v.TrackToParticleIndex(daughter_trk_id);
188 
189  auto const& daughter_part = part_v[daughter_part_index];
190 
191  auto const daughter_edep_index = edep_v.TrackToEdepIndex(daughter_trk_id);
192 
193  if(daughter_edep_index<0) continue;
194 
195  auto const& daughter_edep = edep_v.GetEdepArrayAt(daughter_edep_index);
196 
197  if(!(daughter_edep.size())) continue;
198 
199  // Record first daughter's vtx point
200  double min_dist = sim::kINVALID_DOUBLE;
201  for(auto const& edep : daughter_edep) {
202 
203  double dist = sqrt( pow(edep.pos._x - daughter_part._start_vtx[0],2) +
204  pow(edep.pos._y - daughter_part._start_vtx[1],2) +
205  pow(edep.pos._z - daughter_part._start_vtx[2],2) );
206 
207  if(dist < min_dist) {
208  min_dist = dist;
209  mcs_daughter_vtx[0] = edep.pos._x;
210  mcs_daughter_vtx[1] = edep.pos._y;
211  mcs_daughter_vtx[2] = edep.pos._z;
212  mcs_daughter_vtx[3] = (dist/100. / 2.998e8)*1.e9 + daughter_part._start_vtx[3];
213  }
214 
215  }
216  if(!daughter_stored) {
217  // If daughter is not stored, and shower id energetic enough, attempt to include angle info
218  std::vector<double> shower_dir(3,0);
219  shower_dir[0] = mcshower[mcs_index].Start().Px();
220  shower_dir[1] = mcshower[mcs_index].Start().Py();
221  shower_dir[2] = mcshower[mcs_index].Start().Pz();
222  double magnitude = 0;
223  for(size_t i=0; i<3; ++i)
224  magnitude += pow(shower_dir[i],2);
225 
226  magnitude = sqrt(magnitude);
227 
228  if(magnitude > 1.e-10) {
229  // If enough momentum, include angle info
230  min_dist = sim::kINVALID_DOUBLE;
231 
232  for(auto& v : shower_dir) v /= magnitude;
233 
234  for(auto const& edep : daughter_edep) {
235  std::vector<double> shower_dep_dir(3,0);
236  shower_dep_dir[0] = edep.pos._x - mcshower[mcs_index].Start().X();
237  shower_dep_dir[1] = edep.pos._y - mcshower[mcs_index].Start().Y();
238  shower_dep_dir[2] = edep.pos._z - mcshower[mcs_index].Start().Z();
239 
240  double dist = sqrt( pow(shower_dep_dir[0],2) + pow(shower_dep_dir[1],2) + pow(shower_dep_dir[2],2) );
241  for(auto& v : shower_dep_dir) v /= dist;
242 
243  double angle = acos( shower_dep_dir[0] * shower_dir[0] +
244  shower_dep_dir[1] * shower_dir[1] +
245  shower_dep_dir[2] * shower_dir[2] ) / TMath::Pi() * 180.;
246 
247  if(dist < min_dist && angle < 10) {
248 
249  min_dist = dist;
250  mcs_daughter_vtx[0] = edep.pos._x;
251  mcs_daughter_vtx[1] = edep.pos._y;
252  mcs_daughter_vtx[2] = edep.pos._z;
253  mcs_daughter_vtx[3] = (dist/100. / 2.998e8)*1.e9 + mcshower[mcs_index].Start().T();
254  }
255  }
256  }
257  }
258  break;
259  }
260  // Now take care of momentum & plane charge
261 
262  std::vector<double> mom(3,0);
263  for(auto const& daughter_trk_id : mcshower[mcs_index].DaughterTrackID()) {
264 
265  //auto const daughter_part_index = part_v.TrackToParticleIndex(daughter_trk_id);
266 
267  // for c2: daughter_part is unused
268  //auto const& daughter_part = part_v[daughter_part_index];
269 
270  auto const daughter_edep_index = edep_v.TrackToEdepIndex(daughter_trk_id);
271 
272  if(daughter_edep_index<0) continue;
273 
274  auto const& daughter_edep = edep_v.GetEdepArrayAt(daughter_edep_index);
275 
276  if(!(daughter_edep.size())) continue;
277 
278  //bool first=true; // unused
279  for(auto const& edep : daughter_edep) {
280 
281  // Compute unit vector to this energy deposition
282  mom[0] = edep.pos._x - mcs_daughter_vtx[0];
283  mom[1] = edep.pos._y - mcs_daughter_vtx[1];
284  mom[2] = edep.pos._z - mcs_daughter_vtx[2];
285 
286  // Weight by energy (momentum)
287  double magnitude = sqrt(pow(mom.at(0),2) + pow(mom.at(1),2) + pow(mom.at(2),2));
288 
289  double energy = 0;
290  double npid = 0;
291  for(auto const& pid_energy : edep.deps) {
292  npid++;
293  energy += pid_energy.energy;
294 
295  }
296  energy /= npid;
297  if(magnitude>1.e-10) {
298  mom.at(0) = mom.at(0) * energy / magnitude;
299  mom.at(1) = mom.at(1) * energy / magnitude;
300  mom.at(2) = mom.at(2) * energy / magnitude;
301  mcs_daughter_mom[0] += mom.at(0);
302  mcs_daughter_mom[1] += mom.at(1);
303  mcs_daughter_mom[2] += mom.at(2);
304  }
305  //Determine the direction of the shower right at the start point
306  double E = 0;
307  double N = 0;
308  if(sqrt( pow( edep.pos._x - mcs_daughter_vtx[0],2) +
309  pow( edep.pos._y - mcs_daughter_vtx[1],2) +
310  pow( edep.pos._z - mcs_daughter_vtx[2],2)) < 2.4 && magnitude>1.e-10){
311 
312  mcs_daughter_dir[0] += mom.at(0);
313  mcs_daughter_dir[1] += mom.at(1);
314  mcs_daughter_dir[2] += mom.at(2);
315  E += energy;
316  N += 1;
317  }
318 
319  if(E > 0) E /= N;
320  mcs_daughter_dedxRAD += E;
321 
322  mcs_daughter_mom[3] += energy;
323 
324  // Charge
325  auto const pid = edep.pid;
326  auto q_i = pindex.find(pid);
327  if(q_i != pindex.end())
328  plane_charge[pid.Plane] += (double)(edep.deps[pindex[pid]].charge);
329 
330  }///Looping through the MCShower daughter's energy depositions
331 
332  }///Looping through MCShower daughters
333  mcs_daughter_dedxRAD /= 2.4;
334 
335  for(auto const& daughter_trk_id : mcshower[mcs_index].DaughterTrackID()) {
336 
337  //auto const daughter_part_index = part_v.TrackToParticleIndex(daughter_trk_id);
338 
339  // for c2: daughter_part is unused
340  //auto const& daughter_part = part_v[daughter_part_index];
341 
342  auto const daughter_edep_index = edep_v.TrackToEdepIndex(daughter_trk_id);
343 
344  if(daughter_edep_index<0) continue;
345 
346  auto const& daughter_edep = edep_v.GetEdepArrayAt(daughter_edep_index);
347 
348  if(!(daughter_edep.size())) continue;
349 
350  for(auto const& edep : daughter_edep) {
351 
352  //Defining dEdx
353  //Need to define a plane through the shower start point (x_0, y_0, z_0) with a normal along the momentum vector of the shower
354  //The plane will be defined in the typical way:
355  // a*x + b*y + c*z + d = 0
356  // where, a = dir_x, b = dir_y, c = dir_z, d = - (a*x_0+b*y_0+c*z_0)
357  // then the *signed* distance of any point (x_1, y_1, z_1) from this plane is:
358  // D = (a*x_1 + b*y_1 + c*z_1 + d )/sqrt( pow(a,2) + pow(b,2) + pow(c,2))
359 
360 
361 
362  double p_mag = sqrt( pow(mcs_daughter_dir[0],2) + pow(mcs_daughter_dir[1],2) + pow(mcs_daughter_dir[2],2) );
363  double a = 0, b = 0, c = 0, d = 0;
364  if(p_mag > 1.e-10){
365  a = mcs_daughter_dir[0]/p_mag;
366  b = mcs_daughter_dir[1]/p_mag;
367  c = mcs_daughter_dir[2]/p_mag;
368  d = -1*(a*mcs_daughter_vtx[0] + b*mcs_daughter_vtx[1] + c*mcs_daughter_vtx[2]);
369  }
370  else{mcs_daughter_dedx += 0; continue;}
371  //Radial Distance
372  if( (a*edep.pos._x + b*edep.pos._y + c*edep.pos._z + d)/sqrt( pow(a,2) + pow(b,2) + pow(c,2)) < 2.4 &&
373  (a*edep.pos._x + b*edep.pos._y + c*edep.pos._z + d)/sqrt( pow(a,2) + pow(b,2) + pow(c,2)) > 0){
374 
375  double E = 0;
376  double N = 0;
377 
378  for(auto const& pid_energy : edep.deps) {
379  N += 1;
380  E += pid_energy.energy;
381  }
382 
383  if(N > 0){
384  E /= N;
385  }
386  else{ E = 0;}
387 
388  mcs_daughter_dedx += E;
389 
390  // Charge
391  auto const pid = edep.pid;
392  auto q_i = pindex.find(pid);
393  if(q_i != pindex.end())
394  plane_dqdx[pid.Plane] += (double)(edep.deps[pindex[pid]].charge);
395  }
396  }
397  }
398  mcs_daughter_dedx /= 2.4;
399  plane_dqdx.at(0) /= 2.4;
400  plane_dqdx.at(1) /= 2.4;
401  plane_dqdx.at(2) /= 2.4;
402 
403 
404  }///Looping through MCShowers
405 
406  if(fDebugMode)
407  std::cout << " Found " << mcshower.size() << " MCShowers. Now storing..." << std::endl;
408 
409  // Store plane charge & daughter momentum
410  for(size_t mcs_index=0; mcs_index<mcshower.size(); ++mcs_index) {
411 
412  auto& daughter_vtx = mcs_daughter_vtx_v[mcs_index];
413  auto& daughter_mom = mcs_daughter_mom_v[mcs_index];
414  auto& daughter_dedx = mcs_daughter_dedx_v[mcs_index];
415  auto& daughter_dedxRAD = mcs_daughter_dedxRAD_v[mcs_index];
416  auto& daughter_dir = mcs_daughter_dir_v[mcs_index];
417  auto& plane_charge = plane_charge_v[mcs_index];
418  auto& plane_dqdx = plane_dqdx_v[mcs_index];
419 
420  double magnitude = sqrt(pow(daughter_mom[0],2)+pow(daughter_mom[1],2)+pow(daughter_mom[2],2));
421 
422  if(daughter_mom[3]>1.e-10) {
423  daughter_mom[0] *= daughter_mom[3]/magnitude;
424  daughter_mom[1] *= daughter_mom[3]/magnitude;
425  daughter_mom[2] *= daughter_mom[3]/magnitude;
426  }else
427  for(size_t i=0; i<4; ++i) daughter_mom[i]=0;
428 
429  mcshower.at(mcs_index).DetProfile( MCStep( daughter_vtx, daughter_mom ) );
430  mcshower.at(mcs_index).Charge(plane_charge);
431  mcshower.at(mcs_index).dQdx(plane_dqdx);
432  mcshower.at(mcs_index).dEdx(daughter_dedx);
433  mcshower.at(mcs_index).dEdxRAD(daughter_dedxRAD);
434  mcshower.at(mcs_index).StartDir(daughter_dir);
435 
436  }
437 
438  if(fDebugMode) {
439 
440  for(auto const& prof : mcshower) {
441 
442  std::cout
443 
444  << Form(" Shower particle: PDG=%d : Track ID=%d Start @ (%g,%g,%g,%g) with Momentum (%g,%g,%g,%g)",
445  prof.PdgCode(), prof.TrackID(),
446  prof.Start().X(),prof.Start().Y(),prof.Start().Z(),prof.Start().T(),
447  prof.Start().Px(),prof.Start().Py(),prof.Start().Pz(),prof.Start().E())
448  << std::endl
449  << Form(" Mother particle: PDG=%d : Track ID=%d Start @ (%g,%g,%g,%g) with Momentum (%g,%g,%g,%g)",
450  prof.MotherPdgCode(), prof.MotherTrackID(),
451  prof.MotherStart().X(),prof.MotherStart().Y(),prof.MotherStart().Z(),prof.MotherStart().T(),
452  prof.MotherStart().Px(),prof.MotherStart().Py(),prof.MotherStart().Pz(),prof.MotherStart().E())
453  << std::endl
454  << Form(" Ancestor particle: PDG=%d : Track ID=%d Start @ (%g,%g,%g,%g) with Momentum (%g,%g,%g,%g)",
455  prof.AncestorPdgCode(), prof.AncestorTrackID(),
456  prof.AncestorStart().X(),prof.AncestorStart().Y(),prof.AncestorStart().Z(),prof.AncestorStart().T(),
457  prof.AncestorStart().Px(),prof.AncestorStart().Py(),prof.AncestorStart().Pz(),prof.AncestorStart().E())
458  << std::endl
459  << Form(" ... with %zu daughters: Start @ (%g,%g,%g,%g) with Momentum (%g,%g,%g,%g)",
460  prof.DaughterTrackID().size(),
461  prof.DetProfile().X(),prof.DetProfile().Y(),prof.DetProfile().Z(),prof.DetProfile().T(),
462  prof.DetProfile().Px(),prof.DetProfile().Py(),prof.DetProfile().Pz(),prof.DetProfile().E())
463  << std::endl
464  << " Charge per plane: ";
465  size_t const nPlanes = prof.Charge().size();
466  for(size_t i=0; i<nPlanes; ++i) {
467 
468  std::cout << " | Plane " << i << std::flush;
469  std::cout << " ... Q = " << prof.Charge(i) << std::flush;
470 
471  }
472  std::cout<<std::endl<<std::endl;
473  }
474  }
475  return result;
476  }
const double kINVALID_DOUBLE
Definition: MCLimits.h:10
const MCStep & End() const
Definition: MCShower.h:56
unsigned int TrackID() const
Definition: MCShower.h:53
const std::vector< unsigned int > ShowerMothers() const
MCShowerRecoPart fPartAlg
std::map< geo::PlaneID, size_t > createPlaneIndexMap()
Definition: MCRecoEdep.cxx:26
int PdgCode() const
Definition: MCShower.h:52
process_name E
unsigned int fMinNumDaughters
process_name gaushit a
const std::vector< unsigned int > & DaughterTrackID() const
Definition: MCShower.h:72
simb::Origin_t Origin() const
Definition: MCShower.h:50
int MotherPdgCode() const
Definition: MCShower.h:58
const std::string & AncestorProcess() const
Definition: MCShower.h:66
const MCStep & AncestorStart() const
Definition: MCShower.h:67
const std::string & MotherProcess() const
Definition: MCShower.h:60
const std::vector< unsigned int > & ShowerDaughters(const unsigned int shower_id) const
unsigned int AncestorTrackID() const
Definition: MCShower.h:65
const MCStep & AncestorEnd() const
Definition: MCShower.h:68
const MCStep & Start() const
Definition: MCShower.h:55
const MCStep & MotherEnd() const
Definition: MCShower.h:62
constexpr double dist(const TReal *x, const TReal *y, const unsigned int dimension)
process_name largeant stream1 can override from command line with o or output physics producers generator N
do i e
const unsigned int kINVALID_UINT
Definition: MCLimits.h:14
unsigned int MotherTrackID() const
Definition: MCShower.h:59
const std::string & Process() const
Definition: MCShower.h:54
const MCStep & MotherStart() const
Definition: MCShower.h:61
finds tracks best matching by angle
int AncestorPdgCode() const
Definition: MCShower.h:64
BEGIN_PROLOG could also be cout
void ConstructShower(const MCRecoPart &part_v)
Main function to read-in data and fill variables in this algorithm to reconstruct MC shower...

Member Data Documentation

bool sim::MCShowerRecoAlg::fDebugMode
protected

Definition at line 35 of file MCShowerRecoAlg.h.

unsigned int sim::MCShowerRecoAlg::fMinNumDaughters
protected

Definition at line 37 of file MCShowerRecoAlg.h.

double sim::MCShowerRecoAlg::fMinShowerEnergy
protected

Definition at line 36 of file MCShowerRecoAlg.h.

MCShowerRecoPart sim::MCShowerRecoAlg::fPartAlg
protected

Definition at line 34 of file MCShowerRecoAlg.h.


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