Project
Loading...
Searching...
No Matches
TPCScaler.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_TPC_TPCSCALER
16#define ALICEO2_TPC_TPCSCALER
17
18#include "DataFormatsTPC/Defs.h"
19#include <vector>
20
21class TTree;
22
23namespace o2::tpc
24{
25
26/*
27Class for storing the scalers which are used to calculate an estimate for the mean space-charge density for the last ion drift time
28*/
29
31 float getWeight(float deltaTime) const;
32 bool isValid() const { return (mSamplingTimeMS > 0 && !mWeights.empty()); }
33 float getDurationMS() const { return mSamplingTimeMS * mWeights.size(); }
34 float mSamplingTimeMS = -1;
36 std::vector<float> mWeights;
37
39};
40
42{
43 public:
45 TPCScaler() = default;
46
49
51 int getNValues(o2::tpc::Side side) const { return (side == o2::tpc::Side::A ? mScalerA.size() : mScalerC.size()); }
52
55 void setScaler(const std::vector<float>& values, const o2::tpc::Side side) { (side == o2::tpc::Side::A ? (mScalerA = values) : (mScalerC = values)); }
56
58 void setIonDriftTimeMS(float ionDriftTimeMS) { mIonDriftTimeMS = ionDriftTimeMS; }
59
61 void setRun(int run) { mRun = run; }
62
64 void setFirstTFOrbit(unsigned int firstTFOrbit) { mFirstTFOrbit = firstTFOrbit; }
65
67 void setStartTimeStampMS(double timeStampMS) { mTimeStampMS = timeStampMS; }
68
70 void setIntegrationTimeMS(float integrationTimeMS) { mIntegrationTimeMS = integrationTimeMS; }
71
74 void dumpToFile(const char* file, const char* name);
75
80 void dumpToFile(const char* file, const char* name, double startTimeMS, double endTimeMS);
81
87 void dumpToFileSlices(const char* file, const char* name, int minutesPerObject, float marginMS = 500, float marginCCDBMinutes = 1);
88
91 void loadFromFile(const char* inpf, const char* name);
92
94 void setFromTree(TTree& tpcScalerTree);
95
97 float getScalers(unsigned int idx, o2::tpc::Side side) const { return (side == o2::tpc::Side::A) ? mScalerA[idx] : mScalerC[idx]; }
98
100 const auto& getScalers(o2::tpc::Side side) const { return (side == o2::tpc::Side::A) ? mScalerA : mScalerC; }
101
103 float getIonDriftTimeMS() const { return mIonDriftTimeMS; }
104
106 int getRun() const { return mRun; }
107
109 unsigned int getFirstTFOrbit() const { return mFirstTFOrbit; }
110
112 double getStartTimeStampMS() const { return mTimeStampMS; }
113
115 double getEndTimeStampMS(o2::tpc::Side side) const { return mTimeStampMS + getDurationMS(side); }
116
118 float getIntegrationTimeMS() const { return mIntegrationTimeMS; }
119
121 int getNValuesIonDriftTime() const { return mIonDriftTimeMS / mIntegrationTimeMS + 0.5; }
122
125 float getMeanScaler(double timestamp, o2::tpc::Side side) const;
126
128 float getDurationMS(o2::tpc::Side side) const { return mIntegrationTimeMS * getNValues(side); }
129
133 void clampScalers(float minThreshold, float maxThreshold, Side side);
134
136 void setScalerWeights(const TPCScalerWeights& weights) { mScalerWeights = weights; }
137
139 const auto& getScalerWeights() const { return mScalerWeights; }
140
142 void useWeights(bool useWeights) { mUseWeights = useWeights; }
143
145 bool weightsUsed() const { return mUseWeights; }
146
147 private:
148 float mIonDriftTimeMS{200};
149 int mRun{};
150 unsigned int mFirstTFOrbit{};
151 double mTimeStampMS{};
152 float mIntegrationTimeMS{1.};
153 std::vector<float> mScalerA{};
154 std::vector<float> mScalerC{};
155 TPCScalerWeights mScalerWeights{};
156 bool mUseWeights{false};
157
158 // get index to data for given timestamp
159 int getDataIdx(double timestamp) const { return (timestamp - mTimeStampMS) / mIntegrationTimeMS + 0.5; }
160
161 ClassDefNV(TPCScaler, 2);
162};
163
164} // namespace o2::tpc
165#endif
uint32_t side
Definition RawData.h:0
unsigned int getFirstTFOrbit() const
Definition TPCScaler.h:109
bool weightsUsed() const
return if weights are used
Definition TPCScaler.h:145
int getNValuesIonDriftTime() const
Definition TPCScaler.h:121
void dumpToFileSlices(const char *file, const char *name, int minutesPerObject, float marginMS=500, float marginCCDBMinutes=1)
Definition TPCScaler.cxx:53
TPCScaler()=default
default constructor
float getScalers(unsigned int idx, o2::tpc::Side side) const
Definition TPCScaler.h:97
float getMeanScaler(double timestamp, o2::tpc::Side side) const
void setScaler(const std::vector< float > &values, const o2::tpc::Side side)
Definition TPCScaler.h:55
void setIonDriftTimeMS(float ionDriftTimeMS)
Definition TPCScaler.h:58
TPCScaler & operator=(TPCScaler &&other)=default
default move assignment
int getNValues(o2::tpc::Side side) const
Definition TPCScaler.h:51
float getIntegrationTimeMS() const
Definition TPCScaler.h:118
float getDurationMS(o2::tpc::Side side) const
Definition TPCScaler.h:128
float getIonDriftTimeMS() const
Definition TPCScaler.h:103
double getStartTimeStampMS() const
Definition TPCScaler.h:112
void loadFromFile(const char *inpf, const char *name)
Definition TPCScaler.cxx:95
const auto & getScalers(o2::tpc::Side side) const
Definition TPCScaler.h:100
int getRun() const
Definition TPCScaler.h:106
void clampScalers(float minThreshold, float maxThreshold, Side side)
double getEndTimeStampMS(o2::tpc::Side side) const
Definition TPCScaler.h:115
void dumpToFile(const char *file, const char *name)
Definition TPCScaler.cxx:25
void setFirstTFOrbit(unsigned int firstTFOrbit)
Definition TPCScaler.h:64
void setStartTimeStampMS(double timeStampMS)
Definition TPCScaler.h:67
void setRun(int run)
Definition TPCScaler.h:61
void useWeights(bool useWeights)
enable usage of weights
Definition TPCScaler.h:142
void setFromTree(TTree &tpcScalerTree)
set this object from input tree
const auto & getScalerWeights() const
Definition TPCScaler.h:139
void setIntegrationTimeMS(float integrationTimeMS)
Definition TPCScaler.h:70
void setScalerWeights(const TPCScalerWeights &weights)
setting the weights for the scalers
Definition TPCScaler.h:136
GLsizei const GLuint const GLfloat * weights
Definition glcorearb.h:5475
GLuint const GLchar * name
Definition glcorearb.h:781
GLenum GLsizei GLsizei GLint * values
Definition glcorearb.h:1576
Global TPC definitions and constants.
Definition SimTraits.h:167
Side
TPC readout sidE.
Definition Defs.h:35
@ A
Definition Defs.h:35
float getDurationMS() const
Definition TPCScaler.h:33
float mSamplingTimeMS
sampling of the stored weights
Definition TPCScaler.h:34
std::vector< float > mWeights
stored weights
Definition TPCScaler.h:36
float mFirstTimeStampMS
first timestamp
Definition TPCScaler.h:35
ClassDefNV(TPCScalerWeights, 1)
float getWeight(float deltaTime) const
VectorOfTObjectPtrs other