Line | Branch | Exec | Source |
---|---|---|---|
1 | #pragma once | ||
2 | |||
3 | #include "iguana/algorithms/Algorithm.h" | ||
4 | #include "iguana/algorithms/TypeDefs.h" | ||
5 | |||
6 | namespace iguana::clas12 { | ||
7 | |||
8 | /// @brief_algo Find the sector for all rows in `REC::Particle` | ||
9 | /// | ||
10 | /// @begin_doc_algo{clas12::SectorFinder | Creator} | ||
11 | /// @input_banks{REC::Particle, REC::Track, REC::Calorimeter, REC::Scintillator} | ||
12 | /// @output_banks{%REC::Particle::Sector} | ||
13 | /// @end_doc | ||
14 | /// | ||
15 | /// @begin_doc_config{clas12/SectorFinder} | ||
16 | /// @config_param{bank_charged | string | if not `default`, use this bank for sector finding of charged particles} | ||
17 | /// @config_param{bank_uncharged | string | if not `default`, use this bank for sector finding of neutral particles} | ||
18 | /// @end_doc | ||
19 | /// | ||
20 | /// If `bank_charged` and/or `bank_uncharged` is default, then all of the following banks are needed, in addition to `REC::Particle`: | ||
21 | /// | ||
22 | /// - `REC::Track` | ||
23 | /// - `REC::Calorimeter` | ||
24 | /// - `REC::Scintillator` | ||
25 | /// | ||
26 | /// Otherwise only the bank(s) specified by `bank_charged` and `bank_uncharged` is/are needed, if both of them are non-default. | ||
27 | /// | ||
28 | /// The action function ::GetStandardSector identifies the sector(s) using these banks in a priority order, whereas | ||
29 | /// the action function ::GetSector uses a single bank's data. | ||
30 | /// Note: rows that have been filtered out of `REC::Particle` will still have their sectors determined. | ||
31 | /// | ||
32 | /// @creator_note | ||
33 | class SectorFinder : public Algorithm | ||
34 | { | ||
35 | |||
36 |
7/18✓ Branch 0 (2→3) taken 6 times.
✗ Branch 1 (2→6) not taken.
✗ Branch 2 (5→7) not taken.
✓ Branch 3 (5→8) taken 6 times.
✓ Branch 4 (9→10) taken 6 times.
✗ Branch 5 (9→33) not taken.
✓ Branch 6 (10→11) taken 6 times.
✗ Branch 7 (10→17) not taken.
✓ Branch 8 (17→18) taken 6 times.
✗ Branch 9 (17→43) not taken.
✓ Branch 10 (18→19) taken 6 times.
✗ Branch 11 (18→41) not taken.
✓ Branch 12 (25→26) taken 6 times.
✗ Branch 13 (25→41) not taken.
✗ Branch 14 (33→34) not taken.
✗ Branch 15 (33→40) not taken.
✗ Branch 16 (43→44) not taken.
✗ Branch 17 (43→46) not taken.
|
58 | DEFINE_IGUANA_ALGORITHM(SectorFinder, clas12::SectorFinder) |
37 | |||
38 | public: | ||
39 | |||
40 | void Start(hipo::banklist& banks) override; | ||
41 | void Run(hipo::banklist& banks) const override; | ||
42 | void Stop() override; | ||
43 | |||
44 | /// @action_function{scalar creator} for a given particle with index `pindex_particle`, get its sector from | ||
45 | /// a detector bank's list of `sectors` and `pindices` (both must be ordered in the same way) | ||
46 | /// | ||
47 | /// @note this is done instead of finding the `pindex` in the bank directly, to have an action function | ||
48 | /// | ||
49 | /// @par Example | ||
50 | /// ```cpp | ||
51 | /// | ||
52 | /// //... Initialise algorithms & banks ... | ||
53 | /// | ||
54 | /// //For each event, do: | ||
55 | /// | ||
56 | /// std::vector<int> sectors; | ||
57 | /// std::vector<int> pindices | ||
58 | /// | ||
59 | /// //bank is a hipo::bank object from which we want to get the sectors | ||
60 | /// //for example the bank object related to REC::Calorimeter | ||
61 | /// for(auto const& row : bank.getRowList()) { | ||
62 | /// | ||
63 | /// int det=bank.getInt("detector",row); | ||
64 | /// | ||
65 | /// //NB: you should check you read from an FD detector | ||
66 | /// // e.g. det 7 is the ECAL (see/use `iguana::DetectorType` enum) | ||
67 | /// if(det==7){ | ||
68 | /// sectors.push_back(bank.getInt("sector", row)); | ||
69 | /// pindices.push_back(bank.getShort("pindex", row)); | ||
70 | /// } | ||
71 | /// } | ||
72 | /// | ||
73 | /// //partbank is a hipo::bank object related to REC::Particle | ||
74 | /// //algo_sector_finder is the iguana::clas12::SectorFinder object | ||
75 | /// for(auto const& row : partbank.getRowList()) { | ||
76 | /// int sector = algo_sector_finder.GetSector(sectors, pindices, row); | ||
77 | /// } | ||
78 | /// ``` | ||
79 | /// | ||
80 | /// @see ::GetStandardSector, which calls this method for detectors in a priority order | ||
81 | /// | ||
82 | /// @param sectors list of sectors in a detector bank | ||
83 | /// @param pindices list of pindices in a detector bank | ||
84 | /// @param pindex_particle index in `REC::Particle` bank for which to get sector | ||
85 | /// @returns sector for `pindex_particle` in list, -1 if `pindex_particle` not in inputted list | ||
86 | int GetSector( | ||
87 | std::vector<int> const& sectors, | ||
88 | std::vector<int> const& pindices, | ||
89 | int const& pindex_particle) const; | ||
90 | |||
91 | /// @action_function{scalar creator} for a given particle with index `pindex_particle`, get its sector from | ||
92 | /// using the standard method | ||
93 | /// | ||
94 | /// The following detectors' banks will be searched in order, and once the sector is found for any detector, it is returned: | ||
95 | /// | ||
96 | /// - `REC::Track`, using `sectors_track` and `pindices_track` | ||
97 | /// - `REC::Calorimeter`, using `sectors_cal` and `pindices_cal` | ||
98 | /// - `REC::Scintillator`, using `sectors_scint` and `pindices_scint` | ||
99 | /// | ||
100 | /// @see ::GetSector, which exemplifies using only one bank's lists of `sectors` and `pindices` | ||
101 | /// | ||
102 | /// @param sectors_track list of sectors in `REC::Track` | ||
103 | /// @param pindices_track list of pindices in `REC::Track` | ||
104 | /// @param sectors_cal list of sectors in `REC::Calorimeter` | ||
105 | /// @param pindices_cal list of pindices in `REC::Calorimeter` | ||
106 | /// @param sectors_scint list of sectors in `REC::Scintillator` | ||
107 | /// @param pindices_scint list of pindices in `REC::Scintillator` | ||
108 | /// @param pindex_particle index in `REC::Particle` bank for which to get sector | ||
109 | /// @returns sector for `pindex_particle` in lists, -1 if `pindex_particle` not any of the inputted lists | ||
110 | int GetStandardSector( | ||
111 | std::vector<int> const& sectors_track, | ||
112 | std::vector<int> const& pindices_track, | ||
113 | std::vector<int> const& sectors_cal, | ||
114 | std::vector<int> const& pindices_cal, | ||
115 | std::vector<int> const& sectors_scint, | ||
116 | std::vector<int> const& pindices_scint, | ||
117 | int const& pindex_particle) const; | ||
118 | |||
119 | /// @action_function{vector creator} get sectors for all particles, using the standard method | ||
120 | /// | ||
121 | /// @overloads_scalar | ||
122 | /// | ||
123 | /// @see ::GetSector, which exemplifies using only one bank's lists of `sectors` and `pindices` | ||
124 | /// | ||
125 | /// @param sectors_track list of sectors in `REC::Track` | ||
126 | /// @param pindices_track list of pindices in `REC::Track` | ||
127 | /// @param sectors_cal list of sectors in `REC::Calorimeter` | ||
128 | /// @param pindices_cal list of pindices in `REC::Calorimeter` | ||
129 | /// @param sectors_scint list of sectors in `REC::Scintillator` | ||
130 | /// @param pindices_scint list of pindices in `REC::Scintillator` | ||
131 | /// @param pindices_particle the `REC::Particle` list of `pindices` | ||
132 | /// @returns list of sectors for each particle with `pindex` in `pindices_particle` | ||
133 | std::vector<int> GetStandardSector( | ||
134 | std::vector<int> const& sectors_track, | ||
135 | std::vector<int> const& pindices_track, | ||
136 | std::vector<int> const& sectors_cal, | ||
137 | std::vector<int> const& pindices_cal, | ||
138 | std::vector<int> const& sectors_scint, | ||
139 | std::vector<int> const& pindices_scint, | ||
140 | std::vector<int> const& pindices_particle) const; | ||
141 | |||
142 | /// fill lists of sectors and pindices present in the input bank | ||
143 | /// | ||
144 | /// @note this is not an action function, but here for convenience | ||
145 | /// | ||
146 | /// @param bank bank from which to get lists of sectors and pindices | ||
147 | /// @param sectors list to fill with sectors in the bank | ||
148 | /// @param pindices list to fill with pindices in the bank | ||
149 | void GetListsSectorPindex(hipo::bank const& bank, std::vector<int>& sectors, std::vector<int>& pindices) const; | ||
150 | |||
151 | private: | ||
152 | |||
153 | /// `hipo::banklist` index for the particle bank | ||
154 | hipo::banklist::size_type b_particle; | ||
155 | hipo::banklist::size_type b_calorimeter; | ||
156 | hipo::banklist::size_type b_track; | ||
157 | hipo::banklist::size_type b_scint; | ||
158 | hipo::banklist::size_type b_user_charged; | ||
159 | hipo::banklist::size_type b_user_uncharged; | ||
160 | hipo::banklist::size_type b_result; | ||
161 | bool userSpecifiedBank_charged{false}; | ||
162 | bool userSpecifiedBank_uncharged{true}; | ||
163 | |||
164 | // `b_result` bank item indices | ||
165 | int i_sector; | ||
166 | int i_pindex; | ||
167 | |||
168 | /// Configuration options | ||
169 | std::string o_bankname_charged; | ||
170 | std::string o_bankname_uncharged; | ||
171 | |||
172 | //only want sectors from FD detectors | ||
173 | std::set<int> const listFDDets{ | ||
174 | DetectorType::DC, | ||
175 | DetectorType::ECAL, | ||
176 | DetectorType::FTOF, | ||
177 | DetectorType::HTCC, | ||
178 | DetectorType::LTCC, | ||
179 | DetectorType::RICH | ||
180 | }; | ||
181 | }; | ||
182 | |||
183 | } | ||
184 |