SND@LHC Software
Loading...
Searching...
No Matches
sndPlaneTools.cxx
Go to the documentation of this file.
1#include "sndPlaneTools.h"
2
3#include <vector>
4#include <stdexcept>
5
6#include "TClonesArray.h"
7#include "Scifi.h"
8#include "MuFilter.h"
9#include "sndConfiguration.h"
10#include "sndVetoPlane.h"
11#include "sndScifiPlane.h"
12#include "sndUSPlane.h"
13#include "sndDSPlane.h"
14#include "sndScifiHit.h"
15#include "MuFilterHit.h"
16
17std::vector<snd::analysis_tools::VetoPlane> snd::analysis_tools::FillVeto(const snd::Configuration &configuration, TClonesArray *mufi_hits, MuFilter *mufilter_geometry)
18{
19
20 std::vector<snd::analysis_tools::VetoPlane> veto_planes;
21 const int n_mufi_hits{mufi_hits->GetEntries()};
22
23 const int n_station = configuration.veto_n_stations;
24 std::vector<std::vector<MuFilterHit*>> plane_hits(n_station);
25
26 for (int i{0}; i < n_mufi_hits; ++i) {
27 auto hit = static_cast<MuFilterHit*>(mufi_hits->At(i));
28 if (hit->GetSystem()!=1) continue;
29 int station_id = hit->GetPlane();
30 if (station_id > -1 && station_id < n_station) {
31 plane_hits[station_id].push_back(hit);
32 }
33 else throw std::runtime_error{"Invalid Veto plane"};
34 }
35 for (int st{0}; st < n_station; ++st) {
36 veto_planes.emplace_back(snd::analysis_tools::VetoPlane(plane_hits[st], configuration, mufilter_geometry, st+1));
37 }
38 return veto_planes;
39}
40
41std::vector<snd::analysis_tools::ScifiPlane> snd::analysis_tools::FillScifi(const snd::Configuration &configuration, TClonesArray *sf_hits, Scifi *scifi_geometry)
42{
43
44 std::vector<snd::analysis_tools::ScifiPlane> scifi_planes;
45 const int n_sf_hits{sf_hits->GetEntries()};
46
47 const int max_station = configuration.scifi_n_stations;
48 std::vector<std::vector<sndScifiHit*>> stations_hits(max_station);
49
50 for (int i{0}; i < n_sf_hits; ++i) {
51 auto hit = static_cast<sndScifiHit*>(sf_hits->At(i));
52 int station_id = hit->GetStation()-1;
53
54 if (station_id > -1 && station_id < max_station) {
55 stations_hits[station_id].push_back(hit);
56 }
57 else throw std::runtime_error{"Invalid SciFi station"};
58 }
59 for (int st{0}; st < max_station; ++st) {
60 scifi_planes.emplace_back(snd::analysis_tools::ScifiPlane(stations_hits[st], configuration, scifi_geometry, st+1));
61 }
62 return scifi_planes;
63}
64
65
66std::vector<snd::analysis_tools::USPlane> snd::analysis_tools::FillUS(const snd::Configuration &configuration, TClonesArray *mufi_hits, MuFilter *mufilter_geometry, bool use_small_sipms)
67{
68
69 std::vector<snd::analysis_tools::USPlane> us_planes;
70 const int n_mufi_hits{mufi_hits->GetEntries()};
71
72 const int n_station = configuration.us_n_stations;
73 std::vector<std::vector<MuFilterHit*>> plane_hits(n_station);
74
75 for (int i{0}; i < n_mufi_hits; ++i) {
76 auto hit = static_cast<MuFilterHit*>(mufi_hits->At(i));
77 if (hit->GetSystem()!=2) continue;
78 int station_id = hit->GetPlane();
79 if (station_id > -1 && station_id < n_station) {
80 plane_hits[station_id].push_back(hit);
81 }
82 else throw std::runtime_error{"Invalid US plane"};
83 }
84 for (int st{0}; st < n_station; ++st) {
85 us_planes.emplace_back(snd::analysis_tools::USPlane(plane_hits[st], configuration, mufilter_geometry, st+1, use_small_sipms));
86 }
87 return us_planes;
88}
89
90std::vector<snd::analysis_tools::DSPlane> snd::analysis_tools::FillDS(const snd::Configuration &configuration, TClonesArray *mufi_hits, MuFilter *mufilter_geometry)
91{
92
93 std::vector<snd::analysis_tools::DSPlane> ds_planes;
94 const int n_mufi_hits{mufi_hits->GetEntries()};
95
96 const int n_station = configuration.ds_n_stations;
97 std::vector<std::vector<MuFilterHit*>> plane_hits(n_station);
98
99 for (int i{0}; i < n_mufi_hits; ++i) {
100 auto hit = static_cast<MuFilterHit*>(mufi_hits->At(i));
101 if (hit->GetSystem()!=3) continue;
102 int station_id = hit->GetPlane();
103 if (station_id > -1 && station_id < n_station) {
104 plane_hits[station_id].push_back(hit);
105 }
106 else throw std::runtime_error{"Invalid DS plane"};
107 }
108 for (int st{0}; st < n_station; ++st) {
109 ds_planes.emplace_back(snd::analysis_tools::DSPlane(plane_hits[st], configuration, mufilter_geometry, st+1));
110 }
111 return ds_planes;
112}
int GetPlane()
Definition MuFilterHit.h:42
Definition Scifi.h:20
Int_t GetStation()
Definition sndScifiHit.h:31
std::vector< ScifiPlane > FillScifi(const Configuration &configuration, TClonesArray *sf_hits, Scifi *scifi_geometry)
std::vector< VetoPlane > FillVeto(const Configuration &configuration, TClonesArray *mufi_hits, MuFilter *mufilter_geometry)
std::vector< DSPlane > FillDS(const Configuration &configuration, TClonesArray *mufi_hits, MuFilter *mufilter_geometry)
std::vector< USPlane > FillUS(const Configuration &configuration, TClonesArray *mufi_hits, MuFilter *mufilter_geometry, bool use_small_sipms=false)