Project
Loading...
Searching...
No Matches
AODMcProducerHelpers.cxx
Go to the documentation of this file.
1// Copyright 2023-2099 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
18#include <MathUtils/Utils.h>
19#include <algorithm>
20#define verbosity debug
21
22namespace o2::aodmchelpers
23{
24//==================================================================
26 const std::vector<std::string>& keys,
27 bool anyNotAll)
28{
29 auto check = [&header](const std::string& key) { // Do not format
30 return header.hasInfo(key);
31 };
32 return (anyNotAll ? // Do not format
33 std::any_of(keys.cbegin(), keys.cend(), check)
34 : // Do not format
35 std::all_of(keys.cbegin(), keys.cend(), check));
36}
37//====================================================================
39 int bcId,
40 float time,
42 short generatorId,
43 int sourceId,
44 unsigned int mask)
45{
47 using GenProp = o2::mcgenid::GeneratorProperty;
48 using namespace o2::math_utils;
49
50 int subGenId = getEventInfo(header, GenProp::SUBGENERATORID, -1);
51 int genId = getEventInfo(header, GenProp::GENERATORID,
52 int(generatorId));
53 float weight = getEventInfo(header, Key::weight, 1.f);
54
55 auto encodedGeneratorId = o2::mcgenid::getEncodedGenId(genId,
56 sourceId,
57 subGenId);
58
59 LOG(verbosity) << "Updating MC Collisions table w/bcId=" << bcId;
60 cursor(0,
61 bcId,
62 encodedGeneratorId,
63 truncateFloatFraction(header.GetX(), mask),
64 truncateFloatFraction(header.GetY(), mask),
65 truncateFloatFraction(header.GetZ(), mask),
66 truncateFloatFraction(time, mask),
67 truncateFloatFraction(weight, mask),
68 header.GetB(),
69 getEventInfo(header, Key::planeAngle, header.GetRotZ()));
70 return encodedGeneratorId;
71}
72//--------------------------------------------------------------------
74 int collisionID,
75 short generatorID,
77 HepMCUpdate when)
78{
80
81 if (when == HepMCUpdate::never or
82 (when != HepMCUpdate::always and
83 not hasKeys(header, { // Do not
84 Key::acceptedEvents, // mess with
85 Key::attemptedEvents, // with
86 Key::xSection, // my
87 Key::xSectionError}, // formatting
88 when == HepMCUpdate::anyKey))) {
89 return false;
90 }
91
92 LOG(verbosity) << "Updating HepMC cross-section table w/collisionId "
93 << collisionID;
94 cursor(0,
95 collisionID,
96 generatorID,
97 getEventInfo(header, Key::acceptedEvents, 0),
98 getEventInfo(header, Key::attemptedEvents, 0),
99 getEventInfo(header, Key::xSection, 0.f),
100 getEventInfo(header, Key::xSectionError, 0.f),
101 getEventInfo(header, Key::eventScale, 1.f),
102 getEventInfo(header, Key::mpi, -1),
103 getEventInfo(header, Key::processCode, -1));
104 return true;
105}
106//--------------------------------------------------------------------
108 int collisionID,
109 short generatorID,
110 o2::dataformats::MCEventHeader const& header,
111 HepMCUpdate when)
112{
114
115 if (when == HepMCUpdate::never or
116 (when != HepMCUpdate::always and // Do
117 not hasKeys(header, {Key::pdfParton1Id, // not
118 Key::pdfParton2Id, // mess
119 Key::pdfCode1, // with
120 Key::pdfCode2, // my
121 Key::pdfX1, // formatting
122 Key::pdfX2, // .
123 Key::pdfScale, // It
124 Key::pdfXF1, // is
125 Key::pdfXF2}, // better
126 when == HepMCUpdate::anyKey))) {
127 return false;
128 }
129
130 LOG(verbosity) << "Updating HepMC PDF table w/collisionId " << collisionID;
131 cursor(0,
132 collisionID,
133 generatorID,
134 getEventInfo(header, Key::pdfParton1Id, int(0)),
135 getEventInfo(header, Key::pdfParton2Id, int(0)),
136 getEventInfo(header, Key::pdfCode1, int(0)),
137 getEventInfo(header, Key::pdfCode2, int(0)),
138 getEventInfo(header, Key::pdfX1, 0.f),
139 getEventInfo(header, Key::pdfX2, 0.f),
140 getEventInfo(header, Key::pdfScale, 1.f),
141 getEventInfo(header, Key::pdfXF1, 0.f),
142 getEventInfo(header, Key::pdfXF2, 0.f));
143
144 return true;
145}
146//--------------------------------------------------------------------
148 int collisionID,
149 short generatorID,
150 o2::dataformats::MCEventHeader const& header,
151 HepMCUpdate when)
152{
154
155 if (when == HepMCUpdate::never or
156 (when != HepMCUpdate::always and // clang
157 not hasKeys(header, {Key::nCollHard, // format
158 Key::nPartProjectile, // is
159 Key::nPartTarget, // so
160 Key::nColl, // annoying
161 Key::nCollNNWounded, // .
162 Key::nCollNWoundedN, // It
163 Key::nCollNWoundedNwounded, // messes
164 Key::nSpecProjectileNeutron, // up
165 Key::nSpecTargetNeutron, // the
166 Key::nSpecProjectileProton, // clarity
167 Key::nSpecTargetProton, // of
168 Key::planeAngle, // the
169 "eccentricity", // code
170 Key::sigmaInelNN, // to
171 Key::centrality}, // noavail
172 when == HepMCUpdate::anyKey))) {
173 return false;
174 }
175
176 int specNeutrons = (getEventInfo(header, Key::nSpecProjectileNeutron, -1) +
177 getEventInfo(header, Key::nSpecTargetNeutron, -1));
178 int specProtons = (getEventInfo(header, Key::nSpecProjectileProton, -1) +
179 getEventInfo(header, Key::nSpecTargetProton, -1));
180
181 LOG(verbosity) << "Updating HepMC heavy-ion table w/collisionID "
182 << collisionID;
183 cursor(0,
184 collisionID,
185 generatorID,
186 getEventInfo(header, Key::nCollHard, -1),
187 getEventInfo(header, Key::nPartProjectile, -1),
188 getEventInfo(header, Key::nPartTarget, -1),
189 getEventInfo(header, Key::nColl, -1),
190 getEventInfo(header, Key::nCollNNWounded, -1),
191 getEventInfo(header, Key::nCollNWoundedN, -1),
192 getEventInfo(header, Key::nCollNWoundedNwounded, -1),
193 specNeutrons,
194 specProtons,
195 header.GetB(),
196 getEventInfo(header, Key::planeAngle, header.GetRotZ()),
197 getEventInfo(header, "eccentricity", 0),
198 getEventInfo(header, Key::sigmaInelNN, 0.),
199 getEventInfo(header, Key::centrality, -1));
200
201 return true;
202}
203//--------------------------------------------------------------------
205 const TrackToIndex& toStore,
206 int collisionID,
207 o2::MCTrack const& track,
208 std::vector<MCTrack> const& tracks,
209 uint8_t flags,
210 uint32_t weightMask,
211 uint32_t momentumMask,
212 uint32_t positionMask)
213{
215 using namespace o2::aod::mcparticle::enums;
216 using namespace o2::math_utils;
217 using namespace o2::mcgenstatus;
218
219 auto mapping = [&toStore](int trackNo) {
220 auto iter = toStore.find(trackNo);
221 if (iter == toStore.end()) {
222 return -1;
223 }
224 return iter->second;
225 };
226
227 auto statusCode = track.getStatusCode().fullEncoding;
228 auto hepmc = getHepMCStatusCode(track.getStatusCode());
229 if (not track.isPrimary()) {
230 flags |= ProducedByTransport;
231 statusCode = track.getProcess();
232 }
233 if (MCTrackNavigator::isPhysicalPrimary(track, tracks)) {
234 flags |= PhysicalPrimary;
235 }
236
237 int daughters[2] = {-1, -1};
238 std::vector<int> mothers;
239 int id;
240 if ((id = mapping(track.getMotherTrackId())) >= 0) {
241 mothers.push_back(id);
242 }
243 if ((id = mapping(track.getSecondMotherTrackId())) >= 0) {
244 mothers.push_back(id);
245 }
246 if ((id = mapping(track.getFirstDaughterTrackId())) >= 0) {
247 daughters[0] = id;
248 }
249 if ((id = mapping(track.getLastDaughterTrackId())) >= 0) {
250 daughters[1] = id;
251 } else {
252 daughters[1] = daughters[0];
253 }
254 if (daughters[0] < 0 and daughters[1] >= 0) {
255 LOG(error) << "Problematic daughters: " << daughters[0] << " and "
256 << daughters[1];
257 daughters[0] = daughters[1];
258 }
259 if (daughters[0] > daughters[1]) {
260 std::swap(daughters[0], daughters[1]);
261 }
262
263 float weight = track.getWeight();
264 float pX = float(track.Px());
265 float pY = float(track.Py());
266 float pZ = float(track.Pz());
267 float energy = float(track.GetEnergy());
268 float vX = float(track.Vx());
269 float vY = float(track.Vy());
270 float vZ = float(track.Vz());
271 float time = float(track.T());
272
273 LOG(verbosity) << "McParticle collisionId=" << collisionID << ","
274 << "status=" << statusCode << ","
275 << "hepmc=" << hepmc << ","
276 << "pdg=" << track.GetPdgCode() << ","
277 << "nMothers=" << mothers.size() << ","
278 << "daughters=" << daughters[0] << ","
279 << daughters[1];
280
281 cursor(0,
282 collisionID,
283 track.GetPdgCode(),
284 statusCode,
285 flags,
286 mothers,
287 daughters,
288 truncateFloatFraction(weight, weightMask),
289 truncateFloatFraction(pX, momentumMask),
290 truncateFloatFraction(pY, momentumMask),
291 truncateFloatFraction(pZ, momentumMask),
292 truncateFloatFraction(energy, momentumMask),
293 truncateFloatFraction(vX, positionMask),
294 truncateFloatFraction(vY, positionMask),
295 truncateFloatFraction(vZ, positionMask),
296 truncateFloatFraction(time, positionMask));
297}
298//--------------------------------------------------------------------
299uint32_t updateParticles(const ParticleCursor& cursor,
300 int collisionID,
301 std::vector<MCTrack> const& tracks,
302 TrackToIndex& preselect,
303 uint32_t offset,
304 bool filter,
305 bool background,
306 uint32_t weightMask,
307 uint32_t momentumMask,
308 uint32_t positionMask)
309{
311 using namespace o2::aod::mcparticle::enums;
312 using namespace o2::mcgenstatus;
313
314 // First loop over particles to find out which to store
315 // TrackToIndex toStore(preselect.begin(), preselect.end());
316 //
317 // Guess we need to modifiy the passed in mapping so that MC labels
318 // can be set correctly
319 TrackToIndex& toStore = preselect;
320
321 // Mapping lambda. This maps the track number to the index into
322 // the table exported.
323 auto mapping = [&toStore](int trackNo) {
324 auto iter = toStore.find(trackNo);
325 if (iter == toStore.end()) {
326 return -1;
327 }
328 return iter->second;
329 };
330
331 LOG(verbosity) << "Got a total of " << tracks.size();
332 for (int trackNo = tracks.size() - 1; trackNo >= 0; trackNo--) {
333 auto& track = tracks[trackNo];
334 auto hepmc = getHepMCStatusCode(track.getStatusCode());
335 if (filter) {
336 if (toStore.find(trackNo) == toStore.end() and
337 /* The above test is in-correct. The track may be stored in
338 the list, but with a negative value. In that case, the
339 above test will still check mothers and daughters, and
340 possible store them in the output. This is however how
341 it is done in current `dev` branch, and so to enable
342 comparison on closure, we do this test for now. The
343 correct way it commented out below. */
344 // mapping(trackNo) < 0 and
345 not track.isPrimary() and
346 not MCTrackNavigator::isPhysicalPrimary(track, tracks) and
347 not MCTrackNavigator::isKeepPhysics(track, tracks)) {
348 LOG(verbosity) << "Skipping track " << trackNo << " " << hepmc << " "
349 << mapping(trackNo) << " "
350 << track.isPrimary() << " "
351 << MCTrackNavigator::isPhysicalPrimary(track, tracks)
352 << " "
353 << MCTrackNavigator::isKeepPhysics(track, tracks);
354 continue;
355 }
356 }
357
358 // Store this particle. We mark that putting a 1 in the
359 // `toStore` mapping. This will later on be updated with the
360 // actual index into the table
361 toStore[trackNo] = 1;
362
363 // If we're filtering tracks, then also mark mothers and
364 // daughters(?) to be stored.
365 if (filter) {
366 int id;
367 if ((id = track.getMotherTrackId()) >= 0) {
368 toStore[id] = 1;
369 }
370 if ((id = track.getSecondMotherTrackId()) >= 0) {
371 toStore[id] = 1;
372 }
373 if ((id = track.getFirstDaughterTrackId()) >= 0) {
374 toStore[id] = 1;
375 }
376 if ((id = track.getLastDaughterTrackId()) >= 0) {
377 toStore[id] = 1;
378 }
379 }
380 }
381
382 // Second loop to set indexes. This is needed to be done before
383 // we actually update the table, because a particle may point to a
384 // later particle.
385 LOG(verbosity) << "Will store " << toStore.size() << " particles";
386 size_t index = 0;
387
388 for (size_t trackNo = 0U; trackNo < tracks.size(); trackNo++) {
389 auto storeIt = mapping(trackNo);
390 if (storeIt < 0) {
391 continue;
392 }
393
394 toStore[trackNo] = offset + index;
395 index++;
396 }
397
398 // Make sure we have the right number
399 assert(index == toStore.size());
400 LOG(verbosity) << "Starting index " << offset << ", last index "
401 << offset + index << " "
402 << "Selected " << toStore.size() << " tracks out of "
403 << tracks.size() << " "
404 << "Collision # " << collisionID;
405 // Third loop to actually store the particles in the order given
406 for (size_t trackNo = 0U; trackNo < tracks.size(); trackNo++) {
407 auto storeIt = mapping(trackNo);
408 if (storeIt < 0) {
409 continue;
410 }
411
412 auto& track = tracks[trackNo];
413 auto hepmc = getHepMCStatusCode(track.getStatusCode());
414 uint8_t flags = (background ? FromBackgroundEvent : 0);
415 updateParticle(cursor,
416 toStore,
417 collisionID,
418 track,
419 tracks,
420 flags,
421 weightMask,
422 momentumMask,
423 positionMask);
424 }
425 LOG(verbosity) << "Return new offset " << offset + index;
426 return offset + index;
427}
428} // namespace o2::aodmchelpers
429
430// Local Variables:
431// mode: C++
432// End:
#define verbosity
o2::monitoring::tags::Key Key
General auxilliary methods.
int16_t time
Definition RawEventData.h:4
Utility functions for MC particles.
StringRef key
int getProcess() const
get the production process (id) of this track
Definition MCTrack.h:230
o2::mcgenstatus::MCGenStatusEncoding getStatusCode() const
get generator status code
Definition MCTrack.h:233
Int_t getFirstDaughterTrackId() const
Definition MCTrack.h:77
Double_t Vy() const
Definition MCTrack.h:105
Double_t Vx() const
Definition MCTrack.h:104
bool isPrimary() const
Definition MCTrack.h:75
Int_t getLastDaughterTrackId() const
Definition MCTrack.h:78
_T getWeight() const
return particle weight
Definition MCTrack.h:96
Double_t Vz() const
Definition MCTrack.h:106
Double_t Py() const
Definition MCTrack.h:102
Double_t GetEnergy() const
Definition MCTrack.h:304
Double_t T() const
Definition MCTrack.h:107
Int_t getSecondMotherTrackId() const
Definition MCTrack.h:74
Double_t Px() const
Definition MCTrack.h:101
Double_t Pz() const
Definition MCTrack.h:103
Int_t GetPdgCode() const
Accessors.
Definition MCTrack.h:72
Int_t getMotherTrackId() const
Definition MCTrack.h:73
bool hasInfo(std::string const &key) const
GLuint index
Definition glcorearb.h:781
GLuint GLuint GLfloat weight
Definition glcorearb.h:5477
GLintptr offset
Definition glcorearb.h:660
GLbitfield flags
Definition glcorearb.h:1570
GLint GLint GLint GLint GLint GLint GLint GLbitfield GLenum filter
Definition glcorearb.h:1308
GLint GLuint mask
Definition glcorearb.h:291
GLuint id
Definition glcorearb.h:650
TableCursor< aod::HepMCPdfInfos >::type PdfInfoCursor
bool updateHepMCHeavyIon(const HeavyIonCursor &cursor, int collisionID, short generatorID, o2::dataformats::MCEventHeader const &header, HepMCUpdate when=HepMCUpdate::anyKey)
bool hasKeys(o2::dataformats::MCEventHeader const &header, const std::vector< std::string > &keys, bool anyNotall=true)
TableCursor< aod::StoredMcParticles_001 >::type ParticleCursor
std::unordered_map< int, int > TrackToIndex
short updateMCCollisions(const CollisionCursor &cursor, int bcId, float time, o2::dataformats::MCEventHeader const &header, short generatorId=0, int sourceId=0, unsigned int mask=0xFFFFFFF0)
bool updateHepMCPdfInfo(const PdfInfoCursor &cursor, int collisionID, short generatorID, o2::dataformats::MCEventHeader const &header, HepMCUpdate when=HepMCUpdate::anyKey)
TableCursor< aod::HepMCHeavyIons >::type HeavyIonCursor
TableCursor< aod::HepMCXSections >::type XSectionCursor
uint32_t updateParticles(const ParticleCursor &cursor, int collisionID, std::vector< MCTrack > const &tracks, TrackToIndex &preselect, uint32_t offset=0, bool filter=false, bool background=false, uint32_t weightMask=0xFFFFFFF0, uint32_t momentumMask=0xFFFFFFF0, uint32_t positionMask=0xFFFFFFF0)
TableCursor< aod::McCollisions >::type CollisionCursor
bool updateHepMCXSection(const XSectionCursor &cursor, int collisionID, short generatorID, o2::dataformats::MCEventHeader const &header, HepMCUpdate when=HepMCUpdate::anyKey)
void updateParticle(const ParticleCursor &cursor, const TrackToIndex &toStore, int collisionID, o2::MCTrack const &track, std::vector< MCTrack > const &tracks, uint8_t flags=0, uint32_t weightMask=0xFFFFFFF0, uint32_t momentumMask=0xFFFFFFF0, uint32_t positionMask=0xFFFFFFF0)
const T getEventInfo(o2::dataformats::MCEventHeader const &header, std::string const &key, T const &def)
void check(const std::vector< std::string > &arguments, const std::vector< ConfigParamSpec > &workflowOptions, const std::vector< DeviceSpec > &deviceSpecs, CheckMatrix &matrix)
short getEncodedGenId(int generatorId, int sourceId, int subGeneratorId=-1)
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"