Project
Loading...
Searching...
No Matches
Detector.h
Go to the documentation of this file.
1// Copyright 2019-2020 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
14
15#ifndef ALICEO2_BASE_DETECTOR_H_
16#define ALICEO2_BASE_DETECTOR_H_
17
18#include <map>
19#include <tbb/concurrent_unordered_map.h>
20#include <vector>
21#include <initializer_list>
22#include <memory>
23
24#include "FairDetector.h" // for FairDetector
25#include "FairRootManager.h"
27#include "Rtypes.h" // for Float_t, Int_t, Double_t, Detector::Class, etc
28#include <cxxabi.h>
29#include <typeinfo>
30#include <type_traits>
31#include <string>
34#include <sys/shm.h>
35#include <type_traits>
36#include <unistd.h>
37#include <cassert>
38#include <list>
39#include <mutex>
40#include <thread>
41
42#include <fairmq/FwdDecls.h>
43
44namespace o2::base
45{
46
49class Detector : public FairDetector
50{
51
52 public:
53 Detector(const char* name, Bool_t Active);
54
56 Detector();
57
59 ~Detector() override;
60
61 // Module composition
62 void Material(Int_t imat, const char* name, Float_t a, Float_t z, Float_t dens, Float_t radl, Float_t absl,
63 Float_t* buf = nullptr, Int_t nwbuf = 0);
64
65 void Mixture(Int_t imat, const char* name, Float_t* a, Float_t* z, Float_t dens, Int_t nlmat,
66 Float_t* wmat);
67
68 void Medium(Int_t numed, const char* name, Int_t nmat, Int_t isvol, Int_t ifield, Float_t fieldm,
69 Float_t tmaxfd, Float_t stemax, Float_t deemax, Float_t epsil, Float_t stmin, Float_t* ubuf = nullptr,
70 Int_t nbuf = 0);
71
73 void SpecialCuts(Int_t numed, const std::initializer_list<std::pair<ECut, Float_t>>& parIDValMap);
75 void SpecialCut(Int_t numed, ECut parID, Float_t val);
76
77 void SpecialProcesses(Int_t numed, const std::initializer_list<std::pair<EProc, int>>& parIDValMap);
79 void SpecialProcess(Int_t numed, EProc parID, int val);
80
89 void Matrix(Int_t& nmat, Float_t theta1, Float_t phi1, Float_t theta2, Float_t phi2, Float_t theta3,
90 Float_t phi3) const;
91
92 static void setDensityFactor(Float_t density)
93 {
94 mDensityFactor = density;
95 }
96
98 {
99 return mDensityFactor;
100 }
101
104 void SetSpecialPhysicsCuts() override;
105
107 virtual void addAlignableVolumes() const;
108
110 virtual void fillParallelWorld() const;
111
113 virtual void defineWrapperVolume(Int_t id, Double_t rmin, Double_t rmax, Double_t zspan);
114
116 virtual void setNumberOfWrapperVolumes(Int_t n);
117
118 virtual void defineLayer(Int_t nlay, Double_t phi0, Double_t r, Int_t nladd, Int_t nmod,
119 Double_t lthick = 0., Double_t dthick = 0., UInt_t detType = 0, Int_t buildFlag = 0);
120
121 virtual void defineLayerTurbo(Int_t nlay, Double_t phi0, Double_t r, Int_t nladd, Int_t nmod,
122 Double_t width, Double_t tilt, Double_t lthick = 0., Double_t dthick = 0.,
123 UInt_t detType = 0, Int_t buildFlag = 0);
124
125 // returns global material ID given a "local" material ID for this detector
126 // returns -1 in case local ID not found
127 int getMaterialID(int imat) const
128 {
130 return mgr.getMaterialID(GetName(), imat);
131 }
132
133 // returns global medium ID given a "local" medium ID for this detector
134 // returns -1 in case local ID not found
135 int getMediumID(int imed) const
136 {
138 return mgr.getMediumID(GetName(), imed);
139 }
140
141 // fill the medium index mapping into a standard vector
142 // the vector gets sized properly and will be overridden
143 void getMediumIDMappingAsVector(std::vector<int>& mapping)
144 {
146 mgr.getMediumIDMappingAsVector(GetName(), mapping);
147 }
148
149 // return the name augmented by extension
150 std::string addNameTo(const char* ext) const
151 {
152 std::string s(GetName());
153 return s + ext;
154 }
155
156 // returning the name of the branch (corresponding to probe)
157 // returns zero length string when probe not defined
158 virtual std::string getHitBranchNames(int probe) const = 0;
159
160 // interface to update track indices of data objects
161 // usually called by the Stack, at the end of an event, which might have changed
162 // the track indices due to filtering
163 // FIXME: make private friend of stack?
164 virtual void updateHitTrackIndices(std::map<int, int> const&) = 0;
165
166 // interfaces to attach properly encoded hit information to a FairMQ message
167 // and to decode it
168 virtual void attachHits(fair::mq::Channel&, fair::mq::Parts&) = 0;
169 virtual void fillHitBranch(TTree& tr, fair::mq::Parts& parts, int& index) = 0;
170 virtual void collectHits(int eventID, fair::mq::Parts& parts, int& index) = 0;
171 virtual void mergeHitEntriesAndFlush(int eventID,
172 TTree& target,
173 std::vector<int> const& trackoffsets,
174 std::vector<int> const& nprimaries,
175 std::vector<int> const& subevtsOrdered) = 0;
176
177 // interface needed to merge together hit entries in TBranches (as used by hit merger process)
178 // trackoffsets: a map giving the corresponding trackoffset to be applied to the trackID property when
179 // merging
180 virtual void mergeHitEntries(TTree& origin, TTree& target, std::vector<int> const& trackoffsets, std::vector<int> const& nprimaries, std::vector<int> const& subevtsOrdered) = 0;
181
182 // hook which is called automatically to custom initialize the O2 detectors
183 // all initialization not able to do in constructors should be done here
184 // (typically the case for geometry related stuff, etc)
185 virtual void InitializeO2Detector() = 0;
186
187 // the original FairModule/Detector virtual Initialize function
188 // calls individual customized initializations and makes sure that the mother Initialize
189 // is called as well. Marked final for this reason!
190 void Initialize() final
191 {
193 // make sure the basic initialization is also done
194 FairDetector::Initialize();
195 }
196
197 // a second initialization method for stuff that should be initialized late
198 // (in our case after forking off from the main simulation setup
199 // ... for things that should be setup in each simulation worker separately)
200 virtual void initializeLate() = 0;
201
204 int registerSensitiveVolumeAndGetVolID(std::string const& name);
205
208 int registerSensitiveVolumeAndGetVolID(TGeoVolume const* vol);
209
210 // The GetCollection interface is made final and deprecated since
211 // we no longer support TClonesArrays
212 [[deprecated("Use getHits API on concrete detectors!")]] TClonesArray* GetCollection(int iColl) const final;
213
214 // static and reusable service function to set tracking parameters in relation to field
215 // returns global integration mode (inhomogenety) for the field and the max field value
216 // which is required for media creation
217 static void initFieldTrackingParams(int& mode, float& maxfield);
218
220 static void setDetId2HitBitIndex(std::vector<int> const& v) { Detector::sDetId2HitBitIndex = v; }
221 static std::vector<int> const& getDetId2HitBitIndex() { return Detector::sDetId2HitBitIndex; }
222
223 protected:
224 Detector(const Detector& origin);
225
226 Detector& operator=(const Detector&);
227
228 private:
233 std::map<int, int> mMapMaterial;
234
236 std::map<int, int> mMapMedium;
237
238 static Float_t mDensityFactor;
239 // systematic studies)
240 static std::vector<int> sDetId2HitBitIndex;
241
242 ClassDefOverride(Detector, 1); // Base class for ALICE Modules
243};
244
246inline std::string demangle(const char* name)
247{
248 int status = -4; // some arbitrary value to eliminate compiler warnings
249 std::unique_ptr<char, void (*)(void*)> res{abi::__cxa_demangle(name, nullptr, nullptr, &status), std::free};
250 return (status == 0) ? res.get() : name;
251}
252
253void attachShmMessage(void* hitsptr, fair::mq::Channel& channel, fair::mq::Parts& parts, bool* busy_ptr);
254void* decodeShmCore(fair::mq::Parts& dataparts, int index, bool*& busy);
255
256template <typename T>
257T decodeShmMessage(fair::mq::Parts& dataparts, int index, bool*& busy)
258{
259 return reinterpret_cast<T>(decodeShmCore(dataparts, index, busy));
260}
261
262// this goes into the source
263void attachMessageBufferToParts(fair::mq::Parts& parts, fair::mq::Channel& channel, void* data, TClass* cl);
264
265template <typename Container>
266void attachTMessage(Container const& hits, fair::mq::Channel& channel, fair::mq::Parts& parts)
267{
268 attachMessageBufferToParts(parts, channel, (void*)&hits, TClass::GetClass(typeid(hits)));
269}
270
271void* decodeTMessageCore(fair::mq::Parts& dataparts, int index);
272template <typename T>
273T decodeTMessage(fair::mq::Parts& dataparts, int index)
274{
275 return static_cast<T>(decodeTMessageCore(dataparts, index));
276}
277
278void attachDetIDHeaderMessage(int id, fair::mq::Channel& channel, fair::mq::Parts& parts);
279
280template <typename T>
281TBranch* getOrMakeBranch(TTree& tree, const char* brname, T* ptr)
282{
283 if (auto br = tree.GetBranch(brname)) {
284 br->SetAddress(static_cast<void*>(&ptr));
285 return br;
286 }
287 // otherwise make it
288 return tree.Branch(brname, ptr);
289}
290
291// a trait to determine if we should use shared mem or serialize using TMessage
292template <typename Det>
293struct UseShm {
294 static constexpr bool value = false;
295};
296
297// an implementation helper template which automatically implements
298// common functionality for deriving classes via the CRT pattern
299// (example: it implements the updateHitTrackIndices function and avoids
300// code duplication, while at the same time avoiding virtual function calls)
301template <typename Det>
303{
304 public:
305 // offer same constructors as base
306 using Detector::Detector;
307
308 // default implementation for getHitBranchNames
309 std::string getHitBranchNames(int probe) const override
310 {
311 if (probe == 0) {
312 return addNameTo("Hit");
313 }
314 return std::string(); // empty string as undefined
315 }
316
317 // generic implementation for the updateHitTrackIndices interface
318 // assumes Detectors have a GetHits(int) function that return some iterable
319 // hits which are o2::BaseHits
320 void updateHitTrackIndices(std::map<int, int> const& indexmapping) override
321 {
322 int probe = 0; // some Detectors have multiple hit vectors and we are probing
323 // them via a probe integer until we get a nullptr
324 while (auto hits = static_cast<Det*>(this)->Det::getHits(probe++)) {
325 for (auto& hit : *hits) {
326 auto iter = indexmapping.find(hit.GetTrackID());
327 hit.SetTrackID(iter->second);
328 }
329 }
330 }
331
332 void attachHits(fair::mq::Channel& channel, fair::mq::Parts& parts) override
333 {
334 int probe = 0;
335 // check if there is anything to be attached
336 // at least the first hit index should return non nullptr
337 if (static_cast<Det*>(this)->Det::getHits(0) == nullptr) {
338 return;
339 }
340
341 attachDetIDHeaderMessage(GetDetId(), channel, parts); // the DetId s are universal as they come from o2::detector::DetID
342
343 while (auto hits = static_cast<Det*>(this)->Det::getHits(probe++)) {
344 if (!UseShm<Det>::value || !o2::utils::ShmManager::Instance().isOperational()) {
345 attachTMessage(*hits, channel, parts);
346 } else {
347 // this is the shared mem variant
348 // we will just send the sharedmem ID and the offset inside
349 *mShmBusy[mCurrentBuffer] = true;
350 attachShmMessage((void*)hits, channel, parts, mShmBusy[mCurrentBuffer]);
351 }
352 }
353 }
354
355 // this merges several entries from the TBranch brname from the origin TTree
356 // into a single entry in a target TTree / same branch
357 // (assuming T is typically a vector; merging is simply done by appending)
358 // make this function a (static helper)
359 template <typename T>
360 void mergeAndAdjustHits(std::string const& brname, TTree& origin, TTree& target,
361 std::vector<int> const& trackoffsets, std::vector<int> const& nprimaries, std::vector<int> const& subevtsOrdered)
362 {
363 auto originbr = origin.GetBranch(brname.c_str());
364 if (originbr) {
365 auto targetdata = new T;
366 T* incomingdata = nullptr;
367 originbr->SetAddress(&incomingdata);
368
369 T* filladdress = nullptr;
370 if (origin.GetEntries() == 1) {
371 originbr->GetEntry(0);
372 filladdress = incomingdata;
373 } else {
374 Int_t entries = origin.GetEntries();
375 Int_t nprimTot = 0;
376 for (auto entry = 0; entry < entries; entry++) {
377 nprimTot += nprimaries[entry];
378 }
379 // offset for pimary track index
380 Int_t idelta0 = 0;
381 // offset for secondary track index
382 Int_t idelta1 = nprimTot;
383 for (int entry = entries - 1; entry >= 0; --entry) {
384 // proceed in the order of subevent Ids
385 Int_t index = subevtsOrdered[entry];
386 // numbe of primaries for this event
387 Int_t nprim = nprimaries[index];
388 idelta1 -= nprim;
389 filladdress = targetdata;
390 originbr->GetEntry(index);
391 if (incomingdata) {
392 // fix the trackIDs for this data
393 for (auto& hit : *incomingdata) {
394 const auto oldID = hit.GetTrackID();
395 // offset depends on whether the trackis a primary or secondary
396 Int_t offset = (oldID < nprim) ? idelta0 : idelta1;
397 hit.SetTrackID(oldID + offset);
398 }
399 // this could be further generalized by using a policy for T
400 std::copy(incomingdata->begin(), incomingdata->end(), std::back_inserter(*targetdata));
401 delete incomingdata;
402 incomingdata = nullptr;
403 }
404 // adjust offsets for next subevent
405 idelta0 += nprim;
406 idelta1 += trackoffsets[index];
407 } // subevent loop
408 }
409 // fill target for this event
410 auto targetbr = o2::base::getOrMakeBranch(target, brname.c_str(), &filladdress);
411 targetbr->SetAddress(&filladdress);
412 targetbr->Fill();
413 targetbr->ResetAddress();
414 targetdata->clear();
415 if (incomingdata) {
416 delete incomingdata;
417 incomingdata = nullptr;
418 }
419 delete targetdata;
420 }
421 }
422
423 // this merges several entries from temporary hit buffer into
424 // into a single entry in a target TTree / same branch
425 // (assuming T is typically a vector; merging is simply done by appending)
426 template <typename T, typename L>
427 void mergeAndAdjustHits(std::string const& brname, L& hitbuffervector, TTree& target,
428 std::vector<int> const& trackoffsets, std::vector<int> const& nprimaries,
429 std::vector<int> const& subevtsOrdered)
430 {
431 auto entries = hitbuffervector.size();
432
433 auto targetdata = new T; // used to collect data inside a single container
434 T* filladdress = nullptr; // pointer used for final ROOT IO
435 if (entries == 1) {
436 filladdress = hitbuffervector[0].get();
437 // nothing to do; we can directly do IO from the existing buffer
438 } else {
439 // here we need to do merging and index adjustment
440 int nprimTot = 0;
441 for (auto entry = 0; entry < entries; entry++) {
442 nprimTot += nprimaries[entry];
443 }
444 // offset for pimary track index
445 int idelta0 = 0;
446 // offset for secondary track index
447 int idelta1 = nprimTot;
448 filladdress = targetdata;
449 for (int entry = entries - 1; entry >= 0; --entry) {
450 // proceed in the order of subevent Ids
451 int index = subevtsOrdered[entry];
452 // number of primaries for this event
453 int nprim = nprimaries[index];
454 idelta1 -= nprim;
455
456 // fetch correct data item
457 auto incomingdata = hitbuffervector[index].get();
458 if (incomingdata) {
459 // fix the trackIDs for this data
460 for (auto& hit : *incomingdata) {
461 const auto oldID = hit.GetTrackID();
462 // offset depends on whether the trackis a primary or secondary
463 int offset = (oldID < nprim) ? idelta0 : idelta1;
464 hit.SetTrackID(oldID + offset);
465 }
466 // this could be further generalized by using a policy for T
467 std::copy(incomingdata->begin(), incomingdata->end(), std::back_inserter(*targetdata));
468 }
469 // adjust offsets for next subevent
470 idelta0 += nprim;
471 idelta1 += trackoffsets[index];
472 } // subevent loop
473 }
474 // fill target for this event
475 auto targetbr = o2::base::getOrMakeBranch(target, brname.c_str(), &filladdress);
476 targetbr->SetAddress(&filladdress);
477 targetbr->Fill();
478 targetbr->ResetAddress();
479 targetdata->clear();
480 hitbuffervector.clear();
481 hitbuffervector = L(); // swap with empty vector to release mem
482 delete targetdata;
483 }
484
485 void mergeHitEntries(TTree& origin, TTree& target, std::vector<int> const& trackoffsets, std::vector<int> const& nprimaries, std::vector<int> const& subevtsOrdered) final
486 {
487 // loop over hit containers / different branches
488 // adjust trackID in hits on the go
489 int probe = 0;
490 using Hit_t = decltype(static_cast<Det*>(this)->Det::getHits(probe));
491 std::string name = static_cast<Det*>(this)->getHitBranchNames(probe++);
492 while (name.size() > 0) {
493 mergeAndAdjustHits<typename std::remove_pointer<Hit_t>::type>(name, origin, target, trackoffsets, nprimaries, subevtsOrdered);
494 // next name
495 name = static_cast<Det*>(this)->getHitBranchNames(probe++);
496 }
497 }
498
499 void mergeHitEntriesAndFlush(int eventID, TTree& target, std::vector<int> const& trackoffsets, std::vector<int> const& nprimaries, std::vector<int> const& subevtsOrdered) final
500 {
501 // loop over hit containers / different branches
502 // adjust trackID in hits on the go
503 int probe = 0;
504 using Hit_t = typename std::remove_pointer<decltype(static_cast<Det*>(this)->Det::getHits(0))>::type;
505 // remove buffered event from the hit store
506 using Collector_t = tbb::concurrent_unordered_map<int, std::vector<std::vector<std::unique_ptr<Hit_t>>>>;
507 auto hitbufferPtr = reinterpret_cast<Collector_t*>(mHitCollectorBufferPtr);
508 auto iter = hitbufferPtr->find(eventID);
509 if (iter == hitbufferPtr->end()) {
510 LOG(error) << "No buffered hits available for event " << eventID;
511 return;
512 }
513
514 std::string name = static_cast<Det*>(this)->getHitBranchNames(probe);
515 while (name.size() > 0) {
516 auto& vectorofHitBuffers = (*iter).second[probe];
517 // flushing and buffer removal is done inside here:
518 mergeAndAdjustHits<Hit_t>(name, vectorofHitBuffers, target, trackoffsets, nprimaries, subevtsOrdered);
519 // next name
520 probe++;
521 name = static_cast<Det*>(this)->getHitBranchNames(probe);
522 }
523 }
524
525 public:
529 void collectHits(int eventID, fair::mq::Parts& parts, int& index) override
530 {
531 using Hit_t = typename std::remove_pointer<decltype(static_cast<Det*>(this)->Det::getHits(0))>::type;
532 using Collector_t = tbb::concurrent_unordered_map<int, std::vector<std::vector<std::unique_ptr<Hit_t>>>>;
533 static Collector_t hitcollector; // note: we can't put this as member because
534 // decltype type deduction doesn't seem to work for class members; so we use a static member
535 // and will use some pointer member to communicate this data to other functions
536 mHitCollectorBufferPtr = (char*)&hitcollector;
537
538 int probe = 0;
539 bool* busy = nullptr;
540 using HitPtr_t = decltype(static_cast<Det*>(this)->Det::getHits(probe));
541 std::string name = static_cast<Det*>(this)->getHitBranchNames(probe);
542
543 auto copyToBuffer = [this, eventID](HitPtr_t hitdata, Collector_t& collectbuffer, int probe) {
544 std::vector<std::vector<std::unique_ptr<Hit_t>>>* hitvector = nullptr;
545 {
546 auto eventIter = collectbuffer.find(eventID);
547 if (eventIter == collectbuffer.end()) {
548 // key insertion and traversal are thread-safe with tbb so no need
549 // to protect
550 collectbuffer[eventID] = std::vector<std::vector<std::unique_ptr<Hit_t>>>();
551 }
552 hitvector = &(collectbuffer[eventID]);
553 }
554 if (probe >= hitvector->size()) {
555 hitvector->resize(probe + 1);
556 }
557 // add empty hit bucket to list for this event and probe
558 (*hitvector)[probe].emplace_back(new Hit_t());
559 // copy the data into this bucket
560 *((*hitvector)[probe].back()) = *hitdata;
561 };
562
563 while (name.size() > 0) {
565 // for each branch name we extract/decode hits from the message parts ...
566 auto hitsptr = decodeTMessage<HitPtr_t>(parts, index++);
567 if (hitsptr) {
568 // ... and copy them to the buffer
569 copyToBuffer(hitsptr, hitcollector, probe);
570 delete hitsptr;
571 }
572 } else {
573 // for each branch name we extract/decode hits from the message parts ...
574 auto hitsptr = decodeShmMessage<HitPtr_t>(parts, index++, busy);
575 // ... and copy them to the buffer
576 copyToBuffer(hitsptr, hitcollector, probe);
577 }
578 // next name
579 probe++;
580 name = static_cast<Det*>(this)->getHitBranchNames(probe);
581 }
582 // there is only one busy flag per detector so we need to clear it only
583 // at the end (after all branches have been treated)
584 if (busy) {
585 *busy = false;
586 }
587 }
588
589 void fillHitBranch(TTree& tr, fair::mq::Parts& parts, int& index) override
590 {
591 int probe = 0;
592 bool* busy = nullptr;
593 using Hit_t = decltype(static_cast<Det*>(this)->Det::getHits(probe));
594 std::string name = static_cast<Det*>(this)->getHitBranchNames(probe++);
595 while (name.size() > 0) {
597
598 // for each branch name we extract/decode hits from the message parts ...
599 auto hitsptr = decodeTMessage<Hit_t>(parts, index++);
600 if (hitsptr) {
601 // ... and fill the tree branch
602 auto br = getOrMakeBranch(tr, name.c_str(), hitsptr);
603 br->SetAddress(static_cast<void*>(&hitsptr));
604 br->Fill();
605 br->ResetAddress();
606 delete hitsptr;
607 }
608 } else {
609 // for each branch name we extract/decode hits from the message parts ...
610 auto hitsptr = decodeShmMessage<Hit_t>(parts, index++, busy);
611 // ... and fill the tree branch
612 auto br = getOrMakeBranch(tr, name.c_str(), hitsptr);
613 br->SetAddress(static_cast<void*>(&hitsptr));
614 br->Fill();
615 br->ResetAddress();
616 }
617 // next name
618 name = static_cast<Det*>(this)->getHitBranchNames(probe++);
619 }
620 // there is only one busy flag per detector so we need to clear it only
621 // at the end (after all branches have been treated)
622 if (busy) {
623 *busy = false;
624 }
625 }
626
627 // implementing CloneModule (for G4-MT mode) automatically for each deriving
628 // Detector class "Det"; calls copy constructor of Det
629 FairModule* CloneModule() const final
630 {
631 return new Det(static_cast<const Det&>(*this));
632 }
633
635 {
636 using Hit_t = decltype(static_cast<Det*>(this)->Det::getHits(0));
637 if (UseShm<Det>::value) {
638 for (int buffer = 0; buffer < NHITBUFFERS; ++buffer) {
639 for (auto ptr : mCachedPtr[buffer]) {
640 o2::utils::freeSimVector(static_cast<Hit_t>(ptr));
641 }
642 }
643 }
644 }
645
646 // default implementation for setting hits
647 // always returns false indicating that there is no other
648 // component to assign to apart from i == 0
649 template <typename Hit_t>
650 bool setHits(int i, std::vector<Hit_t>* ptr)
651 {
652 if (i == 0) {
653 static_cast<Det*>(this)->Det::mHits = ptr;
654 }
655 return false;
656 }
657
658 // creating a number of hit buffers (in shared mem) -- to which
659 // detectors can write in round-robin fashion
661 {
662 using VectorHit_t = decltype(static_cast<Det*>(this)->Det::getHits(0));
663 using Hit_t = typename std::remove_pointer<VectorHit_t>::type::value_type;
664 for (int buffer = 0; buffer < NHITBUFFERS; ++buffer) {
665 int probe = 0;
666 bool more{false};
667 do {
668 auto ptr = o2::utils::createSimVector<Hit_t>();
669 more = static_cast<Det*>(this)->Det::setHits(probe, ptr);
670 mCachedPtr[buffer].emplace_back(ptr);
671 probe++;
672 } while (more);
673 }
674 }
675
676 void initializeLate() final
677 {
678 if (!mInitialized) {
679 if (UseShm<Det>::value) {
680 static_cast<Det*>(this)->Det::createHitBuffers();
681 for (int b = 0; b < NHITBUFFERS; ++b) {
682 auto& instance = o2::utils::ShmManager::Instance();
683 mShmBusy[b] = instance.hasSegment() ? (bool*)instance.getmemblock(sizeof(bool)) : new bool;
684 *mShmBusy[b] = false;
685 }
686 }
687 mInitialized = true;
688 mCurrentBuffer = 0;
689 }
690 }
691
692 void BeginEvent() final
693 {
694 if (UseShm<Det>::value) {
696 while (mShmBusy[mCurrentBuffer] != nullptr && *mShmBusy[mCurrentBuffer]) {
697 // this should ideally never happen
698 LOG(info) << " BUSY WAITING SIZE ";
699 sleep(1);
700 }
701
702 using Hit_t = decltype(static_cast<Det*>(this)->Det::getHits(0));
703
704 // now we have to clear the hits before writing again
705 int probe = 0;
706 for (auto bareptr : mCachedPtr[mCurrentBuffer]) {
707 auto hits = static_cast<Hit_t>(bareptr);
708 // assign ..
709 static_cast<Det*>(this)->Det::setHits(probe, hits);
710 hits->clear();
711 probe++;
712 }
713 }
714 }
715
716 ~DetImpl() override
717 {
718 for (int i = 0; i < NHITBUFFERS; ++i) {
719 if (mShmBusy[i]) {
720 auto& instance = o2::utils::ShmManager::Instance();
721 if (instance.hasSegment()) {
722 instance.freememblock(mShmBusy[i]);
723 } else {
724 delete mShmBusy[i];
725 }
726 }
727 }
729 }
730
731 protected:
732 static constexpr int NHITBUFFERS = 3; // number of buffers for hits in order to allow async processing
733 // in the hit merger without blocking nor copying the data
734 // (like done in typical data aquisition systems)
735 bool* mShmBusy[NHITBUFFERS] = {nullptr};
736 std::vector<void*> mCachedPtr[NHITBUFFERS];
737 int mCurrentBuffer = 0; // holding the current buffer information
738 int mInitialized = false;
739
740 char* mHitCollectorBufferPtr = nullptr;
741
743};
744} // namespace o2::base
745
746#endif
int32_t i
uint32_t res
Definition RawData.h:0
TBranch * ptr
void freeHitBuffers()
Definition Detector.h:634
~DetImpl() override
Definition Detector.h:716
void updateHitTrackIndices(std::map< int, int > const &indexmapping) override
Definition Detector.h:320
void mergeHitEntriesAndFlush(int eventID, TTree &target, std::vector< int > const &trackoffsets, std::vector< int > const &nprimaries, std::vector< int > const &subevtsOrdered) final
Definition Detector.h:499
std::vector< void * > mCachedPtr[NHITBUFFERS]
pointer to bool in shared mem indicating of IO busy
Definition Detector.h:736
bool * mShmBusy[NHITBUFFERS]
Definition Detector.h:735
std::string getHitBranchNames(int probe) const override
Definition Detector.h:309
static constexpr int NHITBUFFERS
Definition Detector.h:732
void fillHitBranch(TTree &tr, fair::mq::Parts &parts, int &index) override
Definition Detector.h:589
void createHitBuffers()
Definition Detector.h:660
void BeginEvent() final
Definition Detector.h:692
void attachHits(fair::mq::Channel &channel, fair::mq::Parts &parts) override
Definition Detector.h:332
void mergeAndAdjustHits(std::string const &brname, TTree &origin, TTree &target, std::vector< int > const &trackoffsets, std::vector< int > const &nprimaries, std::vector< int > const &subevtsOrdered)
Definition Detector.h:360
void mergeAndAdjustHits(std::string const &brname, L &hitbuffervector, TTree &target, std::vector< int > const &trackoffsets, std::vector< int > const &nprimaries, std::vector< int > const &subevtsOrdered)
Definition Detector.h:427
char * mHitCollectorBufferPtr
Definition Detector.h:740
bool setHits(int i, std::vector< Hit_t > *ptr)
Definition Detector.h:650
FairModule * CloneModule() const final
Definition Detector.h:629
void initializeLate() final
Definition Detector.h:676
void mergeHitEntries(TTree &origin, TTree &target, std::vector< int > const &trackoffsets, std::vector< int > const &nprimaries, std::vector< int > const &subevtsOrdered) final
Definition Detector.h:485
ClassDefOverride(DetImpl, 0)
pointer to hit (collector) buffer location (strictly internal)
void collectHits(int eventID, fair::mq::Parts &parts, int &index) override
Definition Detector.h:529
virtual std::string getHitBranchNames(int probe) const =0
void SpecialProcess(Int_t numed, EProc parID, int val)
Set process by name and value.
Definition Detector.cxx:98
int registerSensitiveVolumeAndGetVolID(std::string const &name)
Definition Detector.cxx:190
~Detector() override
Default Destructor.
Definition Detector.cxx:87
static std::vector< int > const & getDetId2HitBitIndex()
Definition Detector.h:221
virtual void initializeLate()=0
Detector & operator=(const Detector &)
Definition Detector.cxx:46
virtual void defineLayerTurbo(Int_t nlay, Double_t phi0, Double_t r, Int_t nladd, Int_t nmod, Double_t width, Double_t tilt, Double_t lthick=0., Double_t dthick=0., UInt_t detType=0, Int_t buildFlag=0)
Definition Detector.cxx:117
TClonesArray * GetCollection(int iColl) const final
Definition Detector.cxx:162
virtual void attachHits(fair::mq::Channel &, fair::mq::Parts &)=0
virtual void fillParallelWorld() const
fill parallel geometry with sensitive volumes of detector
Definition Detector.cxx:174
void Initialize() final
Definition Detector.h:190
int getMaterialID(int imat) const
Definition Detector.h:127
void getMediumIDMappingAsVector(std::vector< int > &mapping)
Definition Detector.h:143
static void setDensityFactor(Float_t density)
Definition Detector.h:92
void SpecialCuts(Int_t numed, const std::initializer_list< std::pair< ECut, Float_t > > &parIDValMap)
Custom processes and transport cuts.
Definition Detector.cxx:80
void SpecialCut(Int_t numed, ECut parID, Float_t val)
Set cut by name and value.
Definition Detector.cxx:86
void SpecialProcesses(Int_t numed, const std::initializer_list< std::pair< EProc, int > > &parIDValMap)
Definition Detector.cxx:92
void Matrix(Int_t &nmat, Float_t theta1, Float_t phi1, Float_t theta2, Float_t phi2, Float_t theta3, Float_t phi3) const
Definition Detector.cxx:104
Detector()
Default Constructor.
Definition Detector.cxx:36
virtual void setNumberOfWrapperVolumes(Int_t n)
Books arrays for wrapper volumes.
Definition Detector.cxx:111
virtual void updateHitTrackIndices(std::map< int, int > const &)=0
virtual void mergeHitEntries(TTree &origin, TTree &target, std::vector< int > const &trackoffsets, std::vector< int > const &nprimaries, std::vector< int > const &subevtsOrdered)=0
void Mixture(Int_t imat, const char *name, Float_t *a, Float_t *z, Float_t dens, Int_t nlmat, Float_t *wmat)
Definition Detector.cxx:66
virtual void addAlignableVolumes() const
declare alignable volumes of detector
Definition Detector.cxx:169
void Medium(Int_t numed, const char *name, Int_t nmat, Int_t isvol, Int_t ifield, Float_t fieldm, Float_t tmaxfd, Float_t stemax, Float_t deemax, Float_t epsil, Float_t stmin, Float_t *ubuf=nullptr, Int_t nbuf=0)
Definition Detector.cxx:72
void SetSpecialPhysicsCuts() override
Definition Detector.cxx:122
static Float_t getDensityFactor()
Definition Detector.h:97
int getMediumID(int imed) const
Definition Detector.h:135
virtual void fillHitBranch(TTree &tr, fair::mq::Parts &parts, int &index)=0
static void setDetId2HitBitIndex(std::vector< int > const &v)
set the DetID to HitBitIndex mapping. Succeeds if not already set.
Definition Detector.h:220
static void initFieldTrackingParams(int &mode, float &maxfield)
Definition Detector.cxx:143
virtual void defineWrapperVolume(Int_t id, Double_t rmin, Double_t rmax, Double_t zspan)
Sets per wrapper volume parameters.
Definition Detector.cxx:110
virtual void defineLayer(Int_t nlay, Double_t phi0, Double_t r, Int_t nladd, Int_t nmod, Double_t lthick=0., Double_t dthick=0., UInt_t detType=0, Int_t buildFlag=0)
virtual void mergeHitEntriesAndFlush(int eventID, TTree &target, std::vector< int > const &trackoffsets, std::vector< int > const &nprimaries, std::vector< int > const &subevtsOrdered)=0
void Material(Int_t imat, const char *name, Float_t a, Float_t z, Float_t dens, Float_t radl, Float_t absl, Float_t *buf=nullptr, Int_t nwbuf=0)
Definition Detector.cxx:59
std::string addNameTo(const char *ext) const
Definition Detector.h:150
virtual void InitializeO2Detector()=0
Definition Detector.cxx:98
virtual void collectHits(int eventID, fair::mq::Parts &parts, int &index)=0
static MaterialManager & Instance()
static ShmManager & Instance()
Definition ShmManager.h:61
bool isOperational() const
Definition ShmManager.h:97
GLdouble n
Definition glcorearb.h:1982
GLenum mode
Definition glcorearb.h:266
GLuint buffer
Definition glcorearb.h:655
GLuint entry
Definition glcorearb.h:5735
const GLdouble * v
Definition glcorearb.h:832
GLuint index
Definition glcorearb.h:781
GLuint const GLchar * name
Definition glcorearb.h:781
GLint GLsizei width
Definition glcorearb.h:270
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLsizei const GLfloat * value
Definition glcorearb.h:819
GLenum target
Definition glcorearb.h:1641
GLint GLint GLsizei GLint GLenum GLenum type
Definition glcorearb.h:275
GLboolean * data
Definition glcorearb.h:298
GLintptr offset
Definition glcorearb.h:660
typedef void(APIENTRYP PFNGLCULLFACEPROC)(GLenum mode)
GLuint GLfloat * val
Definition glcorearb.h:1582
GLboolean r
Definition glcorearb.h:1233
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
GLenum GLuint GLenum GLsizei const GLchar * buf
Definition glcorearb.h:2514
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
void attachShmMessage(void *hitsptr, fair::mq::Channel &channel, fair::mq::Parts &parts, bool *busy_ptr)
Definition Detector.cxx:223
ECut
cuts available
EProc
processes available
void * decodeTMessageCore(fair::mq::Parts &dataparts, int index)
Definition Detector.cxx:255
T decodeShmMessage(fair::mq::Parts &dataparts, int index, bool *&busy)
Definition Detector.h:257
TBranch * getOrMakeBranch(TTree &tree, const char *brname, T *ptr)
Definition Detector.h:281
void attachMessageBufferToParts(fair::mq::Parts &parts, fair::mq::Channel &channel, void *data, TClass *cl)
Definition Detector.cxx:207
void attachTMessage(Container const &hits, fair::mq::Channel &channel, fair::mq::Parts &parts)
Definition Detector.h:266
void * decodeShmCore(fair::mq::Parts &dataparts, int index, bool *&busy)
Definition Detector.cxx:240
std::string demangle(const char *name)
utility function to demangle cxx type names
Definition Detector.h:246
void attachDetIDHeaderMessage(int id, fair::mq::Channel &channel, fair::mq::Parts &parts)
Definition Detector.cxx:218
T decodeTMessage(fair::mq::Parts &dataparts, int index)
Definition Detector.h:273
void freeSimVector(std::vector< T > *ptr)
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
std::unique_ptr< TTree > tree((TTree *) flIn.Get(std::string(o2::base::NameConf::CTFTREENAME).c_str()))