Belle II Software development
MillepedeCollectorModule.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#include <alignment/modules/MillepedeCollector/MillepedeCollectorModule.h>
10
11#include <alignment/dataobjects/MilleData.h>
12#include <alignment/GblMultipleScatteringController.h>
13#include <alignment/GlobalDerivatives.h>
14#include <alignment/GlobalLabel.h>
15#include <alignment/GlobalParam.h>
16#include <alignment/GlobalTimeLine.h>
17#include <alignment/Manager.h>
18#include <alignment/reconstruction/AlignableCDCRecoHit.h>
19#include <alignment/reconstruction/AlignablePXDRecoHit.h>
20#include <alignment/reconstruction/AlignableSVDRecoHit.h>
21#include <alignment/reconstruction/AlignableSVDRecoHit2D.h>
22#include <alignment/reconstruction/AlignableBKLMRecoHit.h>
23#include <alignment/reconstruction/AlignableEKLMRecoHit.h>
24#include <analysis/dataobjects/ParticleList.h>
25#include <analysis/utility/ReferenceFrame.h>
26#include <framework/core/FileCatalog.h>
27#include <framework/database/DBObjPtr.h>
28#include <framework/dataobjects/FileMetaData.h>
29#include <framework/datastore/StoreArray.h>
30#include <framework/dbobjects/BeamParameters.h>
31#include <framework/particledb/EvtGenDatabasePDG.h>
32#include <framework/pcore/ProcHandler.h>
33#include <mdst/dbobjects/BeamSpot.h>
34#include <mdst/dataobjects/Track.h>
35#include <tracking/trackFitting/fitter/base/TrackFitter.h>
36#include <tracking/trackFitting/measurementCreator/adder/MeasurementAdder.h>
37
38#include <genfit/FullMeasurement.h>
39#include <genfit/GblFitter.h>
40#include <genfit/KalmanFitterInfo.h>
41#include <genfit/PlanarMeasurement.h>
42#include <genfit/Track.h>
43
44#include <TMath.h>
45#include <TH1F.h>
46#include <TTree.h>
47#include <TDecompSVD.h>
48
49using namespace std;
50using namespace Belle2;
51using namespace alignment;
52
53//-----------------------------------------------------------------
54// Register the Module
55//-----------------------------------------------------------------
56REG_MODULE(MillepedeCollector);
57
58//-----------------------------------------------------------------
59// Implementation
60//-----------------------------------------------------------------
61
63{
65 setDescription("Calibration data collector for Millepede Algorithm");
66
67 // Configure input sample types
68 addParam("tracks", m_tracks, "Names of collections of RecoTracks (already fitted with DAF) for calibration", vector<string>({""}));
69 addParam("particles", m_particles, "Names of particle list of single particles", vector<string>());
70 addParam("vertices", m_vertices,
71 "Name of particle list of (mother) particles with daughters for calibration using vertex constraint", vector<string>());
72 addParam("primaryVertices", m_primaryVertices,
73 "Name of particle list of (mother) particles with daughters for calibration using vertex + IP profile constraint",
74 vector<string>());
75 addParam("twoBodyDecays", m_twoBodyDecays,
76 "Name of particle list of (mother) particles with daughters for calibration using vertex + mass constraint",
77 vector<string>());
78 addParam("primaryTwoBodyDecays", m_primaryTwoBodyDecays,
79 "Name of particle list of (mother) particles with daughters for calibration using vertex + IP profile + kinematics constraint",
80 vector<string>());
81 addParam("primaryMassTwoBodyDecays", m_primaryMassTwoBodyDecays,
82 "Name of particle list of (mother) particles with daughters for calibration using vertex + mass constraint",
83 vector<string>());
84 addParam("primaryMassVertexTwoBodyDecays", m_primaryMassVertexTwoBodyDecays,
85 "Name of particle list of (mother) particles with daughters for calibration using vertex + IP profile + mass constraint",
86 vector<string>());
87
88 addParam("stableParticleWidth", m_stableParticleWidth,
89 "Width (in GeV/c/c) to use for invariant mass constraint for 'stable' particles (like K short). Temporary until proper solution is found.",
90 double(0.002));
91 // Configure output
92 addParam("doublePrecision", m_doublePrecision, "Use double (=true) or single/float (=false) precision for writing binary files",
93 bool(false));
94 addParam("useGblTree", m_useGblTree, "Store GBL trajectories in a tree instead of output to binary files",
95 bool(true));
96 addParam("absFilePaths", m_absFilePaths, "Use absolute paths to remember binary files. Only applies if useGblTree=False",
97 bool(false));
98
99 // Configure global parameters
100 addParam("components", m_components,
101 "Specify which DB objects are calibrated, like ['BeamSpot', 'CDCTimeWalks'] or leave empty to use all components available.",
103 addParam("calibrateVertex", m_calibrateVertex,
104 "For primary vertices / two body decays, beam spot vertex calibration derivatives are added",
105 bool(true));
106 addParam("calibrateKinematics", m_calibrateKinematics,
107 "For primary two body decays, beam spot kinematics calibration derivatives are added",
108 bool(true));
109
110 //Configure GBL fit of individual tracks
111 addParam("externalIterations", m_externalIterations, "Number of external iterations of GBL fitter",
112 int(0));
113 addParam("internalIterations", m_internalIterations, "String defining internal GBL iterations for outlier down-weighting",
114 string(""));
115 addParam("recalcJacobians", m_recalcJacobians, "Up to which external iteration propagation Jacobians should be re-calculated",
116 int(0));
117
118 addParam("minPValue", m_minPValue, "Minimum p-value to write out a (combined) trajectory. Set <0 to write out all.",
119 double(-1.));
120
121 // Configure CDC specific options
122 addParam("fitTrackT0", m_fitTrackT0, "Add local parameter for track T0 fit in GBL",
123 bool(true));
124 addParam("updateCDCWeights", m_updateCDCWeights, "Update L/R weights from previous DAF fit result",
125 bool(true));
126 addParam("minCDCHitWeight", m_minCDCHitWeight, "Minimum (DAF) CDC hit weight for usage by GBL",
127 double(1.0E-6));
128 addParam("minUsedCDCHitFraction", m_minUsedCDCHitFraction, "Minimum used CDC hit fraction to write out a trajectory",
129 double(0.85));
130
131 addParam("hierarchyType", m_hierarchyType, "Type of (VXD only now) hierarchy: 0 = None, 1 = Flat, 2 = Half-Shells, 3 = Full",
132 int(3));
133 addParam("enablePXDHierarchy", m_enablePXDHierarchy, "Enable PXD in hierarchy (flat or full)",
134 bool(true));
135 addParam("enableSVDHierarchy", m_enableSVDHierarchy, "Enable SVD in hierarchy (flat or full)",
136 bool(true));
137
138 addParam("enableWireByWireAlignment", m_enableWireByWireAlignment, "Enable global derivatives for wire-by-wire alignment",
139 bool(false));
140 addParam("enableWireSagging", m_enableWireSagging, "Enable global derivatives for wire sagging",
141 bool(false));
142
143 // Time dependence
144 addParam("events", m_eventNumbers,
145 "List of (event, run, exp) with event numbers at which payloads can change for timedep calibration.",
147 // Time dependence config
148 addParam("timedepConfig", m_timedepConfig,
149 "list{ {list{param1, param2, ...}, list{(ev1, run1, exp1), ...}}, ... }.",
151
152 // Custom mass+width config
153 addParam("customMassConfig", m_customMassConfig,
154 "dict{ list_name: (mass, width), ... } with custom mass and width to use as external measurement.",
156}
157
159{
160 m_eventT0.isOptional();
161
162 if (m_tracks.empty() &&
163 m_particles.empty() &&
164 m_vertices.empty() &&
165 m_primaryVertices.empty() &&
166 m_twoBodyDecays.empty() &&
167 m_primaryTwoBodyDecays.empty() &&
170 B2ERROR("You have to specify either arrays of single tracks or particle lists of single single particles or mothers with vertex constrained daughters.");
171
172 if (!m_tracks.empty()) {
173 for (auto arrayName : m_tracks)
174 continue;
175 // StoreArray<RecoTrack>::required(arrayName);
176 }
177
178 if (!m_particles.empty() || !m_vertices.empty() || !m_primaryVertices.empty()) {
179 // StoreArray<RecoTrack> recoTracks;
180 // StoreArray<Track> tracks;
181 // StoreArray<TrackFitResult> trackFitResults;
182
183 //recoTracks.isRequired();
184 //tracks.isRequired();
185 //trackFitResults.isRequired();
186 }
187
188 for (auto listName : m_particles) {
189 StoreObjPtr<ParticleList> list(listName);
190 //list.isRequired();
191 }
192
193 for (auto listName : m_vertices) {
194 StoreObjPtr<ParticleList> list(listName);
195 //list.isRequired();
196 }
197
198 for (auto listName : m_primaryVertices) {
199 StoreObjPtr<ParticleList> list(listName);
200 //list.isRequired();
201 }
202
203 // Register Mille output
204 registerObject<MilleData>("mille", new MilleData(m_doublePrecision, m_absFilePaths));
205
206 auto gblDataTree = new TTree("GblDataTree", "GblDataTree");
207 gblDataTree->Branch<std::vector<gbl::GblData>>("GblData", &m_currentGblData, 32000, 99);
208 registerObject<TTree>("GblDataTree", gblDataTree);
209
210 registerObject<TH1I>("ndf", new TH1I("ndf", "ndf", 200, 0, 200));
211 registerObject<TH1F>("chi2_per_ndf", new TH1F("chi2_per_ndf", "chi2 divided by ndf", 200, 0., 50.));
212 registerObject<TH1F>("pval", new TH1F("pval", "pval", 100, 0., 1.));
213
214 registerObject<TH1F>("cdc_hit_fraction", new TH1F("cdc_hit_fraction", "cdc_hit_fraction", 100, 0., 1.));
215 registerObject<TH1F>("evt0", new TH1F("evt0", "evt0", 400, -100., 100.));
216
217 // Configure the (VXD) hierarchy before being built
218 if (m_hierarchyType == 0)
219 Belle2::alignment::VXDGlobalParamInterface::s_hierarchyType = VXDGlobalParamInterface::c_None;
220 else if (m_hierarchyType == 1)
221 Belle2::alignment::VXDGlobalParamInterface::s_hierarchyType = VXDGlobalParamInterface::c_Flat;
222 else if (m_hierarchyType == 2)
223 Belle2::alignment::VXDGlobalParamInterface::s_hierarchyType = VXDGlobalParamInterface::c_HalfShells;
224 else if (m_hierarchyType == 3)
225 Belle2::alignment::VXDGlobalParamInterface::s_hierarchyType = VXDGlobalParamInterface::c_Full;
226
229
230 std::vector<EventMetaData> events;
231 for (auto& ev_run_exp : m_eventNumbers) {
232 events.push_back(EventMetaData(std::get<0>(ev_run_exp), std::get<1>(ev_run_exp), std::get<2>(ev_run_exp)));
233 }
234
235 // This will also build the hierarchy for the first time:
236 if (!m_timedepConfig.empty() && m_eventNumbers.empty()) {
237 auto autoEvents = Belle2::alignment::timeline::setupTimedepGlobalLabels(m_timedepConfig);
239 } else if (m_timedepConfig.empty() && !m_eventNumbers.empty()) {
241 } else if (m_timedepConfig.empty() && m_eventNumbers.empty()) {
243 } else {
244 B2ERROR("Cannot set both, event list and timedep config.");
245 }
246
247// Belle2::alignment::GlobalCalibrationManager::getInstance().writeConstraints("constraints.txt");
248
252}
253
255{
257
258 if (!m_useGblTree) {
259 // Open new file on request (at start or after being closed)
260 auto mille = getObjectPtr<MilleData>("mille");
261 if (!mille->isOpen())
262 mille->open(getUniqueMilleName());
263 }
264
265 std::shared_ptr<genfit::GblFitter> gbl(new genfit::GblFitter());
266 double chi2 = -1.;
267 double lostWeight = -1.;
268 int ndf = -1;
269 float evt0 = -9999.;
270
271 for (auto arrayName : m_tracks) {
272 StoreArray<RecoTrack> recoTracks(arrayName);
273 if (!recoTracks.isValid())
274 continue;
275
276 for (auto& recoTrack : recoTracks) {
277
278 if (!fitRecoTrack(recoTrack))
279 continue;
280
281 auto& track = RecoTrackGenfitAccess::getGenfitTrack(recoTrack);
282 if (!track.hasFitStatus())
283 continue;
284 genfit::GblFitStatus* fs = dynamic_cast<genfit::GblFitStatus*>(track.getFitStatus());
285 if (!fs)
286 continue;
287
288 if (!fs->isFittedWithReferenceTrack())
289 continue;
290
291 using namespace gbl;
292 GblTrajectory trajectory(gbl->collectGblPoints(&track, track.getCardinalRep()), fs->hasCurvature());
293
294 trajectory.fit(chi2, ndf, lostWeight);
295 getObjectPtr<TH1I>("ndf")->Fill(ndf);
296 getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
297 getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
298 if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
299 evt0 = m_eventT0->getEventT0();
300 getObjectPtr<TH1F>("evt0")->Fill(evt0);
301 }
302
303 if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(trajectory);
304
305 }
306
307 }
308
309 for (auto listName : m_particles) {
310 StoreObjPtr<ParticleList> list(listName);
311 if (!list.isValid())
312 continue;
313
314 for (unsigned int iParticle = 0; iParticle < list->getListSize(); ++iParticle) {
315 for (auto& track : getParticlesTracks({list->getParticle(iParticle)}, false)) {
316 auto gblfs = dynamic_cast<genfit::GblFitStatus*>(track->getFitStatus());
317
318 gbl::GblTrajectory trajectory(gbl->collectGblPoints(track, track->getCardinalRep()), gblfs->hasCurvature());
319
320 trajectory.fit(chi2, ndf, lostWeight);
321 getObjectPtr<TH1I>("ndf")->Fill(ndf);
322 getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
323 getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
324 if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
325 evt0 = m_eventT0->getEventT0();
326 getObjectPtr<TH1F>("evt0")->Fill(evt0);
327 }
328
329 if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(trajectory);
330
331 }
332 }
333 }
334
335 for (auto listName : m_vertices) {
336 StoreObjPtr<ParticleList> list(listName);
337 if (!list.isValid())
338 continue;
339
340 for (unsigned int iParticle = 0; iParticle < list->getListSize(); ++iParticle) {
341 auto mother = list->getParticle(iParticle);
342 std::vector<std::pair<std::vector<gbl::GblPoint>, TMatrixD> > daughters;
343
344 for (auto& track : getParticlesTracks(mother->getDaughters()))
345 daughters.push_back({
346 gbl->collectGblPoints(track, track->getCardinalRep()),
347 getGlobalToLocalTransform(track->getFittedState()).GetSub(0, 4, 0, 2)
348 });
349
350 if (daughters.size() > 1) {
351 gbl::GblTrajectory combined(daughters);
352
353 combined.fit(chi2, ndf, lostWeight);
354 getObjectPtr<TH1I>("ndf")->Fill(ndf);
355 getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
356 getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
357 if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
358 evt0 = m_eventT0->getEventT0();
359 getObjectPtr<TH1F>("evt0")->Fill(evt0);
360 }
361
362
363 if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(combined);
364
365 B2RESULT("Vertex-constrained fit NDF = " << ndf << " Chi2/NDF = " << chi2 / double(ndf));
366
367 }
368 }
369 }
370
371 for (auto listName : m_primaryVertices) {
372 StoreObjPtr<ParticleList> list(listName);
373 if (!list.isValid())
374 continue;
375
376 for (unsigned int iParticle = 0; iParticle < list->getListSize(); ++iParticle) {
377 auto mother = list->getParticle(iParticle);
378 std::vector<std::pair<std::vector<gbl::GblPoint>, TMatrixD> > daughters;
379
380 TMatrixD extProjection(5, 3);
381 TMatrixD locProjection(3, 5);
382
383 // geometric constaint: 3 common (position) parameters + 3 external (curv., directions) per daughter
384 TMatrixD innerTrafo(5, 3 + 3 * mother->getDaughters().size());
385 unsigned int iCol(3);
386
387 bool first(true);
388 for (auto& track : getParticlesTracks(mother->getDaughters())) {
389 if (first) {
390 // For first trajectory only
391 extProjection = getGlobalToLocalTransform(track->getFittedState()).GetSub(0, 4, 0, 2);
392 locProjection = getLocalToGlobalTransform(track->getFittedState()).GetSub(0, 2, 0, 4);
393 first = false;
394 }
395
396 innerTrafo.Zero();
397 innerTrafo.SetSub(3, 0, getGlobalToLocalTransform(track->getFittedState()).GetSub(3, 4, 0, 2));
398 innerTrafo[0][iCol++] = 1.;
399 innerTrafo[1][iCol++] = 1.;
400 innerTrafo[2][iCol++] = 1.;
401
402
403 daughters.push_back({
404 gbl->collectGblPoints(track, track->getCardinalRep()),
405 innerTrafo
406 });
407 }
408
409 if (daughters.size() > 1) {
410 auto beam = getPrimaryVertexAndCov();
411
412 TMatrixDSym vertexCov(get<TMatrixDSym>(beam));
413 TMatrixDSym vertexPrec(get<TMatrixDSym>(beam).Invert());
414 B2Vector3D vertexResidual = - (B2Vector3D(mother->getVertex()) - get<B2Vector3D>(beam));
415
416 TVectorD extMeasurements(3);
417 extMeasurements[0] = vertexResidual[0];
418 extMeasurements[1] = vertexResidual[1];
419 extMeasurements[2] = vertexResidual[2];
420
421 TMatrixD extDeriv(3, 9);
422 extDeriv.Zero();
423 // beam vertex constraint
424 extDeriv(0, 0) = 1.;
425 extDeriv(1, 1) = 1.;
426 extDeriv(2, 2) = 1.;
427
428 if (m_calibrateVertex) {
429 TMatrixD derivatives(3, 3);
430 derivatives.Zero();
431 derivatives(0, 0) = 1.;
432 derivatives(1, 1) = 1.;
433 derivatives(2, 2) = 1.;
434
435 std::vector<int> labels;
436 GlobalLabel label = GlobalLabel::construct<BeamSpot>(0, 0);
437 labels.push_back(label.setParameterId(1));
438 labels.push_back(label.setParameterId(2));
439 labels.push_back(label.setParameterId(3));
440
441 // Allow to disable BeamSpot externally
442 alignment::GlobalDerivatives globals(labels, derivatives);
443 // Add derivatives for vertex calibration to first point of first trajectory
444 // NOTE: use GlobalDerivatives operators vector<int> and TMatrixD which filter
445 // the derivatives to not pass those with zero labels (useful to get rid of some params)
446 std::vector<int> lab(globals); TMatrixD der(globals);
447
448 // Transformation from local system at (vertex) point to global (vx,vy,vz)
449 // of the (decay) vertex
450 //
451 // d(q/p,u',v',u,v)/d(vy,vy,vz) = dLocal_dExt
452 //
453 //
454 // Note its transpose is its "inverse" in the sense that
455 //
456 // dloc/dext * (dloc/dext)^T = diag(0, 0, 0, 0, 1, 1)
457 //
458 //
459 // N.B. typical dLocal_dExt matrix (5x3):
460 //
461 // | 0 | 1 | 2 |
462 // --------------------------------------------
463 // 0 | 0 0 0
464 // 1 | 0 0 0
465 // 2 | 0 0 0
466 // 3 | -0.02614 -0.9997 0
467 // 4 | 0 0 1
468 //
469 // Therefore one can simplify things by only taking the last two rows/columns in vectors/matrices
470 // and vertex measurement can be expressed as standard 2D measurement in GBL.
471 //
472 TMatrixD dLocal_dExt = extProjection;
473 TMatrixD dExt_dLocal = locProjection;
474
475 TVectorD locRes = dLocal_dExt * extMeasurements;
476 // Do not use inverted covariance - seems to have issues with numeric precision
477 TMatrixD locCov = dLocal_dExt * vertexCov * dExt_dLocal;
478 // Invert here only the 2D sub-matrix (rest is zero due to the foŕm of dLocal_dExt)
479 TMatrixD locPrec = locCov.GetSub(3, 4, 3, 4).Invert();
480 TMatrixDSym locPrec2D(2); locPrec2D.Zero();
481 for (int i = 0; i < 2; ++i)
482 for (int j = 0; j < 2; ++j)
483 locPrec2D(i, j) = locPrec(i, j);
484
485 // Take the 2 last components also for residuals and global derivatives
486 // (in local system of vertex point - defined during fitRecoTrack(..., particle) and using
487 // the (hopefully) updated momentum and position seed after vertex fit by modularAnalysis
488 TVectorD locRes2D = locRes.GetSub(3, 4);
489 TMatrixD locDerivs2D = (extProjection * der).GetSub(3, 4, 0, 2);
490
491 // Attach the primary beamspot vertex position as a measurement at 1st point
492 // of first trajectory (and optionally also the global derivatives for beamspot alignment
493 daughters[0].first[0].addMeasurement(locRes2D, locPrec2D);
494 if (!lab.empty()) {
495 daughters[0].first[0].addGlobals(lab, locDerivs2D);
496 }
497
498 gbl::GblTrajectory combined(daughters);
499 //combined.printTrajectory(100);
500 //combined.printPoints(100);
501
502 combined.fit(chi2, ndf, lostWeight);
503 getObjectPtr<TH1I>("ndf")->Fill(ndf);
504 getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
505 getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
506 if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
507 evt0 = m_eventT0->getEventT0();
508 getObjectPtr<TH1F>("evt0")->Fill(evt0);
509 }
510
511 if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(combined);
512 B2RESULT("Beam vertex constrained fit results NDF = " << ndf << " Chi2/NDF = " << chi2 / double(ndf));
513
514 } else {
515
516 gbl::GblTrajectory combined(daughters, extDeriv, extMeasurements, vertexPrec);
517
518 combined.fit(chi2, ndf, lostWeight);
519 getObjectPtr<TH1I>("ndf")->Fill(ndf);
520 getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
521 getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
522 if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
523 evt0 = m_eventT0->getEventT0();
524 getObjectPtr<TH1F>("evt0")->Fill(evt0);
525 }
526
527 if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(combined);
528
529 B2RESULT("Beam vertex constrained fit results NDF = " << ndf << " Chi2/NDF = " << chi2 / double(ndf));
530
531 }
532 }
533 }
534 }
535
536 for (auto listName : m_twoBodyDecays) {
537 StoreObjPtr<ParticleList> list(listName);
538 if (!list.isValid())
539 continue;
540
541 for (unsigned int iParticle = 0; iParticle < list->getListSize(); ++iParticle) {
542
543 auto mother = list->getParticle(iParticle);
544 auto track12 = getParticlesTracks(mother->getDaughters());
545 if (track12.size() != 2) {
546 B2ERROR("Did not get 2 fitted tracks. Skipping this mother.");
547 continue;
548 }
549
550 auto pdgdb = EvtGenDatabasePDG::Instance();
551 double motherMass = mother->getPDGMass();
552 double motherWidth = pdgdb->GetParticle(mother->getPDGCode())->Width();
553
554 updateMassWidthIfSet(listName, motherMass, motherWidth);
555
556 //TODO: what to take as width for "real" particles? -> make a param for default detector mass resolution??
557 if (motherWidth == 0.) {
558 motherWidth = m_stableParticleWidth * Unit::GeV;
559 B2WARNING("Using artificial width for " << pdgdb->GetParticle(mother->getPDGCode())->GetName() << " : " << motherWidth << " GeV");
560 }
561
562 auto dfdextPlusMinus = getTwoBodyToLocalTransform(*mother, motherMass);
563 std::vector<std::pair<std::vector<gbl::GblPoint>, TMatrixD> > daughters;
564
565 daughters.push_back({gbl->collectGblPoints(track12[0], track12[0]->getCardinalRep()), dfdextPlusMinus.first});
566 daughters.push_back({gbl->collectGblPoints(track12[1], track12[1]->getCardinalRep()), dfdextPlusMinus.second});
567
568 TMatrixDSym massPrec(1); massPrec(0, 0) = 1. / motherWidth / motherWidth;
569 TVectorD massResidual(1); massResidual = - (mother->getMass() - motherMass);
570
571 TVectorD extMeasurements(1);
572 extMeasurements[0] = massResidual[0];
573
574 TMatrixD extDeriv(1, 9);
575 extDeriv.Zero();
576 extDeriv(0, 8) = 1.;
577
578 gbl::GblTrajectory combined(daughters, extDeriv, extMeasurements, massPrec);
579
580 combined.fit(chi2, ndf, lostWeight);
581 //combined.printTrajectory(1000);
582 //combined.printPoints(1000);
583 getObjectPtr<TH1I>("ndf")->Fill(ndf);
584 getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
585 getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
586 if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
587 evt0 = m_eventT0->getEventT0();
588 getObjectPtr<TH1F>("evt0")->Fill(evt0);
589 }
590
591
592 B2RESULT("Mass(PDG) + vertex constrained fit results NDF = " << ndf << " Chi2/NDF = " << chi2 / double(ndf));
593
594 if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(combined);
595
596 }
597 }
598
599 for (auto listName : m_primaryMassTwoBodyDecays) {
600 StoreObjPtr<ParticleList> list(listName);
601 if (!list.isValid())
602 continue;
603
605
606 double motherMass = beam->getMass();
607 double motherWidth = sqrt((beam->getCovHER() + beam->getCovLER())(0, 0));
608
609 updateMassWidthIfSet(listName, motherMass, motherWidth);
610
611 for (unsigned int iParticle = 0; iParticle < list->getListSize(); ++iParticle) {
612
613 auto mother = list->getParticle(iParticle);
614 auto track12 = getParticlesTracks(mother->getDaughters());
615 if (track12.size() != 2) {
616 B2ERROR("Did not get 2 fitted tracks. Skipping this mother.");
617 continue;
618 }
619
620 auto dfdextPlusMinus = getTwoBodyToLocalTransform(*mother, motherMass);
621 std::vector<std::pair<std::vector<gbl::GblPoint>, TMatrixD> > daughters;
622
623 daughters.push_back({gbl->collectGblPoints(track12[0], track12[0]->getCardinalRep()), dfdextPlusMinus.first});
624 daughters.push_back({gbl->collectGblPoints(track12[1], track12[1]->getCardinalRep()), dfdextPlusMinus.second});
625
626 TMatrixDSym massPrec(1); massPrec(0, 0) = 1. / motherWidth / motherWidth;
627 TVectorD massResidual(1); massResidual = - (mother->getMass() - motherMass);
628
629 TVectorD extMeasurements(1);
630 extMeasurements[0] = massResidual[0];
631
632 TMatrixD extDeriv(1, 9);
633 extDeriv.Zero();
634 extDeriv(0, 8) = 1.;
635
636 gbl::GblTrajectory combined(daughters, extDeriv, extMeasurements, massPrec);
637
638 combined.fit(chi2, ndf, lostWeight);
639 getObjectPtr<TH1I>("ndf")->Fill(ndf);
640 getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
641 getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
642 if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
643 evt0 = m_eventT0->getEventT0();
644 getObjectPtr<TH1F>("evt0")->Fill(evt0);
645 }
646
647
648 B2RESULT("Mass constrained fit results NDF = " << ndf << " Chi2/NDF = " << chi2 / double(ndf));
649
650 if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(combined);
651
652 }
653 }
654
655 for (auto listName : m_primaryMassVertexTwoBodyDecays) {
656 StoreObjPtr<ParticleList> list(listName);
657 if (!list.isValid())
658 continue;
659
661
662 double motherMass = beam->getMass();
663 double motherWidth = sqrt((beam->getCovHER() + beam->getCovLER())(0, 0));
664
665 updateMassWidthIfSet(listName, motherMass, motherWidth);
666
667 for (unsigned int iParticle = 0; iParticle < list->getListSize(); ++iParticle) {
668
669 auto mother = list->getParticle(iParticle);
670 auto track12 = getParticlesTracks(mother->getDaughters());
671 if (track12.size() != 2) {
672 B2ERROR("Did not get 2 fitted tracks. Skipping this mother.");
673 continue;
674 }
675
676 auto dfdextPlusMinus = getTwoBodyToLocalTransform(*mother, motherMass);
677 std::vector<std::pair<std::vector<gbl::GblPoint>, TMatrixD> > daughters;
678
679 daughters.push_back({gbl->collectGblPoints(track12[0], track12[0]->getCardinalRep()), dfdextPlusMinus.first});
680 daughters.push_back({gbl->collectGblPoints(track12[1], track12[1]->getCardinalRep()), dfdextPlusMinus.second});
681
682 TMatrixDSym vertexPrec(get<TMatrixDSym>(getPrimaryVertexAndCov()).Invert());
683 B2Vector3D vertexResidual = - (B2Vector3D(mother->getVertex()) - get<B2Vector3D>(getPrimaryVertexAndCov()));
684
685 TMatrixDSym massPrec(1); massPrec(0, 0) = 1. / motherWidth / motherWidth;
686 TVectorD massResidual(1); massResidual = - (mother->getMass() - motherMass);
687
688 TMatrixDSym extPrec(4); extPrec.Zero();
689 extPrec.SetSub(0, 0, vertexPrec);
690 extPrec(3, 3) = massPrec(0, 0);
691
692 TVectorD extMeasurements(4);
693 extMeasurements[0] = vertexResidual[0];
694 extMeasurements[1] = vertexResidual[1];
695 extMeasurements[2] = vertexResidual[2];
696 extMeasurements[3] = massResidual[0];
697
698 TMatrixD extDeriv(4, 9);
699 extDeriv.Zero();
700 extDeriv(0, 0) = 1.;
701 extDeriv(1, 1) = 1.;
702 extDeriv(2, 2) = 1.;
703 extDeriv(3, 8) = 1.;
704
705 gbl::GblTrajectory combined(daughters, extDeriv, extMeasurements, extPrec);
706
707 combined.fit(chi2, ndf, lostWeight);
708 getObjectPtr<TH1I>("ndf")->Fill(ndf);
709 getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
710 getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
711 if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
712 evt0 = m_eventT0->getEventT0();
713 getObjectPtr<TH1F>("evt0")->Fill(evt0);
714 }
715
716
717 if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(combined);
718
719 B2RESULT("Mass + vertex constrained fit results NDF = " << ndf << " Chi2/NDF = " << chi2 / double(ndf));
720
721 }
722 }
723
724 for (auto listName : m_primaryTwoBodyDecays) {
725 B2WARNING("This should NOT be used for production of calibration constants for the real detector (yet)!");
726
727 StoreObjPtr<ParticleList> list(listName);
728 if (!list.isValid())
729 continue;
730
732
733 // For the error of invariant mass M = 2 * sqrt(E_HER * E_LER) (for m_e ~ 0)
734 double M = beam->getMass();
735 double E_HER = beam->getHER().E();
736 double E_LER = beam->getLER().E();
737
738 double pz = beam->getHER().Pz() + beam->getLER().Pz();
739 double E = (beam->getHER() + beam->getLER()).E();
740
741 double motherMass = beam->getMass();
742 double motherWidth = sqrt((E_HER / M) * (E_HER / M) * beam->getCovLER()(0, 0) + (E_LER / M) * (E_LER / M) * beam->getCovHER()(0,
743 0));
744
745 updateMassWidthIfSet(listName, motherMass, motherWidth);
746
747 for (unsigned int iParticle = 0; iParticle < list->getListSize(); ++iParticle) {
748
749 B2WARNING("Two body decays with full kinematic constraint not yet correct - need to resolve strange covariance provided by BeamParameters!");
750
751 auto mother = list->getParticle(iParticle);
752
753 auto track12 = getParticlesTracks(mother->getDaughters());
754 if (track12.size() != 2) {
755 B2ERROR("Did not get exactly 2 fitted tracks. Skipping this mother in list " << listName);
756 continue;
757 }
758
759 auto dfdextPlusMinus = getTwoBodyToLocalTransform(*mother, motherMass);
760 std::vector<std::pair<std::vector<gbl::GblPoint>, TMatrixD> > daughters;
761
762 daughters.push_back({gbl->collectGblPoints(track12[0], track12[0]->getCardinalRep()), dfdextPlusMinus.first});
763 daughters.push_back({gbl->collectGblPoints(track12[1], track12[1]->getCardinalRep()), dfdextPlusMinus.second});
764
765 TMatrixDSym extCov(7); extCov.Zero();
766
767 // 3x3 IP vertex covariance
768 extCov.SetSub(0, 0, get<TMatrixDSym>(getPrimaryVertexAndCov()));
769
770 // 3x3 boost vector covariance
771 //NOTE: BeamSpot return covarince in variables (E, theta_x, theta_y)
772 // We need to transform it to our variables (px, py, pz)
773
774 TMatrixD dBoost_dVect(3, 3);
775 dBoost_dVect(0, 0) = 0.; dBoost_dVect(0, 1) = 1. / pz; dBoost_dVect(0, 2) = 0.;
776 dBoost_dVect(1, 0) = 0.; dBoost_dVect(1, 1) = 0.; dBoost_dVect(1, 2) = 1. / pz;
777 dBoost_dVect(2, 0) = pz / E; dBoost_dVect(2, 1) = 0.; dBoost_dVect(2, 2) = 0.;
778
779 TMatrixD dVect_dBoost(3, 3);
780 dVect_dBoost(0, 0) = 0.; dVect_dBoost(0, 1) = 0.; dVect_dBoost(0, 2) = E / pz;
781 dVect_dBoost(1, 0) = pz; dVect_dBoost(1, 1) = 0.; dVect_dBoost(1, 2) = 0.;
782 dVect_dBoost(2, 0) = 0.; dVect_dBoost(2, 1) = pz; dVect_dBoost(2, 2) = 0.;
783
784 TMatrixD covBoost(3, 3);
785 for (int i = 0; i < 3; ++i) {
786 for (int j = i; j < 3; ++j) {
787 covBoost(j, i) = covBoost(i, j) = (beam->getCovHER() + beam->getCovLER())(i, j);
788 }
789 }
790 //TODO: Temporary fix: if theta_x, theta_y covariance is zero, use arbitrary 10mrad^2
791// if (covBoost(1, 1) == 0.) covBoost(1, 1) = 1.;
792// if (covBoost(2, 2) == 0.) covBoost(2, 2) = 1.;
793 if (covBoost(1, 1) == 0.) covBoost(1, 1) = 1.e-4;
794 if (covBoost(2, 2) == 0.) covBoost(2, 2) = 1.e-4;
795
796 TMatrixD covVect = dBoost_dVect * covBoost * dVect_dBoost;
797
798 extCov.SetSub(3, 3, covVect);
799
800 extCov(6, 6) = motherWidth * motherWidth;
801 auto extPrec = extCov; extPrec.Invert();
802
803 TVectorD extMeasurements(7);
804 extMeasurements[0] = - (B2Vector3D(mother->getVertex()) - get<B2Vector3D>(getPrimaryVertexAndCov()))[0];
805 extMeasurements[1] = - (B2Vector3D(mother->getVertex()) - get<B2Vector3D>(getPrimaryVertexAndCov()))[1];
806 extMeasurements[2] = - (B2Vector3D(mother->getVertex()) - get<B2Vector3D>(getPrimaryVertexAndCov()))[2];
807 extMeasurements[3] = - (B2Vector3D(mother->getMomentum()) - (beam->getHER().Vect() + beam->getLER().Vect()))[0];
808 extMeasurements[4] = - (B2Vector3D(mother->getMomentum()) - (beam->getHER().Vect() + beam->getLER().Vect()))[1];
809 extMeasurements[5] = - (B2Vector3D(mother->getMomentum()) - (beam->getHER().Vect() + beam->getLER().Vect()))[2];
810 extMeasurements[6] = - (mother->getMass() - motherMass);
811
812 B2INFO("mother mass = " << mother->getMass() << " and beam mass = " << beam->getMass());
813
814 TMatrixD extDeriv(7, 9);
815 extDeriv.Zero();
816 // beam vertex constraint
817 extDeriv(0, 0) = 1.;
818 extDeriv(1, 1) = 1.;
819 extDeriv(2, 2) = 1.;
820 // beam kinematics constraint
821 extDeriv(3, 3) = 1.;
822 extDeriv(4, 4) = 1.;
823 extDeriv(5, 5) = 1.;
824 // beam inv. mass constraint
825 extDeriv(6, 8) = 1;
826
828 B2WARNING("Primary vertex+kinematics calibration not (yet?) fully implemented!");
829 B2WARNING("This code is highly experimental and has (un)known issues!");
830
831 // up to d(x,y,z,px,py,pz,theta,phi,M)/d(vx,vy,vz,theta_x,theta_y,E)
832 TMatrixD derivatives(9, 6);
833 std::vector<int> labels;
834 derivatives.Zero();
835
836 if (m_calibrateVertex) {
837 derivatives(0, 0) = 1.;
838 derivatives(1, 1) = 1.;
839 derivatives(2, 2) = 1.;
840 GlobalLabel label = GlobalLabel::construct<BeamSpot>(0, 0);
841 labels.push_back(label.setParameterId(1));
842 labels.push_back(label.setParameterId(2));
843 labels.push_back(label.setParameterId(3));
844 } else {
845 labels.push_back(0);
846 labels.push_back(0);
847 labels.push_back(0);
848 }
849
851 derivatives(3, 3) = mother->getMomentumMagnitude();
852 derivatives(4, 4) = mother->getMomentumMagnitude();
853 derivatives(8, 5) = (beam->getLER().E() + beam->getHER().E()) / beam->getMass();
854
855 GlobalLabel label = GlobalLabel::construct<BeamSpot>(0, 0);
856 labels.push_back(label.setParameterId(4)); //theta_x
857 labels.push_back(label.setParameterId(5)); //theta_y
858 labels.push_back(label.setParameterId(6)); //E
859
860 } else {
861 labels.push_back(0);
862 labels.push_back(0);
863 labels.push_back(0);
864 }
865
866 // Allow to disable BeamSpot externally
867 alignment::GlobalDerivatives globals(labels, derivatives);
868
869 // Add derivatives for vertex calibration to first point of first trajectory
870 // NOTE: use GlobalDerivatives operators vector<int> and TMatrixD which filter
871 // the derivatives to not pass those with zero labels (useful to get rid of some params)
872 std::vector<int> lab(globals); TMatrixD der(globals);
873
874 // I want: dlocal/dext = dlocal/dtwobody * dtwobody/dext = dfdextPlusMinus * dtwobody/dext
875 TMatrixD dTwoBody_dExt(9, 7);
876 dTwoBody_dExt.Zero();
877 // beam vertex constraint
878 dTwoBody_dExt(0, 0) = 1.;
879 dTwoBody_dExt(1, 1) = 1.;
880 dTwoBody_dExt(2, 2) = 1.;
881 // beam kinematics constraint
882 dTwoBody_dExt(3, 3) = 1.;
883 dTwoBody_dExt(4, 4) = 1.;
884 dTwoBody_dExt(5, 5) = 1.;
885 // beam inv. mass constraint
886 dTwoBody_dExt(8, 6) = 1.;
887
888 const TMatrixD dLocal_dExt = dfdextPlusMinus.first * dTwoBody_dExt;
889 TMatrixD dLocal_dExt_T = dLocal_dExt; dLocal_dExt_T.T();
890
891 // The 5x7 transformation matrix d(q/p,u',v',u,v)/d(vx,vy,vz,px,py,pz,M) needs to be "inverted"
892 // to transform the covariance of the beamspot and boost vector of SuperKEKB into the local system
893 // of one GBL point - such that Millepede can align the beamspot (or even beam kinematics) if requested.
894 //
895 // I tested also other methods, but only the Singular Value Decomposition gives nice-enough results,
896 // with almost no code:
897 //
898 TDecompSVD svd(dLocal_dExt_T);
899 TMatrixD dExt_dLocal = svd.Invert().T();
900 //
901 // (dLocal_dExt * dExt_dLocal).Print(); // Check how close we are to unit matrix
902 //
903 // 5x5 matrix is as follows
904 //
905 // | 0 | 1 | 2 | 3 | 4 |
906 // ----------------------------------------------------------------------
907 // 0 | 1 -2.58e-17 6.939e-18 1.571e-17 -1.649e-19
908 // 1 | 1.787e-14 1 5.135e-16 -3.689e-16 -2.316e-18
909 // 2 | -1.776e-15 -7.806e-17 1 5.636e-17 6.193e-18
910 // 3 | -2.453e-15 7.26e-18 2.009e-16 1 -1.14e-16
911 // 4 | -1.689e-14 -9.593e-17 -2.317e-15 -3.396e-17 1
912 //
913 // It took me half a day to find out how to do this with 2 lines of code (3 with the include).
914 // Source: ROOT macro example - actually found at:
915 // <https://root.cern.ch/root/html/tutorials/matrix/solveLinear.C.html>
916 for (int i = 0; i < 7; ++i) {
917 for (int j = 0; j < 5; ++j) {
918 if (fabs(dExt_dLocal(i, j)) < 1.e-6)
919 dExt_dLocal(i, j) = 0.;
920 }
921 }
922 const TVectorD locRes = dLocal_dExt * extMeasurements;
923 const TMatrixD locPrec = dLocal_dExt * extPrec * dExt_dLocal;
924
925 TMatrixDSym locPrecSym(5); locPrecSym.Zero();
926 for (int i = 0; i < 5; ++i) {
927 for (int j = i; j < 5; ++j) {
928 //locPrecSym(j, i) = locPrecSym(i, j) = locPrec(i, j);
929 locPrecSym(j, i) = locPrecSym(i, j) = (fabs(locPrec(i, j)) > 1.e-6) ? locPrec(i, j) : 0.;
930 }
931 }
932
933 daughters[0].first[0].addMeasurement(locRes, locPrecSym);
934 if (!lab.empty())
935 daughters[0].first[0].addGlobals(lab, dfdextPlusMinus.first * der);
936
937 //TODO: Understand this: either find a bug somewhere or improve the parametrization or .... ?
938 // This should be enough, but the parametrization seems to fail for nearly horizontal pairs...
939 //gbl::GblTrajectory combined(daughters);
940 // This should not be needed, it actually seems to make worse Chi2/NDF, but GBL does not fail.
941 // The measurement added just to be able to add the global derivatives (done just above) is redundant
942 // to the external measurement added here:
943 gbl::GblTrajectory combined(daughters, extDeriv, extMeasurements, extPrec);
944 //combined.printTrajectory(1000);
945 //combined.printPoints(1000);
946
947 combined.fit(chi2, ndf, lostWeight);
948 getObjectPtr<TH1I>("ndf")->Fill(ndf);
949 getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
950 getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
951 if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
952 evt0 = m_eventT0->getEventT0();
953 getObjectPtr<TH1F>("evt0")->Fill(evt0);
954 }
955
956
957 B2RESULT("Full kinematic-constrained fit (calibration version) results NDF = " << ndf << " Chi2/NDF = " << chi2 / double(ndf));
958
959 if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(combined);
960
961 } else {
962
963 gbl::GblTrajectory combined(daughters, extDeriv, extMeasurements, extPrec);
964 //combined.printTrajectory(1000);
965 //combined.printPoints(1000);
966
967 combined.fit(chi2, ndf, lostWeight);
968 getObjectPtr<TH1I>("ndf")->Fill(ndf);
969 getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
970 getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
971 if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
972 evt0 = m_eventT0->getEventT0();
973 getObjectPtr<TH1F>("evt0")->Fill(evt0);
974 }
975
976
977 B2RESULT("Full kinematic-constrained fit results NDF = " << ndf << " Chi2/NDF = " << chi2 / double(ndf));
978
979 if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(combined);
980 }
981 }
982 }
983}
984
986{
987 // We close the file at end of run, producing
988 // one file per run (and process id) which is more
989 // convenient than one large binary block.
990 auto mille = getObjectPtr<MilleData>("mille");
991 if (mille->isOpen())
992 mille->close();
993}
994
996{
998
1000 if (!fileMetaData.isValid()) {
1001 B2ERROR("Cannot register binaries in FileCatalog.");
1002 return;
1003 }
1004
1005
1006 const std::vector<string> parents = {fileMetaData->getLfn()};
1007 for (auto binary : getObjectPtr<MilleData>("mille")->getFiles()) {
1008 FileMetaData milleMetaData(*fileMetaData);
1009 // We reset filename to be set directly by the registerFile procedure
1010 milleMetaData.setLfn("");
1011 milleMetaData.setParents(parents);
1012 FileCatalog::Instance().registerFile(binary, milleMetaData);
1013 }
1014
1015}
1016
1017void MillepedeCollectorModule::storeTrajectory(gbl::GblTrajectory& trajectory)
1018{
1019 if (m_useGblTree) {
1020 if (trajectory.isValid())
1021 m_currentGblData = trajectory.getData();
1022 else
1023 m_currentGblData.clear();
1024
1025 if (!m_currentGblData.empty())
1026 getObjectPtr<TTree>("GblDataTree")->Fill();
1027 } else {
1028 getObjectPtr<MilleData>("mille")->fill(trajectory);
1029 }
1030}
1031
1033{
1034 string name = getName();
1035
1036 name += "-e" + to_string(m_evtMetaData->getExperiment());
1037 name += "-r" + to_string(m_evtMetaData->getRun());
1038 name += "-ev" + to_string(m_evtMetaData->getEvent());
1039
1041 name += "-pid" + to_string(ProcHandler::EvtProcID());
1042
1043 name += ".mille";
1044
1045 return name;
1046}
1047
1049{
1050 try {
1051 // For already fitted tracks, try to get fitted (DAF) weights for CDC
1052 if (m_updateCDCWeights && recoTrack.getNumberOfCDCHits() && recoTrack.getTrackFitStatus()
1053 && recoTrack.getTrackFitStatus()->isFitted()) {
1054 double sumCDCWeights = recoTrack.getNumberOfCDCHits(); // start with full weights
1055 // Do the hits synchronisation
1056 auto relatedRecoHitInformation =
1058
1059 for (RecoHitInformation& recoHitInformation : relatedRecoHitInformation) {
1060
1061 if (recoHitInformation.getFlag() == RecoHitInformation::c_pruned) {
1062 B2FATAL("Found pruned point in RecoTrack. Pruned tracks cannot be used in MillepedeCollector.");
1063 }
1064
1065 if (recoHitInformation.getTrackingDetector() != RecoHitInformation::c_CDC) continue;
1066
1067 const genfit::TrackPoint* trackPoint = recoTrack.getCreatedTrackPoint(&recoHitInformation);
1068 if (trackPoint) {
1069 if (not trackPoint->hasFitterInfo(recoTrack.getCardinalRepresentation()))
1070 continue;
1071 auto kalmanFitterInfo = dynamic_cast<genfit::KalmanFitterInfo*>(trackPoint->getFitterInfo());
1072 if (not kalmanFitterInfo) {
1073 continue;
1074 } else {
1075 std::vector<double> weights = kalmanFitterInfo->getWeights();
1076 if (weights.size() == 2) {
1077 if (weights.at(0) > weights.at(1))
1078 recoHitInformation.setRightLeftInformation(RecoHitInformation::c_left);
1079 else if (weights.at(0) < weights.at(1))
1080 recoHitInformation.setRightLeftInformation(RecoHitInformation::c_right);
1081
1082 double weightLR = weights.at(0) + weights.at(1);
1083 if (weightLR < m_minCDCHitWeight) recoHitInformation.setUseInFit(false);
1084 sumCDCWeights += weightLR - 1.; // reduce weight sum if weightLR<1
1085 }
1086 }
1087 }
1088 }
1089
1090 double usedCDCHitFraction = sumCDCWeights / double(recoTrack.getNumberOfCDCHits());
1091 getObjectPtr<TH1F>("cdc_hit_fraction")->Fill(usedCDCHitFraction);
1092 if (usedCDCHitFraction < m_minUsedCDCHitFraction)
1093 return false;
1094 }
1095 } catch (...) {
1096 B2ERROR("Error in checking DAF weights from previous fit to resolve hit ambiguity. Why? Failed fit points in DAF? Skip track to be sure.");
1097 return false;
1098 }
1099
1100 std::shared_ptr<genfit::GblFitter> gbl(new genfit::GblFitter());
1101 gbl->setOptions(m_internalIterations, true, true, m_externalIterations, m_recalcJacobians);
1102 gbl->setTrackSegmentController(new GblMultipleScatteringController);
1103
1104 MeasurementAdder factory("", "", "", "", "");
1105
1106 // We need the store arrays
1112
1113 // Create the genfit::MeasurementFactory
1114 genfit::MeasurementFactory<genfit::AbsMeasurement> genfitMeasurementFactory;
1115
1116 // Add producer for alignable RecoHits to factory
1117 if (pxdHits.isOptional()) {
1118 genfit::MeasurementProducer <RecoHitInformation::UsedPXDHit, AlignablePXDRecoHit>* PXDProducer = new genfit::MeasurementProducer
1120 genfitMeasurementFactory.addProducer(Const::PXD, PXDProducer);
1121 }
1122
1123 if (svdHits.isOptional()) {
1124 genfit::MeasurementProducer <RecoHitInformation::UsedSVDHit, AlignableSVDRecoHit>* SVDProducer = new genfit::MeasurementProducer
1126 genfitMeasurementFactory.addProducer(Const::SVD, SVDProducer);
1127 }
1128
1129 if (cdcHits.isOptional()) {
1130 genfit::MeasurementProducer <RecoHitInformation::UsedCDCHit, AlignableCDCRecoHit>* CDCProducer = new genfit::MeasurementProducer
1132 genfitMeasurementFactory.addProducer(Const::CDC, CDCProducer);
1133 }
1134
1135 if (bklmHits.isOptional()) {
1136 genfit::MeasurementProducer <RecoHitInformation::UsedBKLMHit, AlignableBKLMRecoHit>* BKLMProducer = new genfit::MeasurementProducer
1138 genfitMeasurementFactory.addProducer(Const::BKLM, BKLMProducer);
1139 }
1140
1141 if (eklmHits.isOptional()) {
1142 genfit::MeasurementProducer <RecoHitInformation::UsedEKLMHit, AlignableEKLMRecoHit>* EKLMProducer = new genfit::MeasurementProducer
1144 genfitMeasurementFactory.addProducer(Const::EKLM, EKLMProducer);
1145 }
1146
1147
1148 // Create the measurement creators
1149 std::vector<std::shared_ptr<PXDBaseMeasurementCreator>> pxdMeasurementCreators = { std::shared_ptr<PXDBaseMeasurementCreator>(new PXDCoordinateMeasurementCreator(genfitMeasurementFactory)) };
1150 std::vector<std::shared_ptr<SVDBaseMeasurementCreator>> svdMeasurementCreators = { std::shared_ptr<SVDBaseMeasurementCreator>(new SVDCoordinateMeasurementCreator(genfitMeasurementFactory)) };
1151 // TODO: Create a new MeasurementCreator based on SVDBaseMeasurementCreator (or on SVDCoordinateMeasurementCreator), which does the combination on the fly.
1152
1153 std::vector<std::shared_ptr<CDCBaseMeasurementCreator>> cdcMeasurementCreators = { std::shared_ptr<CDCBaseMeasurementCreator>(new CDCCoordinateMeasurementCreator(genfitMeasurementFactory)) };
1154 std::vector<std::shared_ptr<BKLMBaseMeasurementCreator>> bklmMeasurementCreators = { std::shared_ptr<BKLMBaseMeasurementCreator>(new BKLMCoordinateMeasurementCreator(genfitMeasurementFactory)) };
1155 std::vector<std::shared_ptr<EKLMBaseMeasurementCreator>> eklmMeasurementCreators = { std::shared_ptr<EKLMBaseMeasurementCreator>(new EKLMCoordinateMeasurementCreator(genfitMeasurementFactory)) };
1156
1157 // TODO: Or put it in here and leave the svdMeasurementCreators empty.
1158 std::vector<std::shared_ptr<BaseMeasurementCreator>> additionalMeasurementCreators = {};
1159 factory.resetMeasurementCreators(pxdMeasurementCreators, svdMeasurementCreators, cdcMeasurementCreators, bklmMeasurementCreators,
1160 eklmMeasurementCreators, additionalMeasurementCreators);
1161 factory.addMeasurements(recoTrack);
1162
1163 auto& gfTrack = RecoTrackGenfitAccess::getGenfitTrack(recoTrack);
1164
1165 int currentPdgCode = TrackFitter::createCorrectPDGCodeForChargedStable(Const::muon, recoTrack);
1166 if (particle)
1167 currentPdgCode = particle->getPDGCode();
1168
1169 genfit::AbsTrackRep* trackRep = RecoTrackGenfitAccess::createOrReturnRKTrackRep(recoTrack, currentPdgCode);
1170 gfTrack.setCardinalRep(gfTrack.getIdForRep(trackRep));
1171
1172 if (particle) {
1173 B2Vector3D vertexPos = particle->getVertex();
1174 B2Vector3D vertexMom = particle->getMomentum();
1175 gfTrack.setStateSeed(vertexPos, vertexMom);
1176
1177 genfit::StateOnPlane vertexSOP(gfTrack.getCardinalRep());
1178 B2Vector3D vertexRPhiDir(vertexPos[0], vertexPos[1], 0);
1179 B2Vector3D vertexZDir(0, 0, vertexPos[2]);
1180 //FIXME: This causes problem to current GBL version in genfit -> needs update of GBL to re-enable
1181 // genfit::SharedPlanePtr vertexPlane(new genfit::DetPlane(vertexPos, vertexRPhiDir, vertexZDir));
1182 //This works instead fine:
1183 genfit::SharedPlanePtr vertexPlane(new genfit::DetPlane(vertexPos, vertexMom));
1184
1185 vertexSOP.setPlane(vertexPlane);
1186 vertexSOP.setPosMom(vertexPos, vertexMom);
1187 TMatrixDSym vertexCov(5);
1188 vertexCov.UnitMatrix();
1189 // By using negative covariance no measurement is added to GBL. But this first point
1190 // is then used as additional point in trajectory at the assumed point of its fitted vertex
1191 vertexCov *= -1.;
1192 genfit::MeasuredStateOnPlane mop(vertexSOP, vertexCov);
1193 genfit::FullMeasurement* vertex = new genfit::FullMeasurement(mop, Const::IR);
1194 gfTrack.insertMeasurement(vertex, 0);
1195 }
1196
1197 try {
1198 for (unsigned int i = 0; i < gfTrack.getNumPoints() - 1; ++i) {
1199 //if (gfTrack.getPointWithMeasurement(i)->getNumRawMeasurements() != 1)
1200 // continue;
1201 genfit::PlanarMeasurement* planarMeas1 = dynamic_cast<genfit::PlanarMeasurement*>(gfTrack.getPointWithMeasurement(
1202 i)->getRawMeasurement(0));
1203 genfit::PlanarMeasurement* planarMeas2 = dynamic_cast<genfit::PlanarMeasurement*>(gfTrack.getPointWithMeasurement(
1204 i + 1)->getRawMeasurement(0));
1205
1206 if (planarMeas1 != NULL && planarMeas2 != NULL &&
1207 planarMeas1->getDetId() == planarMeas2->getDetId() &&
1208 planarMeas1->getPlaneId() != -1 && // -1 is default plane id
1209 planarMeas1->getPlaneId() == planarMeas2->getPlaneId()) {
1210 Belle2::AlignableSVDRecoHit* hit1 = dynamic_cast<Belle2::AlignableSVDRecoHit*>(planarMeas1);
1211 Belle2::AlignableSVDRecoHit* hit2 = dynamic_cast<Belle2::AlignableSVDRecoHit*>(planarMeas2);
1212 if (hit1 && hit2) {
1213 Belle2::AlignableSVDRecoHit* hitU(NULL);
1214 Belle2::AlignableSVDRecoHit* hitV(NULL);
1215 // We have to decide U/V now (else AlignableSVDRecoHit2D could throw FATAL)
1216 if (hit1->isU() && !hit2->isU()) {
1217 hitU = hit1;
1218 hitV = hit2;
1219 } else if (!hit1->isU() && hit2->isU()) {
1220 hitU = hit2;
1221 hitV = hit1;
1222 } else {
1223 continue;
1224 }
1226 // insert measurement before point i (increases number of correct point to i+1)
1227 gfTrack.insertMeasurement(hit, i);
1228 // now delete current point (at its original place, we have the new 2D recohit)
1229 gfTrack.deletePoint(i + 1);
1230 gfTrack.deletePoint(i + 1);
1231 }
1232 }
1233 }
1234 } catch (std::exception& e) {
1235 B2ERROR(e.what());
1236 B2ERROR("SVD Cluster combination failed. This is symptomatic of pruned tracks. MillepedeCollector cannot process pruned tracks.");
1237 return false;
1238 }
1239
1240 try {
1241 gbl->processTrackWithRep(&gfTrack, gfTrack.getCardinalRep(), true);
1242 } catch (genfit::Exception& e) {
1243 B2ERROR(e.what());
1244 return false;
1245 } catch (...) {
1246 B2ERROR("GBL fit failed.");
1247 return false;
1248 }
1249
1250 return true;
1251}
1252
1253std::vector< genfit::Track* > MillepedeCollectorModule::getParticlesTracks(std::vector<Particle*> particles, bool addVertexPoint)
1254{
1255 std::vector< genfit::Track* > tracks;
1256 for (auto particle : particles) {
1257 auto belle2Track = particle->getTrack();
1258 if (!belle2Track) {
1259 B2WARNING("No Belle2::Track for particle (particle->X");
1260 continue;
1261 }
1262// auto trackFitResult = belle2Track->getTrackFitResult(Const::chargedStableSet.find(abs(particle->getPDGCode())));
1263// if (!trackFitResult) {
1264// B2INFO("No track fit result for track");
1265// continue;
1266// }
1267// auto recoTrack = trackFitResult->getRelatedFrom<RecoTrack>();
1268 auto recoTrack = belle2Track->getRelatedTo<RecoTrack>();
1269
1270 if (!recoTrack) {
1271 B2WARNING("No related RecoTrack for Belle2::Track (particle->Track->X)");
1272 continue;
1273 }
1274
1275 // If any track fails, fail completely
1276 if (!fitRecoTrack(*recoTrack, (addVertexPoint) ? particle : nullptr))
1277 return {};
1278
1279 auto& track = RecoTrackGenfitAccess::getGenfitTrack(*recoTrack);
1280
1281 if (!track.hasFitStatus()) {
1282 B2WARNING("Track has no fit status");
1283 continue;
1284 }
1285 genfit::GblFitStatus* fs = dynamic_cast<genfit::GblFitStatus*>(track.getFitStatus());
1286 if (!fs) {
1287 B2WARNING("Track FitStatus is not GblFitStatus.");
1288 continue;
1289 }
1290 if (!fs->isFittedWithReferenceTrack()) {
1291 B2WARNING("Track is not fitted with reference track.");
1292 continue;
1293 }
1294
1295 tracks.push_back(&track);
1296 }
1297
1298 return tracks;
1299}
1300
1302 double motherMass)
1303{
1304 std::vector<TMatrixD> result;
1305
1306 double px = mother.getPx();
1307 double py = mother.getPy();
1308 double pz = mother.getPz();
1309 double pt = sqrt(px * px + py * py);
1310 double p = mother.getMomentumMagnitude();
1311 double M = motherMass;
1312 double m = mother.getDaughter(0)->getPDGMass();
1313
1314 if (mother.getNDaughters() != 2
1315 || m != mother.getDaughter(1)->getPDGMass()) B2FATAL("Only two same-mass daughters (V0->f+f- decays) allowed.");
1316
1317 // Rotation matrix from mother reference system to lab system
1318 TMatrixD mother2lab(3, 3);
1319 mother2lab(0, 0) = px * pz / pt / p; mother2lab(0, 1) = - py / pt; mother2lab(0, 2) = px / p;
1320 mother2lab(1, 0) = py * pz / pt / p; mother2lab(1, 1) = px / pt; mother2lab(1, 2) = py / p;
1321 mother2lab(2, 0) = - pt / p; mother2lab(2, 1) = 0; mother2lab(2, 2) = pz / p;
1322 ROOT::Math::Rotation3D lab2mother;
1323 lab2mother.SetRotationMatrix(mother2lab); lab2mother.Invert();
1324
1325 // Need to rotate and boost daughters' momenta to know which goes forward (+sign in decay model)
1326 // and to get the angles theta, phi of the decaying daughter system in mothers' reference frame
1327 RestFrame boostedFrame(&mother);
1328 ROOT::Math::PxPyPzEVector fourVector1 = mother.getDaughter(0)->get4Vector();
1329 ROOT::Math::PxPyPzEVector fourVector2 = mother.getDaughter(1)->get4Vector();
1330
1331 auto mom1 = lab2mother * boostedFrame.getMomentum(fourVector1).Vect();
1332 auto mom2 = lab2mother * boostedFrame.getMomentum(fourVector2).Vect();
1333 // One momentum has opposite direction (otherwise should be same in CMS of mother), but which?
1334 double sign = 1.;
1335 auto avgMom = 0.5 * (mom1 - mom2);
1336 if (avgMom.Z() < 0.) {
1337 avgMom *= -1.;
1338 // switch meaning of plus/minus trajectories
1339 sign = -1.;
1340 }
1341
1342 double theta = atan2(avgMom.rho(), avgMom.Z());
1343 double phi = atan2(avgMom.Y(), avgMom.X());
1344 if (phi < 0.) phi += 2. * TMath::Pi();
1345
1346 double alpha = M / 2. / m;
1347 double c1 = m * sqrt(alpha * alpha - 1.);
1348 double c2 = 0.5 * sqrt((alpha * alpha - 1.) / alpha / alpha * (p * p + M * M));
1349
1350 double p3 = p * p * p;
1351 double pt3 = pt * pt * pt;
1352
1353
1354 for (auto& track : getParticlesTracks(mother.getDaughters())) {
1355
1356
1357 TMatrixD R = mother2lab;
1358 B2Vector3D P(sign * c1 * sin(theta) * cos(phi),
1359 sign * c1 * sin(theta) * sin(phi),
1360 p / 2. + sign * c2 * cos(theta));
1361
1362 TMatrixD dRdpx(3, 3);
1363 dRdpx(0, 0) = - pz * (pow(px, 4.) - pow(py, 4.) - py * py * pz * pz) / pt3 / p3;
1364 dRdpx(0, 1) = px * py / pt3;
1365 dRdpx(0, 2) = (py * py + pz * pz) / p3;
1366
1367 dRdpx(1, 0) = - px * py * pz * (2. * px * px + 2. * py * py + pz * pz) / pt3 / p3;
1368 dRdpx(1, 1) = - py * py / pt3;
1369 dRdpx(1, 2) = px * py / p3;
1370
1371 dRdpx(2, 0) = - px * pz * pz / pt / p3;
1372 dRdpx(2, 1) = 0.;
1373 dRdpx(2, 2) = - px * pz / p3;
1374
1375 TMatrixD dRdpy(3, 3);
1376 dRdpy(0, 0) = - px * py * pz * (2. * px * px + 2. * py * py + pz * pz) / pt3 / p3;
1377 dRdpy(0, 1) = - px * px / pt3;
1378 dRdpy(0, 2) = px * pz / p3;
1379
1380 dRdpy(1, 0) = - pz * (- pow(px, 4.) - px * px * pz * pz + pow(py, 4.)) / pt3 / p3;
1381 dRdpy(1, 1) = px * py / pt3;
1382 dRdpy(1, 2) = (px * px + pz * pz) / p3;
1383
1384 dRdpy(2, 0) = - py * pz * pz / pt / p3;
1385 dRdpy(2, 1) = 0.;
1386 dRdpy(2, 2) = - py * pz / p3;
1387
1388 TMatrixD dRdpz(3, 3);
1389 dRdpz(0, 0) = px * pt / p3;
1390 dRdpz(0, 1) = 0.;
1391 dRdpz(0, 2) = - px * pz / p3;
1392
1393 dRdpz(1, 0) = py * pt / p3;
1394 dRdpz(1, 1) = 0.;
1395 dRdpz(1, 2) = py * pz / p3;
1396
1397 dRdpz(2, 0) = pz * pt / p3;
1398 dRdpz(2, 1) = 0.;
1399 dRdpz(2, 2) = (px * px + py * py) / p3;
1400
1401 auto K = 1. / 2. / p + sign * cos(theta) * m * m * (M * M / 4. / m / m - 1.) / M / M / sqrt(m * m * (M * M / 4. / m / m - 1.) *
1402 (M * M + p * p) / M / M);
1403
1404 B2Vector3D dpdpx = dRdpx * P + R * K * px * B2Vector3D(0., 0., 1.);
1405 B2Vector3D dpdpy = dRdpy * P + R * K * py * B2Vector3D(0., 0., 1.);
1406 B2Vector3D dpdpz = dRdpz * P + R * K * pz * B2Vector3D(0., 0., 1.);
1407
1408 B2Vector3D dpdtheta = R * B2Vector3D(sign * c1 * cos(theta) * cos(phi),
1409 sign * c1 * cos(theta) * sin(phi),
1410 sign * c2 * (- sin(theta)));
1411
1412
1413 B2Vector3D dpdphi = R * B2Vector3D(sign * c1 * sin(theta) * (- sin(phi)),
1414 sign * c1 * sin(theta) * cos(phi),
1415 0.);
1416
1417 double dc1dM = m * M / (2. * sqrt(M * M - 4. * m * m));
1418 double dc2dM = M * (4. * m * m * p * p + pow(M, 4)) / (2 * M * M * M * sqrt((M * M - 4. * m * m) * (p * p + M * M)));
1419
1420 B2Vector3D dpdM = R * B2Vector3D(sign * sin(theta) * cos(phi) * dc1dM,
1421 sign * sin(theta) * sin(phi) * dc1dM,
1422 sign * cos(theta) * dc2dM);
1423
1424 TMatrixD dpdz(3, 6);
1425 dpdz(0, 0) = dpdpx(0); dpdz(0, 1) = dpdpy(0); dpdz(0, 2) = dpdpz(0); dpdz(0, 3) = dpdtheta(0); dpdz(0, 4) = dpdphi(0);
1426 dpdz(0, 5) = dpdM(0);
1427 dpdz(1, 0) = dpdpx(1); dpdz(1, 1) = dpdpy(1); dpdz(1, 2) = dpdpz(1); dpdz(1, 3) = dpdtheta(1); dpdz(1, 4) = dpdphi(1);
1428 dpdz(1, 5) = dpdM(1);
1429 dpdz(2, 0) = dpdpx(2); dpdz(2, 1) = dpdpy(2); dpdz(2, 2) = dpdpz(2); dpdz(2, 3) = dpdtheta(2); dpdz(2, 4) = dpdphi(2);
1430 dpdz(2, 5) = dpdM(2);
1431
1432 TMatrixD dqdv = getGlobalToLocalTransform(track->getFittedState()).GetSub(0, 4, 0, 2);
1433 TMatrixD dqdp = getGlobalToLocalTransform(track->getFittedState()).GetSub(0, 4, 3, 5);
1434 TMatrixD dfdvz(5, 9);
1435 dfdvz.SetSub(0, 0, dqdv);
1436 dfdvz.SetSub(0, 3, dqdp * dpdz);
1437
1438 result.push_back(dfdvz);
1439
1440 // switch sign for second trajectory
1441 sign *= -1.;
1442 }
1443
1444 return {result[0], result[1]};
1445}
1446
1447TMatrixD MillepedeCollectorModule::getGlobalToLocalTransform(const genfit::MeasuredStateOnPlane& msop)
1448{
1449 auto state = msop;
1450 const B2Vector3D& U(state.getPlane()->getU());
1451 const B2Vector3D& V(state.getPlane()->getV());
1452 const B2Vector3D& O(state.getPlane()->getO());
1453 const B2Vector3D& W(state.getPlane()->getNormal());
1454
1455 const double* state5 = state.getState().GetMatrixArray();
1456
1457 double spu = 1.;
1458
1459 const TVectorD& auxInfo = state.getAuxInfo();
1460 if (auxInfo.GetNrows() == 2
1461 || auxInfo.GetNrows() == 1) // backwards compatibility with old RKTrackRep
1462 spu = state.getAuxInfo()(0);
1463
1464 TVectorD state7(7);
1465
1466 state7[0] = O.X() + state5[3] * U.X() + state5[4] * V.X(); // x
1467 state7[1] = O.Y() + state5[3] * U.Y() + state5[4] * V.Y(); // y
1468 state7[2] = O.Z() + state5[3] * U.Z() + state5[4] * V.Z(); // z
1469
1470 state7[3] = spu * (W.X() + state5[1] * U.X() + state5[2] * V.X()); // a_x
1471 state7[4] = spu * (W.Y() + state5[1] * U.Y() + state5[2] * V.Y()); // a_y
1472 state7[5] = spu * (W.Z() + state5[1] * U.Z() + state5[2] * V.Z()); // a_z
1473
1474 // normalize dir
1475 double norm = 1. / sqrt(state7[3] * state7[3] + state7[4] * state7[4] + state7[5] * state7[5]);
1476 for (unsigned int i = 3; i < 6; ++i) state7[i] *= norm;
1477
1478 state7[6] = state5[0]; // q/p
1479
1480 const double AtU = state7[3] * U.X() + state7[4] * U.Y() + state7[5] * U.Z();
1481 const double AtV = state7[3] * V.X() + state7[4] * V.Y() + state7[5] * V.Z();
1482 const double AtW = state7[3] * W.X() + state7[4] * W.Y() + state7[5] * W.Z();
1483
1484 // J_Mp matrix is d(q/p,u',v',u,v) / d(x,y,z,px,py,pz) (in is 6x6)
1485
1486 const double qop = state7[6];
1487 const double p = state.getCharge() / qop; // momentum
1488
1489 TMatrixD J_Mp_6x5(6, 5);
1490 J_Mp_6x5.Zero();
1491
1492 //d(u)/d(x,y,z)
1493 J_Mp_6x5(0, 3) = U.X(); // [0][3]
1494 J_Mp_6x5(1, 3) = U.Y(); // [1][3]
1495 J_Mp_6x5(2, 3) = U.Z(); // [2][3]
1496 //d(v)/d(x,y,z)
1497 J_Mp_6x5(0, 4) = V.X(); // [0][4]
1498 J_Mp_6x5(1, 4) = V.Y(); // [1][4]
1499 J_Mp_6x5(2, 4) = V.Z(); // [2][4]
1500
1501 // d(q/p)/d(px,py,pz)
1502 double fact = (-1.) * qop / p;
1503 J_Mp_6x5(3, 0) = fact * state7[3]; // [3][0]
1504 J_Mp_6x5(4, 0) = fact * state7[4]; // [4][0]
1505 J_Mp_6x5(5, 0) = fact * state7[5]; // [5][0]
1506 // d(u')/d(px,py,pz)
1507 fact = 1. / (p * AtW * AtW);
1508 J_Mp_6x5(3, 1) = fact * (U.X() * AtW - W.X() * AtU); // [3][1]
1509 J_Mp_6x5(4, 1) = fact * (U.Y() * AtW - W.Y() * AtU); // [4][1]
1510 J_Mp_6x5(5, 1) = fact * (U.Z() * AtW - W.Z() * AtU); // [5][1]
1511 // d(v')/d(px,py,pz)
1512 J_Mp_6x5(3, 2) = fact * (V.X() * AtW - W.X() * AtV); // [3][2]
1513 J_Mp_6x5(4, 2) = fact * (V.Y() * AtW - W.Y() * AtV); // [4][2]
1514 J_Mp_6x5(5, 2) = fact * (V.Z() * AtW - W.Z() * AtV); // [5][2]
1515
1516 return J_Mp_6x5.T();
1517}
1518
1519TMatrixD MillepedeCollectorModule::getLocalToGlobalTransform(const genfit::MeasuredStateOnPlane& msop)
1520{
1521 auto state = msop;
1522 // get vectors and aux variables
1523 const B2Vector3D& U(state.getPlane()->getU());
1524 const B2Vector3D& V(state.getPlane()->getV());
1525 const B2Vector3D& W(state.getPlane()->getNormal());
1526
1527 const TVectorD& state5(state.getState());
1528 double spu = 1.;
1529
1530 const TVectorD& auxInfo = state.getAuxInfo();
1531 if (auxInfo.GetNrows() == 2
1532 || auxInfo.GetNrows() == 1) // backwards compatibility with old RKTrackRep
1533 spu = state.getAuxInfo()(0);
1534
1535 TVectorD pTilde(3);
1536 pTilde[0] = spu * (W.X() + state5(1) * U.X() + state5(2) * V.X()); // a_x
1537 pTilde[1] = spu * (W.Y() + state5(1) * U.Y() + state5(2) * V.Y()); // a_y
1538 pTilde[2] = spu * (W.Z() + state5(1) * U.Z() + state5(2) * V.Z()); // a_z
1539
1540 const double pTildeMag = sqrt(pTilde[0] * pTilde[0] + pTilde[1] * pTilde[1] + pTilde[2] * pTilde[2]);
1541 const double pTildeMag2 = pTildeMag * pTildeMag;
1542
1543 const double utpTildeOverpTildeMag2 = (U.X() * pTilde[0] + U.Y() * pTilde[1] + U.Z() * pTilde[2]) / pTildeMag2;
1544 const double vtpTildeOverpTildeMag2 = (V.X() * pTilde[0] + V.Y() * pTilde[1] + V.Z() * pTilde[2]) / pTildeMag2;
1545
1546 //J_pM matrix is d(x,y,z,px,py,pz) / d(q/p,u',v',u,v) (out is 6x6)
1547
1548 const double qop = state5(0);
1549 const double p = state.getCharge() / qop; // momentum
1550
1551 TMatrixD J_pM_5x6(5, 6);
1552 J_pM_5x6.Zero();
1553
1554 // d(px,py,pz)/d(q/p)
1555 double fact = -1. * p / (pTildeMag * qop);
1556 J_pM_5x6(0, 3) = fact * pTilde[0]; // [0][3]
1557 J_pM_5x6(0, 4) = fact * pTilde[1]; // [0][4]
1558 J_pM_5x6(0, 5) = fact * pTilde[2]; // [0][5]
1559 // d(px,py,pz)/d(u')
1560 fact = p * spu / pTildeMag;
1561 J_pM_5x6(1, 3) = fact * (U.X() - pTilde[0] * utpTildeOverpTildeMag2); // [1][3]
1562 J_pM_5x6(1, 4) = fact * (U.Y() - pTilde[1] * utpTildeOverpTildeMag2); // [1][4]
1563 J_pM_5x6(1, 5) = fact * (U.Z() - pTilde[2] * utpTildeOverpTildeMag2); // [1][5]
1564 // d(px,py,pz)/d(v')
1565 J_pM_5x6(2, 3) = fact * (V.X() - pTilde[0] * vtpTildeOverpTildeMag2); // [2][3]
1566 J_pM_5x6(2, 4) = fact * (V.Y() - pTilde[1] * vtpTildeOverpTildeMag2); // [2][4]
1567 J_pM_5x6(2, 5) = fact * (V.Z() - pTilde[2] * vtpTildeOverpTildeMag2); // [2][5]
1568 // d(x,y,z)/d(u)
1569 J_pM_5x6(3, 0) = U.X(); // [3][0]
1570 J_pM_5x6(3, 1) = U.Y(); // [3][1]
1571 J_pM_5x6(3, 2) = U.Z(); // [3][2]
1572 // d(x,y,z)/d(v)
1573 J_pM_5x6(4, 0) = V.X(); // [4][0]
1574 J_pM_5x6(4, 1) = V.Y(); // [4][1]
1575 J_pM_5x6(4, 2) = V.Z(); // [4][2]
1576
1577 return J_pM_5x6.T();
1578
1579}
1580
1581tuple<B2Vector3D, TMatrixDSym> MillepedeCollectorModule::getPrimaryVertexAndCov() const
1582{
1583 DBObjPtr<BeamSpot> beam;
1584 return {beam->getIPPosition(), beam->getSizeCovMatrix()};
1585}
1586
1587void MillepedeCollectorModule::updateMassWidthIfSet(string listName, double& mass, double& width)
1588{
1589 if (m_customMassConfig.find(listName) != m_customMassConfig.end()) {
1590 auto massWidth = m_customMassConfig.at(listName);
1591 mass = std::get<0>(massWidth);
1592 width = std::get<1>(massWidth);
1593 }
1594}
1595
R E
internal precision of FFTW codelets
double R
typedef autogenerated by FFTW
#define K(x)
macro autogenerated by FFTW
This class is used to transfer CDC information to the track fit and Millepede.
static bool s_enableWireSaggingGlobalDerivative
Static enabling(true) or disabling(false) addition of global derivative for wire sagging coefficient ...
static bool s_enableWireByWireAlignmentGlobalDerivatives
Static enabling(true) or disabling(false) addition of global derivatives for wire-by-wire alignment.
static bool s_enableTrackT0LocalDerivative
Static enabling(true) or disabling(false) addition of local derivative for track T0.
This class is used to transfer PXD information to the track fit.
This class is used to transfer SVD information to the track fit.
This class is used to transfer SVD information to the track fit.
DataType Z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:435
DataType X() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:431
DataType Y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:433
Calibration collector module base class.
static const ChargedStable muon
muon particle
Definition: Const.h:660
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
@ c_Persistent
Object is available during entire execution time.
Definition: DataStore.h:60
Store event, run, and experiment numbers.
Definition: EventMetaData.h:33
static EvtGenDatabasePDG * Instance()
Instance method that loads the EvtGen table.
static FileCatalog & Instance()
Static method to get a reference to the FileCatalog instance.
Definition: FileCatalog.cc:23
virtual bool registerFile(const std::string &fileName, FileMetaData &metaData, const std::string &oldLFN="")
Register a file in the (local) file catalog.
Definition: FileCatalog.cc:90
Metadata information about a file.
Definition: FileMetaData.h:29
void setLfn(const std::string &lfn)
Setter for LFN.
Definition: FileMetaData.h:139
void setParents(const std::vector< std::string > &parents)
Parents setter.
Definition: FileMetaData.h:173
TrackSegmentController for use with GblFitter in Belle2.
Class to convert to/from global labels for Millepede II to/from detector & parameter identificators.
Definition: GlobalLabel.h:41
Algorithm class to translate the added detector hits (e.g.
void resetMeasurementCreators(const std::vector< std::shared_ptr< PXDBaseMeasurementCreator > > &pxdMeasurementCreators, const std::vector< std::shared_ptr< SVDBaseMeasurementCreator > > &svdMeasurementCreators, const std::vector< std::shared_ptr< CDCBaseMeasurementCreator > > &cdcMeasurementCreators, const std::vector< std::shared_ptr< BKLMBaseMeasurementCreator > > &bklmMeasurementCreators, const std::vector< std::shared_ptr< EKLMBaseMeasurementCreator > > &eklmMeasurementCreators, const std::vector< std::shared_ptr< BaseMeasurementCreator > > &additionalMeasurementCreators)
If you want to use non-default settings for the store arrays, you can create your own instances of th...
bool addMeasurements(RecoTrack &recoTrack) const
After you have filled the internal storage with measurement creators (either by providing your own or...
Mergeable class holding list of so far opened mille binaries and providing the binaries.
Definition: MilleData.h:24
bool m_updateCDCWeights
Update L/R weights from previous DAF fit result?
std::vector< std::string > m_twoBodyDecays
Name of particle list with mothers of daughters to be used with vertex + mass constraint in calibrati...
TMatrixD getLocalToGlobalTransform(const genfit::MeasuredStateOnPlane &msop)
Compute the transformation matrix d(x,y,z,px,py,pz)/d(q/p,u',v',u,v) from state at first track point ...
std::vector< std::string > m_tracks
Names of arrays with single RecoTracks fitted by GBL.
MillepedeCollectorModule()
Constructor: Sets the description, the properties and the parameters of the module.
StoreObjPtr< EventT0 > m_eventT0
Optional input for EventT0.
std::vector< std::string > m_components
Whether to use VXD alignment hierarchy.
double m_minCDCHitWeight
Minimum CDC hit weight.
std::vector< std::string > m_primaryMassTwoBodyDecays
Name of particle list with mothers of daughters to be used with vertex + IP profile + mass constraint...
double m_minPValue
Minimum p.value for output.
std::string getUniqueMilleName()
Make a name for mille binary (encodes module name + starting exp, run and event + process id)
std::vector< std::tuple< int, int, int > > m_eventNumbers
List of event meta data entries at which payloads can change for timedep calibration.
bool m_absFilePaths
Use absolute path to locate binary files in MilleData.
std::vector< std::string > m_vertices
Name of particle list with mothers of daughters to be used with vertex constraint in calibration.
bool m_fitTrackT0
Add local parameter for track T0 fit in GBL (local derivative)
std::vector< genfit::Track * > getParticlesTracks(std::vector< Particle * > particles, bool addVertexPoint=true)
Get all usable tracks for particles.
std::vector< gbl::GblData > m_currentGblData
Current vector of GBL data from trajectory to be stored in a tree.
bool fitRecoTrack(RecoTrack &recoTrack, Particle *particle=nullptr)
Fit given RecoTrack with GBL.
bool m_enableWireSagging
Enable global derivatives for wire sagging.
TMatrixD getGlobalToLocalTransform(const genfit::MeasuredStateOnPlane &msop)
Compute the transformation matrix d(q/p,u',v',u,v)/d(x,y,z,px,py,pz) from state at first track point ...
std::tuple< B2Vector3D, TMatrixDSym > getPrimaryVertexAndCov() const
Get the primary vertex position estimation and its size from BeamSpot.
std::map< std::string, std::tuple< double, double > > m_customMassConfig
Map of list_name -> (mass, width) for custom mass and width setting.
int m_recalcJacobians
Up to which external iteration propagation Jacobians should be re-calculated.
bool m_useGblTree
Whether to use TTree to accumulate GBL data instead of binary files.
void storeTrajectory(gbl::GblTrajectory &trajectory)
Write down a GBL trajectory (to TTree or binary file)
bool m_doublePrecision
Use double (instead of single/float) precision for binary files.
virtual void collect() override
Data collection.
StoreObjPtr< EventMetaData > m_evtMetaData
Required object pointer to EventMetaData.
std::pair< TMatrixD, TMatrixD > getTwoBodyToLocalTransform(Particle &mother, double motherMass)
Compute the transformation matrices d(q/p,u'v',u,v)/d(vx,vy,vz,px,py,pz,theta,phi,...
bool m_enablePXDHierarchy
enable PXD hierarchy
virtual void closeRun() override
Only for closing mille binaries after each run.
bool m_calibrateKinematics
Add derivatives for beam spot kinematics calibration for primary vertices.
double m_minUsedCDCHitFraction
Minimum CDC used hit fraction.
void updateMassWidthIfSet(std::string listName, double &mass, double &width)
Update mass and width of the particle (mother in list) with user custom-defined values.
virtual void prepare() override
Prepration.
bool m_enableSVDHierarchy
enable SVD hierarchy
std::string m_internalIterations
String defining internal GBL iterations for outlier down-weighting.
std::vector< std::tuple< std::vector< int >, std::vector< std::tuple< int, int, int > > > > m_timedepConfig
Config for time dependence: list( tuple( list( param1, param2, ... ), list( (ev, run,...
std::vector< std::string > m_particles
Names of particle list with single particles.
virtual void finish() override
Register mille binaries in file catalog.
std::vector< std::string > m_primaryVertices
Name of particle list with mothers of daughters to be used with vertex + IP profile (+ optional calib...
int m_externalIterations
Number of external iterations of GBL fitter.
int m_hierarchyType
Type of alignment hierarchy (for VXD only for now): 0 = None, 1 = Flat (only constraints,...
bool m_enableWireByWireAlignment
Enable global derivatives for wire-by-wire alignment.
double m_stableParticleWidth
Width (in GeV/c/c) to use for invariant mass constraint for 'stable' particles (like K short).
std::vector< std::string > m_primaryTwoBodyDecays
Name of particle list with mothers of daughters to be used with vertex + IP profile (+ optional calib...
std::vector< std::string > m_primaryMassVertexTwoBodyDecays
Name of particle list with mothers of daughters to be used with vertex + IP profile + mass constraint...
bool m_calibrateVertex
Add derivatives for beam spot vertex calibration for primary vertices.
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
const std::string & getName() const
Returns the name of the module.
Definition: Module.h:186
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
Class to store reconstructed particles.
Definition: Particle.h:76
double getPx() const
Returns x component of momentum.
Definition: Particle.h:607
double getPz() const
Returns z component of momentum.
Definition: Particle.h:625
double getPy() const
Returns y component of momentum.
Definition: Particle.h:616
unsigned getNDaughters(void) const
Returns number of daughter particles.
Definition: Particle.h:747
std::vector< Belle2::Particle * > getDaughters() const
Returns a vector of pointers to daughter particles.
Definition: Particle.cc:668
double getPDGMass(void) const
Returns uncertainty on the invariant mass (requires valid momentum error matrix)
Definition: Particle.cc:635
ROOT::Math::PxPyPzEVector get4Vector() const
Returns Lorentz vector.
Definition: Particle.h:567
double getMomentumMagnitude() const
Returns momentum magnitude.
Definition: Particle.h:589
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
Definition: Particle.cc:662
static int EvtProcID()
Return ID of the current process.
Definition: ProcHandler.cc:248
static bool parallelProcessingUsed()
Returns true if multiple processes have been spawned, false in single-core mode.
Definition: ProcHandler.cc:226
This class stores additional information to every CDC/SVD/PXD hit stored in a RecoTrack.
CDCHit UsedCDCHit
Define, use of CDC hits as CDC hits (for symmetry).
EKLMAlignmentHit UsedEKLMHit
Define, use of EKLMHit2d as EKLM hits.
KLMHit2d UsedBKLMHit
Define, use of KLMHit2d as BKLM hits.
PXDCluster UsedPXDHit
Define, use of clusters or true hits for PXD.
SVDCluster UsedSVDHit
Define, use of clusters or true hits for SVD.
static genfit::Track & getGenfitTrack(RecoTrack &recoTrack)
Give access to the RecoTrack's genfit::Track.
Definition: RecoTrack.cc:404
static genfit::AbsTrackRep * createOrReturnRKTrackRep(RecoTrack &recoTrack, int PDGcode)
Checks if a TrackRap for the PDG id of the RecoTrack (and its charge conjugate) does already exit and...
Definition: RecoTrack.cc:409
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
genfit::AbsTrackRep * getCardinalRepresentation() const
Get a pointer to the cardinal track representation. You are not allowed to modify or delete it!
Definition: RecoTrack.h:631
unsigned int getNumberOfCDCHits() const
Return the number of cdc hits.
Definition: RecoTrack.h:427
const std::string & getStoreArrayNameOfRecoHitInformation() const
Name of the store array of the reco hit information.
Definition: RecoTrack.h:747
const genfit::TrackPoint * getCreatedTrackPoint(const RecoHitInformation *recoHitInformation) const
Get a pointer to the TrackPoint that was created from this hit.
Definition: RecoTrack.cc:230
const genfit::FitStatus * getTrackFitStatus(const genfit::AbsTrackRep *representation=nullptr) const
Return the track fit status for the given representation or for the cardinal one. You are not allowed...
Definition: RecoTrack.h:621
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
Rest frame of a particle.
virtual ROOT::Math::PxPyPzEVector getMomentum(const ROOT::Math::PxPyPzEVector &vector) const override
Get Lorentz vector in rest frame System.
bool isU() const
Is the coordinate u or v?
Definition: SVDRecoHit.h:91
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
bool isValid() const
Check whether the array was registered.
Definition: StoreArray.h:288
TClonesArray * getPtr() const
Raw access to the underlying TClonesArray.
Definition: StoreArray.h:311
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:95
bool isValid() const
Check whether the object was created.
Definition: StoreObjPtr.h:110
static int createCorrectPDGCodeForChargedStable(const Const::ChargedStable &particleType, const RecoTrack &recoTrack)
Helper function to multiply the PDG code of a charged stable with the charge of the reco track (if ne...
Definition: TrackFitter.cc:24
static const double GeV
Standard of [energy, momentum, mass].
Definition: Unit.h:51
void writeConstraints(std::string txtFilename)
Write-out complete hierarchy to a text file.
Definition: Manager.cc:169
void initialize(const std::vector< std::string > &components={}, const std::vector< EventMetaData > &timeSlices={})
Initialize the manager with given configuration (from MillepedeCollector)
Definition: Manager.cc:52
void preCollect(const EventMetaData &emd)
Notice manager of a coming event (from MillepedeCollector)
Definition: Manager.cc:98
static GlobalCalibrationManager & getInstance()
Get instance of the Manager auto& gcm = GlobalCalibrationManager::getInstance();.
Definition: Manager.cc:27
Class for easier manipulation with global derivatives (and their labels)
static bool s_enablePXD
Enable PXD in hierarchy?
Definition: GlobalParam.h:85
static bool s_enableSVD
Enable SVD in hierarchy?
Definition: GlobalParam.h:87
static E_VXDHierarchyType s_hierarchyType
What type of hierarchy to use for VXD?
Definition: GlobalParam.h:83
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:649
B2Vector3< double > B2Vector3D
typedef for common usage with double
Definition: B2Vector3.h:516
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
CoordinateMeasurementCreator< RecoHitInformation::UsedSVDHit, Const::SVD > SVDCoordinateMeasurementCreator
Hit to reco hit measurement creator for the SVD.
CoordinateMeasurementCreator< RecoHitInformation::UsedPXDHit, Const::PXD > PXDCoordinateMeasurementCreator
Hit to reco hit measurement creator for the PXD.
CoordinateMeasurementCreator< RecoHitInformation::UsedBKLMHit, Const::BKLM > BKLMCoordinateMeasurementCreator
Hit to reco hit measurement creator for the BKLM.
CoordinateMeasurementCreator< RecoHitInformation::UsedCDCHit, Const::CDC > CDCCoordinateMeasurementCreator
Needed for templating.
CoordinateMeasurementCreator< RecoHitInformation::UsedEKLMHit, Const::EKLM > EKLMCoordinateMeasurementCreator
Hit to reco hit measurement creator for the EKLM.
Abstract base class for different kinds of events.
STL namespace.