9 #include "cdc/modules/cdcCalibrationCollector/CDCT0CalibrationCollector.h"
11 #include <framework/datastore/RelationArray.h>
13 #include <genfit/TrackPoint.h>
14 #include <genfit/KalmanFitterInfo.h>
15 #include <genfit/MeasuredStateOnPlane.h>
16 #include <Math/ProbFuncMathCore.h>
18 #include <cdc/dataobjects/WireID.h>
19 #include <cdc/geometry/CDCGeometryPar.h>
30 setDescription(
"Collector module for cdc T0 calibration");
31 setPropertyFlags(c_ParallelProcessingCertified);
32 addParam(
"RecoTracksColName", m_recoTrackArrayName,
"Name of collection hold genfit::Track", std::string(
""));
33 addParam(
"BField", m_BField,
"If true -> #Params ==5 else #params ==4 for calculate P-Val",
false);
34 addParam(
"EventT0Extraction", m_EventT0Extraction,
"use event t0 extract t0 or not",
false);
35 addParam(
"MinimumPt", m_MinimumPt,
"Tracks whose Pt lower than this value will not be recoreded", 0.);
36 addParam(
"PvalCut", m_PvalCut,
"Tracks which Pval small than this will not be recoreded", 0.);
37 addParam(
"NDFCut", m_ndfCut,
"Tracks which NDF small than this will not be recoreded", 0.);
38 addParam(
"Xmin", m_xmin,
"minimum dift length", 0.1);
41 CDCT0CalibrationCollectorModule::~CDCT0CalibrationCollectorModule()
45 void CDCT0CalibrationCollectorModule::prepare()
47 m_Tracks.isRequired(m_trackArrayName);
48 m_RecoTracks.isRequired(m_recoTrackArrayName);
49 m_TrackFitResults.isRequired(m_trackFitResultArrayName);
50 m_CDCHits.isRequired(m_cdcHitArrayName);
51 RelationArray relRecoTrackTrack(m_RecoTracks, m_Tracks, m_relRecoTrackTrackName);
53 m_relRecoTrackTrackName = relRecoTrackTrack.
getName();
55 auto m_hNDF =
new TH1D(
"hNDF",
"NDF of fitted track;NDF;Tracks", 71, -1, 70);
56 auto m_hPval =
new TH1D(
"hPval",
"p-values of tracks;pVal;Tracks", 1000, 0, 1);
58 registerObject<TH1D>(
"hNDF", m_hNDF);
59 registerObject<TH1D>(
"hPval", m_hPval);
62 for (
int i = 0; i < 56; ++i) {
65 halfCSize[i] = M_PI * R / nW;
68 auto m_hTotal =
new TH1F(
"hTotal",
"hTotal", 30, -15, 15);
69 registerObject<TH1F>(
"hTotal", m_hTotal);
73 for (
int il = 0; il < 56; il++) {
75 for (
int ic = 0; ic < nW; ++ic) {
76 m_h1[il][ic] =
new TH1F(Form(
"hdT_lay%d_cell%d", il, ic), Form(
"Lay%d_cell%d", il, ic), 30, -15, 15);
77 registerObject<TH1F>(Form(
"hdT_lay%d_cell%d", il, ic), m_h1[il][ic]);
81 for (
int ib = 0; ib < 300; ++ib) {
82 m_hT0b[ib] =
new TH1F(Form(
"hT0b%d", ib), Form(
"boardID_%d", ib), 100, -20, 20);
83 registerObject<TH1F>(Form(
"hT0b%d", ib), m_hT0b[ib]);
88 void CDCT0CalibrationCollectorModule::collect()
90 const RelationArray relTrackTrack(m_RecoTracks, m_Tracks, m_relRecoTrackTrackName);
95 const int nTr = m_RecoTracks.getEntries();
97 for (
int i = 0; i < nTr; ++i) {
99 if (track->getDirtyFlag())
continue;
101 if (!fs || !fs->
isFitted())
continue;
106 B2DEBUG(99,
"No relation found");
113 B2WARNING(
"track was fitted but Relation not found");
125 Pval = std::max(0., ROOT::Math::chisquared_cdf_c(Chi2, ndf));
128 getObjectPtr<TH1D>(
"hPval")->Fill(Pval);
129 getObjectPtr<TH1D>(
"hNDF")->Fill(ndf);
131 if (Pval < m_PvalCut || ndf < m_ndfCut)
continue;
133 if (fitresult->
getMomentum().Perp() < m_MinimumPt)
continue;
135 if (m_EventT0Extraction) {
137 if (m_eventTimeStoreObject.isValid() && m_eventTimeStoreObject->hasEventT0()) {
138 evtT0 = m_eventTimeStoreObject->getEventT0();
144 B2DEBUG(99,
"start collect hit");
148 const genfit::TrackPoint* tp = track->getCreatedTrackPoint(track->getRecoHitInformation(hit));
149 int lay = hit->getICLayer();
150 int IWire = hit->getIWire();
151 unsigned short tdc = hit->getTDCCount();
152 unsigned short adc = hit->getADCCount();
153 WireID wireid(lay, IWire);
155 if (!kfi) {B2DEBUG(199,
"No Fitter Info: Layer " << hit->getICLayer());
continue;}
156 for (
unsigned int iMeas = 0; iMeas < kfi->getNumMeasurements(); ++iMeas) {
160 const TVector3 pocaOnWire = mop.getPlane()->getO();
162 const TVector3 pocaMom = mop.getMom();
163 double alpha = cdcgeo.
getAlpha(pocaOnWire, pocaMom) ;
164 double theta = cdcgeo.
getTheta(pocaMom);
168 double xmax = halfCSize[lay] - 0.1;
169 if ((fabs(x_u) < m_xmin) || (fabs(x_u) > xmax))
continue;
176 if (fabs(alpha) > M_PI / 2) {
184 double t_fit = cdcgeo.
getDriftTime(std::abs(x_u), lay, lr, alpha , theta);
189 dt_flight = mop.getTime();
190 if (dt_flight < 50) {
196 double z = pocaOnWire.Z();
198 double z_prop = z - m_backWirePos.Z();
199 B2DEBUG(99,
"z_prop = " << z_prop <<
" |z " << z <<
" |back wire poss: " << m_backWirePos.Z());
211 getObjectPtr<TH1F>(Form(
"hdT_lay%d_cell%d", lay, IWire))->Fill(t - t_fit);
212 getObjectPtr<TH1F>(Form(
"hT0b%d", boardID))->Fill(t - t_fit);
213 getObjectPtr<TH1F>(
"hTotal")->Fill(t - t_fit);
221 void CDCT0CalibrationCollectorModule::finish()
Class containing the result of the unpacker in raw data and the result of the digitizer in simulation...
The Class for CDC Geometry Parameters.
double getTimeWalk(const WireID &wID, unsigned short adcCount) const
Returns time-walk.
unsigned short getBoardID(const WireID &wID) const
Returns frontend board id. corresponding to the wire id.
const TVector3 wireBackwardPosition(int layerId, int cellId, EWirePosition set=c_Base) const
Returns the backward position of the input sense wire.
double getAlpha(const TVector3 &posOnWire, const TVector3 &momentum) const
Returns track incident angle in rphi plane (alpha in rad.).
double getTdcBinWidth() const
Return TDC bin width (nsec).
double getOutgoingAlpha(const double alpha) const
Converts incoming- to outgoing-alpha.
float getT0(const WireID &wireID) const
Returns t0 parameter of the specified sense wire.
unsigned short getOutgoingLR(const unsigned short lr, const double alpha) const
Converts incoming-lr to outgoing-lr.
double getOutgoingTheta(const double alpha, const double theta) const
Converts incoming- to outgoing-theta.
unsigned nWiresInLayer(int layerId) const
Returns wire numbers in a layer.
double getTheta(const TVector3 &momentum) const
Returns track incident angle (theta in rad.).
double getPropSpeedInv(const unsigned int layerID) const
Get the inversel of propagation speed in the sense wire.
double getDriftTime(double dist, unsigned short layer, unsigned short lr, double alpha, double theta) const
Return the drift time to the sense wire.
double senseWireR(int layerId) const
Returns radius of sense wire in each layer.
Collect hit information for cdc calibration with CAF.
Calibration collector module base class.
This is the Reconstruction Event-Data Model Track.
Low-level class to create/modify relations between StoreArrays.
const std::string & getName() const
Return name under which the object is saved in the DataStore.
Values of the result of a track fit with a given particle hypothesis.
TVector3 getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
Class that bundles various TrackFitResults.
const TrackFitResult * getTrackFitResult(const Const::ChargedStable &chargedStable) const
Access to TrackFitResults.
Class to identify a wire inside the CDC.
Class where important numbers and properties of a fit can be stored.
bool isFitConverged(bool inAllPoints=true) const
Did the fit converge (in all Points or only partially)?
virtual double getPVal() const
Get the p value of the fit.
double getChi2() const
Get chi^2 of the fit.
bool isFitted() const
Has the track been fitted?
double getNdf() const
Get the degrees of freedom of the fit.
Collects information needed and produced by a AbsKalmanFitter implementations and is specific to one ...
const MeasuredStateOnPlane & getFittedState(bool biased=true) const override
Get unbiased or biased (default) smoothed state.
std::vector< double > getWeights() const
Get weights of measurements.
#StateOnPlane with additional covariance matrix.
Object containing AbsMeasurement and AbsFitterInfo objects.
KalmanFitterInfo * getKalmanFitterInfo(const AbsTrackRep *rep=nullptr) const
Helper to avoid casting.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.
Defines for I/O streams used for error and debug printing.