forked from AliceO2Group/O2Physics
-
Notifications
You must be signed in to change notification settings - Fork 0
Expand file tree
/
Copy pathft0CorrectedTable.cxx
More file actions
120 lines (109 loc) · 5.05 KB
/
ft0CorrectedTable.cxx
File metadata and controls
120 lines (109 loc) · 5.05 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
// All rights not expressly granted are reserved.
//
// This software is distributed under the terms of the GNU General Public
// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
//
// In applying this license CERN does not waive the privileges and immunities
// granted to it by virtue of its status as an Intergovernmental Organization
// or submit itself to any jurisdiction.
#include "Common/Core/CollisionTypeHelper.h"
#include "Common/DataModel/EventSelection.h"
#include "Common/DataModel/FT0Corrected.h"
#include <CCDB/BasicCCDBManager.h>
#include <CommonConstants/PhysicsConstants.h>
#include <DataFormatsFT0/Digit.h>
#include <DataFormatsParameters/GRPLHCIFData.h>
#include <Framework/AnalysisDataModel.h>
#include <Framework/AnalysisHelpers.h>
#include <Framework/AnalysisTask.h>
#include <Framework/Configurable.h>
#include <Framework/HistogramRegistry.h>
#include <Framework/HistogramSpec.h>
#include <Framework/InitContext.h>
#include <Framework/OutputObjHeader.h>
#include <Framework/runDataProcessing.h>
#include <TRandom3.h>
#include <bitset>
#include <chrono>
#include <cstdint>
#include <string>
using namespace o2;
using namespace o2::framework;
using namespace o2::aod;
struct ft0CorrectedTable {
// Configurables
Configurable<float> resoFT0A{"resoFT0A", 20.f, "FT0A resolution in ps for the MC override"};
Configurable<float> resoFT0C{"resoFT0C", 20.f, "FT0C resolution in ps for the MC override"};
Configurable<bool> addHistograms{"addHistograms", false, "Add QA histograms"};
Configurable<int> cfgCollisionSystem{"collisionSystem", -2, "Collision system: -2 (use cfg values), -1 (autoset), 0 (pp), 1 (PbPb), 2 (XeXe), 3 (pPb)"};
Configurable<std::string> cfgUrl{"ccdb-url", "http://alice-ccdb.cern.ch", "url of the ccdb repository"};
Configurable<std::string> cfgPathGrpLhcIf{"ccdb-path-grplhcif", "GLO/Config/GRPLHCIF", "Path on the CCDB for the GRPLHCIF object"};
Configurable<int64_t> cfgTimestamp{"ccdb-timestamp", -1, "timestamp of the object"};
Service<o2::ccdb::BasicCCDBManager> ccdb;
// Producer
Produces<o2::aod::FT0sCorrected> table;
using BCsWithMatchings = soa::Join<aod::BCs, aod::Run3MatchedToBCSparse>;
using CollisionEvSel = soa::Join<aod::Collisions, aod::EvSels>::iterator;
static constexpr float invLightSpeedCm2NS = 1.f / o2::constants::physics::LightSpeedCm2NS;
HistogramRegistry histos{"Histos", {}, OutputObjHandlingPolicy::AnalysisObject};
void init(o2::framework::InitContext&)
{
ccdb->setURL(cfgUrl);
ccdb->setTimestamp(cfgTimestamp);
ccdb->setCaching(true);
ccdb->setLocalObjectValidityChecking();
// Not later than now objects
ccdb->setCreatedNotAfter(std::chrono::duration_cast<std::chrono::milliseconds>(std::chrono::system_clock::now().time_since_epoch()).count());
if (!addHistograms) {
return;
}
histos.add("t0A", "t0A", kTH1D, {{1000, -1, 1, "t0A (ns)"}});
histos.add("t0C", "t0C", kTH1D, {{1000, -1, 1, "t0C (ns)"}});
histos.add("t0AC", "t0AC", kTH1D, {{1000, -1000, 1000, "t0AC (ns)"}});
histos.add("deltat0AC", "deltat0AC", kTH1D, {{1000, -1, 1, "#Deltat0AC (ns)"}});
histos.add("deltat0ACps", "deltat0ACps", kTH1D, {{1000, -1000, 1000, "#Deltat0AC (ps)"}});
}
void process(soa::Join<aod::Collisions, aod::EvSels> const& collisions,
BCsWithMatchings const&,
aod::FT0s const&)
{
table.reserve(collisions.size());
float t0A = 1e10f;
float t0C = 1e10f;
for (const auto& collision : collisions) {
t0A = 1e10f;
t0C = 1e10f;
const float vertexPV = collision.posZ();
const float vertex_corr = vertexPV * invLightSpeedCm2NS;
constexpr float dummyTime = 30.; // Due to HW limitations time can be only within range (-25,25) ns, dummy time is around 32 ns
if (collision.has_foundFT0()) {
const auto& ft0 = collision.foundFT0();
const std::bitset<8>& triggers = ft0.triggerMask();
const bool ora = triggers[o2::ft0::Triggers::bitA];
const bool orc = triggers[o2::ft0::Triggers::bitC];
LOGF(debug, "triggers OrA %i OrC %i ", ora, orc);
LOGF(debug, " T0A = %f, T0C %f, vertex_corr %f", ft0.timeA(), ft0.timeC(), vertex_corr);
if (ora && ft0.timeA() < dummyTime) {
t0A = ft0.timeA() + vertex_corr;
}
if (orc && ft0.timeC() < dummyTime) {
t0C = ft0.timeC() - vertex_corr;
}
}
LOGF(debug, " T0 collision time T0A = %f, T0C = %f", t0A, t0C);
if (addHistograms) {
histos.fill(HIST("t0A"), t0A);
histos.fill(HIST("t0C"), t0C);
if (t0A < 1e10f && t0C < 1e10f) {
histos.fill(HIST("t0AC"), (t0A + t0C) * 0.5f);
histos.fill(HIST("deltat0AC"), (t0A - t0C) * 0.5f);
histos.fill(HIST("deltat0ACps"), (t0A - t0C) * 500.f);
}
}
table(t0A, t0C);
}
}
};
WorkflowSpec defineDataProcessing(ConfigContext const& cfgc) { return WorkflowSpec{adaptAnalysisTask<ft0CorrectedTable>(cfgc)}; }