summaryrefslogtreecommitdiff
path: root/io/root.C
blob: 43f4e730cab71a31b6c7e62b85873b3d60851372 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
#include <TVector3.h>
#include <vector>
#include <TTree.h>
#include <string>

struct Photon {
  double t;
  TVector3 pos;
  TVector3 dir;
  TVector3 pol;
  double wavelength; // nm
  int history;
  int last_hit_triangle;
};

struct MC {
  std::string particle;
  TVector3 gen_pos;
  TVector3 gen_dir;
  double gen_total_e;

  std::vector<Photon> photon_start;
  std::vector<Photon> photon_stop;

};

struct Channel {
  Channel() : channel_id(-1), t(-9999.0), q(-9999.0) { };
  int channel_id;
  double t;
  double q;
};
  
struct Event {
  int event_id;
  MC mc;
  int nhit;
  std::vector<Channel> channel;

};

void fill_photons(Event *ev, bool start, 
		  unsigned int nphotons, float *pos, float *dir,
		  float *pol, float *wavelength, float *t0,
		  int *histories=0, int *last_hit_triangle=0)
{
  std::vector<Photon> &photons = start ? ev->mc.photon_start : ev->mc.photon_stop;
  photons.resize(nphotons);
  
  for (unsigned int i=0; i < nphotons; i++) {
    Photon &photon = photons[i];
    photon.t = t0[i];
    photon.pos.SetXYZ(pos[3*i], pos[3*i + 1], pos[3*i + 2]);
    photon.dir.SetXYZ(dir[3*i], dir[3*i + 1], dir[3*i + 2]);
    photon.pol.SetXYZ(pol[3*i], pol[3*i + 1], pol[3*i + 2]);
    photon.wavelength = wavelength[i];
    if (histories)
      photon.history = histories[i];
    else
      photon.history = 0;

    if (last_hit_triangle)
      photon.last_hit_triangle = last_hit_triangle[i];
    else
      photon.last_hit_triangle = -1;
  }
}
		  

void fill_hits(Event *ev, unsigned int nchannels, float *channel_times)
{
  ev->channel.resize(0);
  ev->nhit = 0;
  Channel ch;
  for (unsigned int i=0; i < nchannels; i++) {
    if (channel_times[i] < 1e8) {
      ev->nhit++;
      ch.channel_id = i;
      ch.t = channel_times[i];
      ch.q = 1.0;
      ev->channel.push_back(ch);
    }
  }
}


#ifdef __MAKECINT__
#pragma link C++ class vector<Photon>;
#pragma link C++ class vector<Channel>;
#endif