Project
Loading...
Searching...
No Matches
TRDDigitWriterSpec.cxx
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
12#ifndef STEER_DIGITIZERWORKFLOW_SRC_TRDDIGITWRITERSPEC_H_
13#define STEER_DIGITIZERWORKFLOW_SRC_TRDDIGITWRITERSPEC_H_
14
17#include "Framework/InputSpec.h"
23
24namespace o2
25{
26namespace trd
27{
28
29template <typename T>
30using BranchDefinition = framework::MakeRootTreeWriterSpec::BranchDefinition<T>;
31
32o2::framework::DataProcessorSpec getTRDDigitWriterSpec(bool mctruth, bool inpFromDigitizer)
33{
37
38 // the callback to be set as hook for custom action when the writer is closed
39 auto finishWriting = [](TFile* outputfile, TTree* outputtree) {
40 const auto* brArr = outputtree->GetListOfBranches();
41 int64_t nent = 0;
42 for (const auto* brc : *brArr) {
43 int64_t n = ((const TBranch*)brc)->GetEntries();
44 if (nent && (nent != n)) {
45 LOG(error) << "Branches have different number of entries";
46 }
47 nent = n;
48 }
49 outputtree->SetEntries(nent);
50 outputtree->Write();
51 outputfile->Close();
52 };
53
54 // custom handler for labels:
55 // essentially transform the input container (as registered in the original branch definition) to the special output format for labels
56 auto customlabelhandler = [](TBranch& branch, std::vector<char> const& labeldata, DataRef const& ref) {
57 // make the actual output object by adopting/casting the buffer
58 // into a split format
59 o2::dataformats::IOMCTruthContainerView outputcontainer(labeldata);
60 auto ptr = &outputcontainer;
62 br->Fill();
63 br->ResetAddress();
64 };
65
66 auto labelsdef = BranchDefinition<std::vector<char>>{InputSpec{"labelinput", "TRD", "LABELS"},
67 "TRDMCLabels", "labels-branch-name",
68 // this branch definition is disabled if MC labels are not processed
69 (mctruth ? 1 : 0),
70 customlabelhandler};
71
72 return MakeRootTreeWriterSpec("TRDDigitWriter",
73 "trddigits.root",
74 "o2sim",
75 // setting a custom callback for closing the writer
77 BranchDefinition<std::vector<o2::trd::Digit>>{InputSpec{"input", "TRD", "DIGITS", (inpFromDigitizer ? 1u : 0u)}, "TRDDigit"},
78 BranchDefinition<std::vector<o2::trd::TriggerRecord>>{InputSpec{"trinput", "TRD", "TRKTRGRD", (inpFromDigitizer ? 1u : 0u)}, "TriggerRecord"},
79 std::move(labelsdef))();
80}
81
82} // end namespace trd
83} // end namespace o2
84
85#endif /* STEER_DIGITIZERWORKFLOW_SRC_TRDDIGITWRITERSPEC_H_ */
A const (ready only) version of MCTruthContainer.
A special IO container - splitting a given vector to enable ROOT IO.
Configurable generator for RootTreeWriter processor spec.
TBranch * ptr
Generate a processor spec for the RootTreeWriter utility.
static TBranch * remapBranch(TBranch &branchRef, T **newdata)
GLdouble n
Definition glcorearb.h:1982
GLint ref
Definition glcorearb.h:291
o2::framework::DataProcessorSpec getTRDDigitWriterSpec(bool mctruth=true, bool inpFromDigitizer=true)
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"