Class Reference for E1039 Core & Analysis Software
PHG4DSTReader.cc
Go to the documentation of this file.
1 // $Id: PHG4DSTReader.cc,v 1.11 2015/01/06 02:52:07 jinhuang Exp $
2 
11 #include "PHG4DSTReader.h"
12 
15 
16 #include <g4main/PHG4InEvent.h>
17 //#include <g4main/PHG4Particle.h>
18 #include <g4main/PHG4HitEval.h>
19 #include <g4main/PHG4Particlev2.h>
20 #include <g4main/PHG4VtxPointv1.h>
21 
22 #include <fun4all/PHTFileServer.h>
24 
25 #include <phool/getClass.h>
26 
27 #include <TClonesArray.h>
28 #include <TTree.h>
29 
30 #include <boost/foreach.hpp>
31 #include <map>
32 #include <set>
33 #include <cassert>
34 
35 #include<sstream>
36 
37 #include <memory>
38 
39 
40 using namespace std;
41 
45 
46 PHG4DSTReader::PHG4DSTReader(const string &filename) :
47  SubsysReco("PHG4DSTReader"), nblocks(0), _event(0), //
48  _out_file_name(filename), /*_file(nullptr), */_T(nullptr), //
49  _save_particle(true), _load_all_particle(false), _load_active_particle(
50  true), _save_vertex(true)
51 {
52  // TODO Auto-generated constructor stub
53 
54 }
55 
57 {
58  cout << "PHG4DSTReader::destructor - Clean ups" << endl;
59 
60  if (_T)
61  {
62  _T->ResetBranchAddresses();
63  }
64 
65  _records.clear();
66 }
67 
68 int
70 {
71 
72  const int arr_size = 100;
73 
74 // BOOST_FOREACH(string nodenam, _node_postfix)
75  for (vector<string>::const_iterator it = _node_postfix.begin();
76  it != _node_postfix.end(); ++it)
77  {
78  const char * class_name = hit_type::Class()->GetName();
79 
80  const string & nodenam = *it;
81 
82  string hname = "G4HIT_" + nodenam;
83 // _node_name.push_back(hname);
84  cout << "PHG4DSTReader::Init - saving hits from node: " << hname << " - "
85  << class_name << endl;
86 
87  record rec;
88  rec._cnt = 0;
89  rec._name = hname;
90  rec._arr = std::make_shared<TClonesArray>(class_name, arr_size);
91  rec._arr_ptr = rec._arr.get();
92  rec._type = record::typ_hit;
93 
94  _records.push_back(rec);
95 
96  nblocks++;
97  }
98 
99  if (_save_particle)
100  {
101  //save particles
102 
103  const char * class_name = part_type::Class()->GetName();
104 
105  cout << "PHG4DSTReader::Init - saving Particles node:" << class_name
106  << endl;
107 
108  record rec;
109  rec._cnt = 0;
110  rec._name = "PHG4Particle";
111  rec._arr = std::make_shared<TClonesArray>(class_name, arr_size);
112  rec._arr_ptr = rec._arr.get();
113  rec._type = record::typ_part;
114 
115  _records.push_back(rec);
116 
117  nblocks++;
118  }
119 
120  if (_save_vertex)
121  {
122  //save particles
123 
124  const char * class_name = vertex_type::Class()->GetName();
125 
126  cout << "PHG4DSTReader::Init - saving vertex for particles under node:"
127  << class_name << endl;
128 
129  record rec;
130  rec._cnt = 0;
131  rec._name = "PHG4VtxPoint";
132  rec._arr = std::make_shared<TClonesArray>(class_name, arr_size);
133  rec._arr_ptr = rec._arr.get();
135 
136  _records.push_back(rec);
137 
138  nblocks++;
139  }
140 
141  cout << "PHG4DSTReader::Init - requested " << nblocks << " nodes" << endl;
142 
143  build_tree();
144 
145  return 0;
146 }
147 
148 void
150 {
151  cout << "PHG4DSTReader::build_tree - output to " << _out_file_name << endl;
152 
153  static const int BUFFER_SIZE = 32000;
154 
155  // open TFile
156  PHTFileServer::get().open(_out_file_name, "RECREATE");
157 
158  _T = new TTree("T", "PHG4DSTReader");
159 
160  nblocks = 0;
161  for (records_t::iterator it = _records.begin(); it != _records.end(); ++it)
162  {
163  record & rec = *it;
164 
165  cout << "PHG4DSTReader::build_tree - Add " << rec._name << endl;
166 
167  const string name_cnt = "n_" + rec._name;
168  const string name_cnt_desc = name_cnt + "/I";
169  _T->Branch(name_cnt.c_str(), &(rec._cnt), name_cnt_desc.c_str(),
170  BUFFER_SIZE);
171  _T->Branch(rec._name.c_str(), &(rec._arr_ptr), BUFFER_SIZE, 99);
172 
173  nblocks++;
174  }
175 
176  cout << "PHG4DSTReader::build_tree - added " << nblocks << " nodes" << endl;
177 
178  _T->SetAutoSave(16000);
179 }
180 
181 int
183 {
184 
185 // const double significand = _event / TMath::Power(10, (int) (log10(_event)));
186 //
187 // if (fmod(significand, 1.0) == 0 && significand <= 10)
188 // cout << "PHG4DSTReader::process_event - " << _event << endl;
189  _event++;
190 
191  //clean ups
192  _particle_set.clear();
193  _vertex_set.clear();
194 
196  PHG4TruthInfoContainer>(topNode, "G4TruthInfo");
197  if (!truthInfoList)
198  {
199  if (_event < 2)
200  cout
201  << "PHG4DSTReader::process_event - Error - can not find node G4TruthInfo. Quit processing!"
202  << endl;
204  }
205 
206  for (records_t::iterator it = _records.begin(); it != _records.end(); ++it)
207  {
208  record & rec = *it;
209 
210  rec._cnt = 0;
211  assert(rec._arr.get() == rec._arr_ptr);
212  assert(rec._arr.get());
213  rec._arr->Clear();
214 
215  if (rec._type == record::typ_hit)
216  {
217  if (Verbosity() >= 2)
218  cout << "PHG4DSTReader::process_event - processing " << rec._name
219  << endl;
220 
221  PHG4HitContainer *hits = findNode::getClass<PHG4HitContainer>(topNode,
222  rec._name);
223  if (!hits)
224  {
225  if (_event < 2)
226  cout
227  << "PHG4DSTReader::process_event - Error - can not find node "
228  << rec._name << endl;
229 
230  }
231  else
232  {
233  PHG4HitContainer::ConstRange hit_range = hits->getHits();
234 
235  if (Verbosity() >= 2)
236  cout << "PHG4DSTReader::process_event - processing "
237  << rec._name << " and received " << hits->size() << " hits"
238  << endl;
239 
240  for (PHG4HitContainer::ConstIterator hit_iter = hit_range.first;
241  hit_iter != hit_range.second; hit_iter++)
242  {
243  PHG4Hit * hit = hit_iter->second;
244 
245 // hit_type * hit = dynamic_cast<hit_type *>(hit_raw);
246 //
247  assert(hit);
248 
249  new ((*(rec._arr.get()))[rec._cnt]) hit_type(*hit);
250 
251  hit_type * new_hit =
252  dynamic_cast<hit_type *>(rec._arr.get()->At(rec._cnt));
253  assert(new_hit);
254 
255 // for (int i = 0; i < 2; i++)
256 // {
257 // new_hit->set_x(i, hit->get_x(i));
258 // new_hit->set_y(i, hit->get_y(i));
259 // new_hit->set_z(i, hit->get_z(i));
260 // new_hit->set_t(i, hit->get_t(i));
261 // }
262 // new_hit->set_edep(hit->get_edep());
263 // new_hit->set_layer(hit->get_layer());
264 // new_hit->set_hit_id(hit->get_hit_id());
265 // new_hit->set_trkid(hit->get_trkid());
266 
267 // *new_hit = (*hit);
268 
270  _particle_set.insert(hit->get_trkid());
271 
272  if (Verbosity() >= 2)
273  cout << "PHG4DSTReader::process_event - processing "
274  << rec._name << " and hit " << hit->get_hit_id()
275  << " with track id " << hit->get_trkid() << endl;
276 
277  rec._cnt++;
278  }
279  } // if (!hits)
280  }
281  else if (rec._type == record::typ_part)
282  {
283 
284  map<int, PHG4Particle *>::const_iterator particle_iter;
285 
286  if (_load_all_particle)
287  {
288  static bool once = true;
289  if (once)
290  {
291  cout
292  << "PHG4DSTReader::process_event - will load all particle from G4TruthInfo"
293  << endl;
294 
295  once = false;
296  }
297 
298  for (particle_iter = truthInfoList->GetMap().begin();
299  particle_iter != truthInfoList->GetMap().end();
300  particle_iter++)
301  {
302 
303  _particle_set.insert(particle_iter->first);
304 
305  }
306 
307  } // if (_load_all_particle)
308  else
309  {
310  static bool once = true;
311  if (once)
312  {
313  cout
314  << "PHG4DSTReader::process_event - will load primary particle from G4TruthInfo"
315  << endl;
316 
317  once = false;
318  }
319 
320  PHG4TruthInfoContainer::ConstRange primary_range =
321  truthInfoList->GetPrimaryParticleRange();
322 
323  for (PHG4TruthInfoContainer::ConstIterator particle_iter =
324  primary_range.first; particle_iter != primary_range.second;
325  ++particle_iter)
326  {
327 
328  _particle_set.insert(particle_iter->first);
329 
330  } // if (_load_all_particle) else
331  }
332  for (PartSet_t::const_iterator i = _particle_set.begin();
333  i != _particle_set.end(); i++)
334  {
335  particle_iter = truthInfoList->GetMap().find(*i);
336  if (particle_iter == truthInfoList->GetMap().end())
337  {
338 
339  cout
340  << "PHG4DSTReader::process_event - ERROR - can not find particle/track ID "
341  << *i << " in G4TruthInfo" << endl;
342 
343  continue;
344  }
345 
346  PHG4Particle * part = particle_iter->second;
347  add_particle(rec, part);
348  } // for(PartSet_t::const_iterator i = _particle_set.begin();i!=_particle_set.end();i++)
349 
350  } // if (rec._type == record::typ_part)
351  else if (rec._type == record::typ_vertex)
352  {
353  static bool once = true;
354  if (once)
355  {
356  cout
357  << "PHG4DSTReader::process_event - will load vertex from G4TruthInfo"
358  << endl;
359 
360  once = false;
361  }
362 
363  for (PartSet_t::const_iterator i = _vertex_set.begin();
364  i != _vertex_set.end(); ++i)
365  {
366 
367  PHG4VtxPoint * v = truthInfoList->GetVtx(*i);
368  if (!v)
369  {
370  cout
371  << "PHG4DSTReader::process_event - ERROR - can not find vertex ID "
372  << *i << " in G4TruthInfo" << endl;
373 
374  continue;
375  }
376 
377  new ((*(rec._arr.get()))[rec._cnt]) vertex_type();
378 
379  if (Verbosity() >= 2)
380  cout << "PHG4DSTReader::process_event - saving vertex id "
381  << *i << endl;
382 
383  vertex_type * new_v =
384  static_cast<vertex_type *>(rec._arr.get()->At(rec._cnt));
385  assert(new_v);
386 
387  new_v->set_x(v->get_x());
388  new_v->set_y(v->get_y());
389  new_v->set_z(v->get_z());
390  new_v->set_t(v->get_t());
391  new_v->set_id(v->get_id());
392 
393  rec._cnt++;
394  } // else if (rec._type == record::typ_vertex)
395 
396  } // if (_load_all_particle)
397 
398  } // for (records_t::iterator it = _records.begin(); it != _records.end(); ++it)
399 
400  if (_T)
401  _T->Fill();
402 
403  return 0;
404 } // for (records_t::iterator it = _records.begin(); it != _records.end(); ++it)
405 
406 void
408 {
409 
410  assert(part);
411 
412 // if (Verbosity() >= 2)
413 // cout << "PHG4DSTReader::process_event - Particle type is "
414 // << p_raw->GetName() << " - " << p_raw->ClassName()
415 // << " get_track_id = " << p_raw->get_track_id() << endl;
416 
417 // part_type * part = dynamic_cast<part_type *>(p_raw);
418 
419 // assert(part);
420 
421  new ((*(rec._arr.get()))[rec._cnt]) part_type();
422 
423  part_type * new_part = static_cast<part_type *>(rec._arr.get()->At(rec._cnt));
424  assert(new_part);
425 
426  new_part->set_track_id(part->get_track_id());
427  new_part->set_vtx_id(part->get_vtx_id());
428  new_part->set_parent_id(part->get_parent_id());
429  new_part->set_primary_id(part->get_primary_id());
430  new_part->set_name(part->get_name());
431  new_part->set_pid(part->get_pid());
432  new_part->set_px(part->get_px());
433  new_part->set_py(part->get_py());
434  new_part->set_pz(part->get_pz());
435  new_part->set_e(part->get_e());
436 
437  _vertex_set.insert(part->get_vtx_id());
438 
439  rec._cnt++;
440 }
441 
442 int
444 {
445  cout << "PHG4DSTReader::End - Clean ups" << endl;
446 
447  if (_T)
448  {
450  _T->Write();
451  _T->ResetBranchAddresses();
452  }
453 
454  _records.clear();
455 
456  return 0;
457 }
PHG4VtxPointv1 vertex_type
PHG4HitEval hit_type
PHG4Particlev2 part_type
TFile clean handling.
virtual int Verbosity() const
Gets the verbosity of this module.
Definition: Fun4AllBase.h:64
virtual ~PHG4DSTReader()
bool _save_vertex
save vertex for particles?
int Init(PHCompositeNode *)
full initialization
bool _load_active_particle
load all particle that produced a saved hit
std::vector< std::string > _node_postfix
void add_particle(record &rec, PHG4Particle *part)
add a particle and associated vertex if _save_vertex
int process_event(PHCompositeNode *)
event processing method
PartSet_t _vertex_set
PartSet_t _particle_set
bool _load_all_particle
load all particle in truth info module?
PHG4DSTReader(const std::string &filename)
records_t _records
bool _save_particle
master switch to save particles
int End(PHCompositeNode *)
end of run method
std::string _out_file_name
Map::const_iterator ConstIterator
unsigned int size(void) const
ConstRange getHits(const unsigned int detid) const
return all hits matching a given detid
std::pair< ConstIterator, ConstIterator > ConstRange
PHG4HitEval for evaluating PHG4Hitv1 in CINT readable evaluation trees.
Definition: PHG4HitEval.h:20
virtual PHG4HitDefs::keytype get_hit_id() const
Definition: PHG4Hit.h:36
virtual int get_trkid() const
Definition: PHG4Hit.h:41
virtual int get_track_id() const
Definition: PHG4Particle.h:21
virtual int get_pid() const
Definition: PHG4Particle.h:14
virtual int get_parent_id() const
Definition: PHG4Particle.h:23
virtual double get_px() const
Definition: PHG4Particle.h:16
virtual int get_primary_id() const
Definition: PHG4Particle.h:24
virtual int get_vtx_id() const
Definition: PHG4Particle.h:22
virtual std::string get_name() const
Definition: PHG4Particle.h:15
virtual double get_py() const
Definition: PHG4Particle.h:17
virtual double get_e() const
Definition: PHG4Particle.h:19
virtual double get_pz() const
Definition: PHG4Particle.h:18
void set_pid(const int i)
void set_pz(const double x)
void set_px(const double x)
void set_name(const std::string &name)
void set_py(const double x)
void set_primary_id(const int i)
void set_parent_id(const int i)
void set_e(const double e)
void set_track_id(const int i)
void set_vtx_id(const int i)
std::pair< ConstIterator, ConstIterator > ConstRange
PHG4VtxPoint * GetVtx(const int vtxid)
Map::const_iterator ConstIterator
const Map & GetMap() const
Get the Particle Map storage.
virtual double get_y() const
Definition: PHG4VtxPoint.h:21
virtual double get_x() const
Definition: PHG4VtxPoint.h:20
virtual int get_id() const
Definition: PHG4VtxPoint.h:24
virtual double get_z() const
Definition: PHG4VtxPoint.h:22
virtual double get_t() const
Definition: PHG4VtxPoint.h:23
void set_y(const double r)
void set_t(const double r)
void set_z(const double r)
void set_x(const double r)
void set_id(const int i)
static PHTFileServer & get(void)
return reference to class singleton
Definition: PHTFileServer.h:36
void open(const std::string &filename, const std::string &type="RECREATE")
open a SafeTFile. If filename is not found in the map, create a new TFile and append to the map; incr...
bool cd(const std::string &filename)
change to directory of TFile matching filename
T * getClass(PHCompositeNode *top, const std::string &name)
Definition: getClass.h:18
TClonesArray * _arr_ptr