Belle II Software  release-08-01-10
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 
49 using namespace std;
50 using namespace Belle2;
51 using namespace alignment;
52 
53 //-----------------------------------------------------------------
54 // Register the Module
55 //-----------------------------------------------------------------
56 REG_MODULE(MillepedeCollector);
57 
58 //-----------------------------------------------------------------
59 // Implementation
60 //-----------------------------------------------------------------
61 
62 MillepedeCollectorModule::MillepedeCollectorModule() : CalibrationCollectorModule()
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.",
102  m_components);
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() &&
168  m_primaryMassTwoBodyDecays.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  bool first(true);
384  for (auto& track : getParticlesTracks(mother->getDaughters())) {
385  if (first) {
386  // For first trajectory only
387  extProjection = getGlobalToLocalTransform(track->getFittedState()).GetSub(0, 4, 0, 2);
388  locProjection = getLocalToGlobalTransform(track->getFittedState()).GetSub(0, 2, 0, 4);
389  first = false;
390  }
391  daughters.push_back({
392  gbl->collectGblPoints(track, track->getCardinalRep()),
393  getGlobalToLocalTransform(track->getFittedState()).GetSub(0, 4, 0, 2)
394  });
395  }
396 
397  if (daughters.size() > 1) {
398  auto beam = getPrimaryVertexAndCov();
399 
400  TMatrixDSym vertexCov(get<TMatrixDSym>(beam));
401  TMatrixDSym vertexPrec(get<TMatrixDSym>(beam).Invert());
402  B2Vector3D vertexResidual = - (B2Vector3D(mother->getVertex()) - get<B2Vector3D>(beam));
403 
404  TVectorD extMeasurements(3);
405  extMeasurements[0] = vertexResidual[0];
406  extMeasurements[1] = vertexResidual[1];
407  extMeasurements[2] = vertexResidual[2];
408 
409  TMatrixD extDeriv(3, 3);
410  extDeriv.Zero();
411  // beam vertex constraint
412  extDeriv(0, 0) = 1.;
413  extDeriv(1, 1) = 1.;
414  extDeriv(2, 2) = 1.;
415 
416  if (m_calibrateVertex) {
417  TMatrixD derivatives(3, 3);
418  derivatives.Zero();
419  derivatives(0, 0) = 1.;
420  derivatives(1, 1) = 1.;
421  derivatives(2, 2) = 1.;
422 
423  std::vector<int> labels;
424  GlobalLabel label = GlobalLabel::construct<BeamSpot>(0, 0);
425  labels.push_back(label.setParameterId(1));
426  labels.push_back(label.setParameterId(2));
427  labels.push_back(label.setParameterId(3));
428 
429  // Allow to disable BeamSpot externally
430  alignment::GlobalDerivatives globals(labels, derivatives);
431  // Add derivatives for vertex calibration to first point of first trajectory
432  // NOTE: use GlobalDerivatives operators vector<int> and TMatrixD which filter
433  // the derivatives to not pass those with zero labels (usefull to get rid of some params)
434  std::vector<int> lab(globals); TMatrixD der(globals);
435 
436  // Transformation from local system at (vertex) point to global (vx,vy,vz)
437  // of the (decay) vertex
438  //
439  // d(q/p,u',v',u,v)/d(vy,vy,vz) = dLocal_dExt
440  //
441  //
442  // Note its transpose is its "inverse" in the sense that
443  //
444  // dloc/dext * (dloc/dext)^T = diag(0, 0, 0, 0, 1, 1)
445  //
446  //
447  // N.B. typical dLocal_dExt matrix (5x3):
448  //
449  // | 0 | 1 | 2 |
450  // --------------------------------------------
451  // 0 | 0 0 0
452  // 1 | 0 0 0
453  // 2 | 0 0 0
454  // 3 | -0.02614 -0.9997 0
455  // 4 | 0 0 1
456  //
457  // Therefore one can simplify things by only taking the last two rows/columns in vectors/matrices
458  // and vertex measurement can be expressed as standard 2D measurement in GBL.
459  //
460  TMatrixD dLocal_dExt = extProjection;
461  TMatrixD dExt_dLocal = locProjection;
462 
463  TVectorD locRes = dLocal_dExt * extMeasurements;
464  // Do not use inverted covariance - seems to have issues with numeric precision
465  TMatrixD locCov = dLocal_dExt * vertexCov * dExt_dLocal;
466  // Invert here only the 2D sub-matrix (rest is zero due to the foÅ•m of dLocal_dExt)
467  TMatrixD locPrec = locCov.GetSub(3, 4, 3, 4).Invert();
468  TMatrixDSym locPrec2D(2); locPrec2D.Zero();
469  for (int i = 0; i < 2; ++i)
470  for (int j = 0; j < 2; ++j)
471  locPrec2D(i, j) = locPrec(i, j);
472 
473  // Take the 2 last components also for residuals and global derivatives
474  // (in local system of vertex point - defined during fitRecoTrack(..., particle) and using
475  // the (hopefully) updated momentum and position seed after vertex fit by modularAnalysis
476  TVectorD locRes2D = locRes.GetSub(3, 4);
477  TMatrixD locDerivs2D = (extProjection * der).GetSub(3, 4, 0, 2);
478 
479  // Attach the primary beamspot vertex position as a measurement at 1st point
480  // of first trajectory (and optionaly also the global derivatives for beamspot alignment
481  daughters[0].first[0].addMeasurement(locRes2D, locPrec2D);
482  if (!lab.empty()) {
483  daughters[0].first[0].addGlobals(lab, locDerivs2D);
484  }
485 
486  gbl::GblTrajectory combined(daughters);
487  //combined.printTrajectory(100);
488  //combined.printPoints(100);
489 
490  combined.fit(chi2, ndf, lostWeight);
491  getObjectPtr<TH1I>("ndf")->Fill(ndf);
492  getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
493  getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
494  if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
495  evt0 = m_eventT0->getEventT0();
496  getObjectPtr<TH1F>("evt0")->Fill(evt0);
497  }
498 
499  if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(combined);
500  B2RESULT("Beam vertex constrained fit results NDF = " << ndf << " Chi2/NDF = " << chi2 / double(ndf));
501 
502  } else {
503 
504  gbl::GblTrajectory combined(daughters, extDeriv, extMeasurements, vertexPrec);
505 
506  combined.fit(chi2, ndf, lostWeight);
507  getObjectPtr<TH1I>("ndf")->Fill(ndf);
508  getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
509  getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
510  if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
511  evt0 = m_eventT0->getEventT0();
512  getObjectPtr<TH1F>("evt0")->Fill(evt0);
513  }
514 
515  if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(combined);
516 
517  B2RESULT("Beam vertex constrained fit results NDF = " << ndf << " Chi2/NDF = " << chi2 / double(ndf));
518 
519  }
520  }
521  }
522  }
523 
524  for (auto listName : m_twoBodyDecays) {
525  StoreObjPtr<ParticleList> list(listName);
526  if (!list.isValid())
527  continue;
528 
529  for (unsigned int iParticle = 0; iParticle < list->getListSize(); ++iParticle) {
530 
531  auto mother = list->getParticle(iParticle);
532  auto track12 = getParticlesTracks(mother->getDaughters());
533  if (track12.size() != 2) {
534  B2ERROR("Did not get 2 fitted tracks. Skipping this mother.");
535  continue;
536  }
537 
538  auto pdgdb = EvtGenDatabasePDG::Instance();
539  double motherMass = mother->getPDGMass();
540  double motherWidth = pdgdb->GetParticle(mother->getPDGCode())->Width();
541 
542  updateMassWidthIfSet(listName, motherMass, motherWidth);
543 
544  //TODO: what to take as width for "real" particles? -> make a param for default detector mass resolution??
545  if (motherWidth == 0.) {
546  motherWidth = m_stableParticleWidth * Unit::GeV;
547  B2WARNING("Using artificial width for " << pdgdb->GetParticle(mother->getPDGCode())->GetName() << " : " << motherWidth << " GeV");
548  }
549 
550  auto dfdextPlusMinus = getTwoBodyToLocalTransform(*mother, motherMass);
551  std::vector<std::pair<std::vector<gbl::GblPoint>, TMatrixD> > daughters;
552 
553  daughters.push_back({gbl->collectGblPoints(track12[0], track12[0]->getCardinalRep()), dfdextPlusMinus.first});
554  daughters.push_back({gbl->collectGblPoints(track12[1], track12[1]->getCardinalRep()), dfdextPlusMinus.second});
555 
556  TMatrixDSym massPrec(1); massPrec(0, 0) = 1. / motherWidth / motherWidth;
557  TVectorD massResidual(1); massResidual = - (mother->getMass() - motherMass);
558 
559  TVectorD extMeasurements(1);
560  extMeasurements[0] = massResidual[0];
561 
562  TMatrixD extDeriv(1, 9);
563  extDeriv.Zero();
564  extDeriv(0, 8) = 1.;
565 
566  gbl::GblTrajectory combined(daughters, extDeriv, extMeasurements, massPrec);
567 
568  combined.fit(chi2, ndf, lostWeight);
569  //combined.printTrajectory(1000);
570  //combined.printPoints(1000);
571  getObjectPtr<TH1I>("ndf")->Fill(ndf);
572  getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
573  getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
574  if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
575  evt0 = m_eventT0->getEventT0();
576  getObjectPtr<TH1F>("evt0")->Fill(evt0);
577  }
578 
579 
580  B2RESULT("Mass(PDG) + vertex constrained fit results NDF = " << ndf << " Chi2/NDF = " << chi2 / double(ndf));
581 
582  if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(combined);
583 
584  }
585  }
586 
587  for (auto listName : m_primaryMassTwoBodyDecays) {
588  StoreObjPtr<ParticleList> list(listName);
589  if (!list.isValid())
590  continue;
591 
593 
594  double motherMass = beam->getMass();
595  double motherWidth = sqrt((beam->getCovHER() + beam->getCovLER())(0, 0));
596 
597  updateMassWidthIfSet(listName, motherMass, motherWidth);
598 
599  for (unsigned int iParticle = 0; iParticle < list->getListSize(); ++iParticle) {
600 
601  auto mother = list->getParticle(iParticle);
602  auto track12 = getParticlesTracks(mother->getDaughters());
603  if (track12.size() != 2) {
604  B2ERROR("Did not get 2 fitted tracks. Skipping this mother.");
605  continue;
606  }
607 
608  auto dfdextPlusMinus = getTwoBodyToLocalTransform(*mother, motherMass);
609  std::vector<std::pair<std::vector<gbl::GblPoint>, TMatrixD> > daughters;
610 
611  daughters.push_back({gbl->collectGblPoints(track12[0], track12[0]->getCardinalRep()), dfdextPlusMinus.first});
612  daughters.push_back({gbl->collectGblPoints(track12[1], track12[1]->getCardinalRep()), dfdextPlusMinus.second});
613 
614  TMatrixDSym massPrec(1); massPrec(0, 0) = 1. / motherWidth / motherWidth;
615  TVectorD massResidual(1); massResidual = - (mother->getMass() - motherMass);
616 
617  TVectorD extMeasurements(1);
618  extMeasurements[0] = massResidual[0];
619 
620  TMatrixD extDeriv(1, 9);
621  extDeriv.Zero();
622  extDeriv(0, 8) = 1.;
623 
624  gbl::GblTrajectory combined(daughters, extDeriv, extMeasurements, massPrec);
625 
626  combined.fit(chi2, ndf, lostWeight);
627  getObjectPtr<TH1I>("ndf")->Fill(ndf);
628  getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
629  getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
630  if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
631  evt0 = m_eventT0->getEventT0();
632  getObjectPtr<TH1F>("evt0")->Fill(evt0);
633  }
634 
635 
636  B2RESULT("Mass constrained fit results NDF = " << ndf << " Chi2/NDF = " << chi2 / double(ndf));
637 
638  if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(combined);
639 
640  }
641  }
642 
643  for (auto listName : m_primaryMassVertexTwoBodyDecays) {
644  StoreObjPtr<ParticleList> list(listName);
645  if (!list.isValid())
646  continue;
647 
649 
650  double motherMass = beam->getMass();
651  double motherWidth = sqrt((beam->getCovHER() + beam->getCovLER())(0, 0));
652 
653  updateMassWidthIfSet(listName, motherMass, motherWidth);
654 
655  for (unsigned int iParticle = 0; iParticle < list->getListSize(); ++iParticle) {
656 
657  auto mother = list->getParticle(iParticle);
658  auto track12 = getParticlesTracks(mother->getDaughters());
659  if (track12.size() != 2) {
660  B2ERROR("Did not get 2 fitted tracks. Skipping this mother.");
661  continue;
662  }
663 
664  auto dfdextPlusMinus = getTwoBodyToLocalTransform(*mother, motherMass);
665  std::vector<std::pair<std::vector<gbl::GblPoint>, TMatrixD> > daughters;
666 
667  daughters.push_back({gbl->collectGblPoints(track12[0], track12[0]->getCardinalRep()), dfdextPlusMinus.first});
668  daughters.push_back({gbl->collectGblPoints(track12[1], track12[1]->getCardinalRep()), dfdextPlusMinus.second});
669 
670  TMatrixDSym vertexPrec(get<TMatrixDSym>(getPrimaryVertexAndCov()).Invert());
671  B2Vector3D vertexResidual = - (B2Vector3D(mother->getVertex()) - get<B2Vector3D>(getPrimaryVertexAndCov()));
672 
673  TMatrixDSym massPrec(1); massPrec(0, 0) = 1. / motherWidth / motherWidth;
674  TVectorD massResidual(1); massResidual = - (mother->getMass() - motherMass);
675 
676  TMatrixDSym extPrec(4); extPrec.Zero();
677  extPrec.SetSub(0, 0, vertexPrec);
678  extPrec(3, 3) = massPrec(0, 0);
679 
680  TVectorD extMeasurements(4);
681  extMeasurements[0] = vertexResidual[0];
682  extMeasurements[1] = vertexResidual[1];
683  extMeasurements[2] = vertexResidual[2];
684  extMeasurements[3] = massResidual[0];
685 
686  TMatrixD extDeriv(4, 9);
687  extDeriv.Zero();
688  extDeriv(0, 0) = 1.;
689  extDeriv(1, 1) = 1.;
690  extDeriv(2, 2) = 1.;
691  extDeriv(3, 8) = 1.;
692 
693  gbl::GblTrajectory combined(daughters, extDeriv, extMeasurements, extPrec);
694 
695  combined.fit(chi2, ndf, lostWeight);
696  getObjectPtr<TH1I>("ndf")->Fill(ndf);
697  getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
698  getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
699  if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
700  evt0 = m_eventT0->getEventT0();
701  getObjectPtr<TH1F>("evt0")->Fill(evt0);
702  }
703 
704 
705  if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(combined);
706 
707  B2RESULT("Mass + vertex constrained fit results NDF = " << ndf << " Chi2/NDF = " << chi2 / double(ndf));
708 
709  }
710  }
711 
712  for (auto listName : m_primaryTwoBodyDecays) {
713  B2WARNING("This should NOT be used for production of calibration constants for the real detector (yet)!");
714 
715  StoreObjPtr<ParticleList> list(listName);
716  if (!list.isValid())
717  continue;
718 
720 
721  // For the error of invariant mass M = 2 * sqrt(E_HER * E_LER) (for m_e ~ 0)
722  double M = beam->getMass();
723  double E_HER = beam->getHER().E();
724  double E_LER = beam->getLER().E();
725 
726  double pz = beam->getHER().Pz() + beam->getLER().Pz();
727  double E = (beam->getHER() + beam->getLER()).E();
728 
729  double motherMass = beam->getMass();
730  double motherWidth = sqrt((E_HER / M) * (E_HER / M) * beam->getCovLER()(0, 0) + (E_LER / M) * (E_LER / M) * beam->getCovHER()(0,
731  0));
732 
733  updateMassWidthIfSet(listName, motherMass, motherWidth);
734 
735  for (unsigned int iParticle = 0; iParticle < list->getListSize(); ++iParticle) {
736 
737  B2WARNING("Two body decays with full kinematic constraint not yet correct - need to resolve strange covariance provided by BeamParameters!");
738 
739  auto mother = list->getParticle(iParticle);
740 
741  auto track12 = getParticlesTracks(mother->getDaughters());
742  if (track12.size() != 2) {
743  B2ERROR("Did not get exactly 2 fitted tracks. Skipping this mother in list " << listName);
744  continue;
745  }
746 
747  auto dfdextPlusMinus = getTwoBodyToLocalTransform(*mother, motherMass);
748  std::vector<std::pair<std::vector<gbl::GblPoint>, TMatrixD> > daughters;
749 
750  daughters.push_back({gbl->collectGblPoints(track12[0], track12[0]->getCardinalRep()), dfdextPlusMinus.first});
751  daughters.push_back({gbl->collectGblPoints(track12[1], track12[1]->getCardinalRep()), dfdextPlusMinus.second});
752 
753  TMatrixDSym extCov(7); extCov.Zero();
754 
755  // 3x3 IP vertex covariance
756  extCov.SetSub(0, 0, get<TMatrixDSym>(getPrimaryVertexAndCov()));
757 
758  // 3x3 boost vector covariance
759  //NOTE: BeamSpot return covarince in variables (E, theta_x, theta_y)
760  // We need to transform it to our variables (px, py, pz)
761 
762  TMatrixD dBoost_dVect(3, 3);
763  dBoost_dVect(0, 0) = 0.; dBoost_dVect(0, 1) = 1. / pz; dBoost_dVect(0, 2) = 0.;
764  dBoost_dVect(1, 0) = 0.; dBoost_dVect(1, 1) = 0.; dBoost_dVect(1, 2) = 1. / pz;
765  dBoost_dVect(2, 0) = pz / E; dBoost_dVect(2, 1) = 0.; dBoost_dVect(2, 2) = 0.;
766 
767  TMatrixD dVect_dBoost(3, 3);
768  dVect_dBoost(0, 0) = 0.; dVect_dBoost(0, 1) = 0.; dVect_dBoost(0, 2) = E / pz;
769  dVect_dBoost(1, 0) = pz; dVect_dBoost(1, 1) = 0.; dVect_dBoost(1, 2) = 0.;
770  dVect_dBoost(2, 0) = 0.; dVect_dBoost(2, 1) = pz; dVect_dBoost(2, 2) = 0.;
771 
772  TMatrixD covBoost(3, 3);
773  for (int i = 0; i < 3; ++i) {
774  for (int j = i; j < 3; ++j) {
775  covBoost(j, i) = covBoost(i, j) = (beam->getCovHER() + beam->getCovLER())(i, j);
776  }
777  }
778  //TODO: Temporary fix: if theta_x, theta_y covariance is zero, use arbitrary 10mrad^2
779 // if (covBoost(1, 1) == 0.) covBoost(1, 1) = 1.;
780 // if (covBoost(2, 2) == 0.) covBoost(2, 2) = 1.;
781  if (covBoost(1, 1) == 0.) covBoost(1, 1) = 1.e-4;
782  if (covBoost(2, 2) == 0.) covBoost(2, 2) = 1.e-4;
783 
784  TMatrixD covVect = dBoost_dVect * covBoost * dVect_dBoost;
785 
786  extCov.SetSub(3, 3, covVect);
787 
788  extCov(6, 6) = motherWidth * motherWidth;
789  auto extPrec = extCov; extPrec.Invert();
790 
791  TVectorD extMeasurements(7);
792  extMeasurements[0] = - (B2Vector3D(mother->getVertex()) - get<B2Vector3D>(getPrimaryVertexAndCov()))[0];
793  extMeasurements[1] = - (B2Vector3D(mother->getVertex()) - get<B2Vector3D>(getPrimaryVertexAndCov()))[1];
794  extMeasurements[2] = - (B2Vector3D(mother->getVertex()) - get<B2Vector3D>(getPrimaryVertexAndCov()))[2];
795  extMeasurements[3] = - (B2Vector3D(mother->getMomentum()) - (beam->getHER().Vect() + beam->getLER().Vect()))[0];
796  extMeasurements[4] = - (B2Vector3D(mother->getMomentum()) - (beam->getHER().Vect() + beam->getLER().Vect()))[1];
797  extMeasurements[5] = - (B2Vector3D(mother->getMomentum()) - (beam->getHER().Vect() + beam->getLER().Vect()))[2];
798  extMeasurements[6] = - (mother->getMass() - motherMass);
799 
800  B2INFO("mother mass = " << mother->getMass() << " and beam mass = " << beam->getMass());
801 
802  TMatrixD extDeriv(7, 9);
803  extDeriv.Zero();
804  // beam vertex constraint
805  extDeriv(0, 0) = 1.;
806  extDeriv(1, 1) = 1.;
807  extDeriv(2, 2) = 1.;
808  // beam kinematics constraint
809  extDeriv(3, 3) = 1.;
810  extDeriv(4, 4) = 1.;
811  extDeriv(5, 5) = 1.;
812  // beam inv. mass constraint
813  extDeriv(6, 8) = 1;
814 
816  B2WARNING("Primary vertex+kinematics calibration not (yet?) fully implemented!");
817  B2WARNING("This code is highly experimental and has (un)known issues!");
818 
819  // up to d(x,y,z,px,py,pz,theta,phi,M)/d(vx,vy,vz,theta_x,theta_y,E)
820  TMatrixD derivatives(9, 6);
821  std::vector<int> labels;
822  derivatives.Zero();
823 
824  if (m_calibrateVertex) {
825  derivatives(0, 0) = 1.;
826  derivatives(1, 1) = 1.;
827  derivatives(2, 2) = 1.;
828  GlobalLabel label = GlobalLabel::construct<BeamSpot>(0, 0);
829  labels.push_back(label.setParameterId(1));
830  labels.push_back(label.setParameterId(2));
831  labels.push_back(label.setParameterId(3));
832  } else {
833  labels.push_back(0);
834  labels.push_back(0);
835  labels.push_back(0);
836  }
837 
838  if (m_calibrateKinematics) {
839  derivatives(3, 3) = mother->getMomentumMagnitude();
840  derivatives(4, 4) = mother->getMomentumMagnitude();
841  derivatives(8, 5) = (beam->getLER().E() + beam->getHER().E()) / beam->getMass();
842 
843  GlobalLabel label = GlobalLabel::construct<BeamSpot>(0, 0);
844  labels.push_back(label.setParameterId(4)); //theta_x
845  labels.push_back(label.setParameterId(5)); //theta_y
846  labels.push_back(label.setParameterId(6)); //E
847 
848  } else {
849  labels.push_back(0);
850  labels.push_back(0);
851  labels.push_back(0);
852  }
853 
854  // Allow to disable BeamSpot externally
855  alignment::GlobalDerivatives globals(labels, derivatives);
856 
857  // Add derivatives for vertex calibration to first point of first trajectory
858  // NOTE: use GlobalDerivatives operators vector<int> and TMatrixD which filter
859  // the derivatives to not pass those with zero labels (usefull to get rid of some params)
860  std::vector<int> lab(globals); TMatrixD der(globals);
861 
862  // I want: dlocal/dext = dlocal/dtwobody * dtwobody/dext = dfdextPlusMinus * dtwobody/dext
863  TMatrixD dTwoBody_dExt(9, 7);
864  dTwoBody_dExt.Zero();
865  // beam vertex constraint
866  dTwoBody_dExt(0, 0) = 1.;
867  dTwoBody_dExt(1, 1) = 1.;
868  dTwoBody_dExt(2, 2) = 1.;
869  // beam kinematics constraint
870  dTwoBody_dExt(3, 3) = 1.;
871  dTwoBody_dExt(4, 4) = 1.;
872  dTwoBody_dExt(5, 5) = 1.;
873  // beam inv. mass constraint
874  dTwoBody_dExt(8, 6) = 1.;
875 
876  const TMatrixD dLocal_dExt = dfdextPlusMinus.first * dTwoBody_dExt;
877  TMatrixD dLocal_dExt_T = dLocal_dExt; dLocal_dExt_T.T();
878 
879  // The 5x7 transformation matrix d(q/p,u',v',u,v)/d(vx,vy,vz,px,py,pz,M) needs to be "inverted"
880  // to transform the covariance of the beamspot and boost vector of SuperKEKB into the local system
881  // of one GBL point - such that Millepede can align the beamspot (or even beam kinematics) if requested.
882  //
883  // I tested also other methods, but only the Singular Value Decomposition gives nice-enough results,
884  // with almost no code:
885  //
886  TDecompSVD svd(dLocal_dExt_T);
887  TMatrixD dExt_dLocal = svd.Invert().T();
888  //
889  // (dLocal_dExt * dExt_dLocal).Print(); // Check how close we are to unit matrix
890  //
891  // 5x5 matrix is as follows
892  //
893  // | 0 | 1 | 2 | 3 | 4 |
894  // ----------------------------------------------------------------------
895  // 0 | 1 -2.58e-17 6.939e-18 1.571e-17 -1.649e-19
896  // 1 | 1.787e-14 1 5.135e-16 -3.689e-16 -2.316e-18
897  // 2 | -1.776e-15 -7.806e-17 1 5.636e-17 6.193e-18
898  // 3 | -2.453e-15 7.26e-18 2.009e-16 1 -1.14e-16
899  // 4 | -1.689e-14 -9.593e-17 -2.317e-15 -3.396e-17 1
900  //
901  // It took me half a day to find out how to do this with 2 lines of code (3 with the include).
902  // Source: ROOT macro example - actually found at:
903  // <https://root.cern.ch/root/html/tutorials/matrix/solveLinear.C.html>
904  for (int i = 0; i < 7; ++i) {
905  for (int j = 0; j < 5; ++j) {
906  if (fabs(dExt_dLocal(i, j)) < 1.e-6)
907  dExt_dLocal(i, j) = 0.;
908  }
909  }
910  const TVectorD locRes = dLocal_dExt * extMeasurements;
911  const TMatrixD locPrec = dLocal_dExt * extPrec * dExt_dLocal;
912 
913  TMatrixDSym locPrecSym(5); locPrecSym.Zero();
914  for (int i = 0; i < 5; ++i) {
915  for (int j = i; j < 5; ++j) {
916  //locPrecSym(j, i) = locPrecSym(i, j) = locPrec(i, j);
917  locPrecSym(j, i) = locPrecSym(i, j) = (fabs(locPrec(i, j)) > 1.e-6) ? locPrec(i, j) : 0.;
918  }
919  }
920 
921  daughters[0].first[0].addMeasurement(locRes, locPrecSym);
922  if (!lab.empty())
923  daughters[0].first[0].addGlobals(lab, dfdextPlusMinus.first * der);
924 
925  //TODO: Understand this: either find a bug somewhere or improve the parametrization or .... ?
926  // This should be enough, but the parametrization seems to fail for nearly horizontal pairs...
927  //gbl::GblTrajectory combined(daughters);
928  // This should not be needed, it actually seems to make worse Chi2/NDF, but GBL does not fail.
929  // The measurement added just to be able to add the global derivatives (done just above) is redundant
930  // to the external measurement added here:
931  gbl::GblTrajectory combined(daughters, extDeriv, extMeasurements, extPrec);
932  //combined.printTrajectory(1000);
933  //combined.printPoints(1000);
934 
935  combined.fit(chi2, ndf, lostWeight);
936  getObjectPtr<TH1I>("ndf")->Fill(ndf);
937  getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
938  getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
939  if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
940  evt0 = m_eventT0->getEventT0();
941  getObjectPtr<TH1F>("evt0")->Fill(evt0);
942  }
943 
944 
945  B2RESULT("Full kinematic-constrained fit (calibration version) results NDF = " << ndf << " Chi2/NDF = " << chi2 / double(ndf));
946 
947  if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(combined);
948 
949  } else {
950 
951  gbl::GblTrajectory combined(daughters, extDeriv, extMeasurements, extPrec);
952  //combined.printTrajectory(1000);
953  //combined.printPoints(1000);
954 
955  combined.fit(chi2, ndf, lostWeight);
956  getObjectPtr<TH1I>("ndf")->Fill(ndf);
957  getObjectPtr<TH1F>("chi2_per_ndf")->Fill(chi2 / double(ndf));
958  getObjectPtr<TH1F>("pval")->Fill(TMath::Prob(chi2, ndf));
959  if (m_eventT0.isValid() && m_eventT0->hasEventT0()) {
960  evt0 = m_eventT0->getEventT0();
961  getObjectPtr<TH1F>("evt0")->Fill(evt0);
962  }
963 
964 
965  B2RESULT("Full kinematic-constrained fit results NDF = " << ndf << " Chi2/NDF = " << chi2 / double(ndf));
966 
967  if (TMath::Prob(chi2, ndf) > m_minPValue) storeTrajectory(combined);
968  }
969  }
970  }
971 }
972 
974 {
975  // We close the file at end of run, producing
976  // one file per run (and process id) which is more
977  // convenient than one large binary block.
978  auto mille = getObjectPtr<MilleData>("mille");
979  if (mille->isOpen())
980  mille->close();
981 }
982 
984 {
986 
988  if (!fileMetaData.isValid()) {
989  B2ERROR("Cannot register binaries in FileCatalog.");
990  return;
991  }
992 
993 
994  const std::vector<string> parents = {fileMetaData->getLfn()};
995  for (auto binary : getObjectPtr<MilleData>("mille")->getFiles()) {
996  FileMetaData milleMetaData(*fileMetaData);
997  // We reset filename to be set directly by the registerFile procedure
998  milleMetaData.setLfn("");
999  milleMetaData.setParents(parents);
1000  FileCatalog::Instance().registerFile(binary, milleMetaData);
1001  }
1002 
1003 }
1004 
1006 {
1007  if (m_useGblTree) {
1008  if (trajectory.isValid())
1009  m_currentGblData = trajectory.getData();
1010  else
1011  m_currentGblData.clear();
1012 
1013  if (!m_currentGblData.empty())
1014  getObjectPtr<TTree>("GblDataTree")->Fill();
1015  } else {
1016  getObjectPtr<MilleData>("mille")->fill(trajectory);
1017  }
1018 }
1019 
1021 {
1022  string name = getName();
1023 
1024  name += "-e" + to_string(m_evtMetaData->getExperiment());
1025  name += "-r" + to_string(m_evtMetaData->getRun());
1026  name += "-ev" + to_string(m_evtMetaData->getEvent());
1027 
1029  name += "-pid" + to_string(ProcHandler::EvtProcID());
1030 
1031  name += ".mille";
1032 
1033  return name;
1034 }
1035 
1037 {
1038  try {
1039  // For already fitted tracks, try to get fitted (DAF) weights for CDC
1040  if (m_updateCDCWeights && recoTrack.getNumberOfCDCHits() && recoTrack.getTrackFitStatus()
1041  && recoTrack.getTrackFitStatus()->isFitted()) {
1042  double sumCDCWeights = recoTrack.getNumberOfCDCHits(); // start with full weights
1043  // Do the hits synchronisation
1044  auto relatedRecoHitInformation =
1046 
1047  for (RecoHitInformation& recoHitInformation : relatedRecoHitInformation) {
1048 
1049  if (recoHitInformation.getFlag() == RecoHitInformation::c_pruned) {
1050  B2FATAL("Found pruned point in RecoTrack. Pruned tracks cannot be used in MillepedeCollector.");
1051  }
1052 
1053  if (recoHitInformation.getTrackingDetector() != RecoHitInformation::c_CDC) continue;
1054 
1055  const genfit::TrackPoint* trackPoint = recoTrack.getCreatedTrackPoint(&recoHitInformation);
1056  if (trackPoint) {
1057  if (not trackPoint->hasFitterInfo(recoTrack.getCardinalRepresentation()))
1058  continue;
1059  auto kalmanFitterInfo = dynamic_cast<genfit::KalmanFitterInfo*>(trackPoint->getFitterInfo());
1060  if (not kalmanFitterInfo) {
1061  continue;
1062  } else {
1063  std::vector<double> weights = kalmanFitterInfo->getWeights();
1064  if (weights.size() == 2) {
1065  if (weights.at(0) > weights.at(1))
1066  recoHitInformation.setRightLeftInformation(RecoHitInformation::c_left);
1067  else if (weights.at(0) < weights.at(1))
1068  recoHitInformation.setRightLeftInformation(RecoHitInformation::c_right);
1069 
1070  double weightLR = weights.at(0) + weights.at(1);
1071  if (weightLR < m_minCDCHitWeight) recoHitInformation.setUseInFit(false);
1072  sumCDCWeights += weightLR - 1.; // reduce weight sum if weightLR<1
1073  }
1074  }
1075  }
1076  }
1077 
1078  double usedCDCHitFraction = sumCDCWeights / double(recoTrack.getNumberOfCDCHits());
1079  getObjectPtr<TH1F>("cdc_hit_fraction")->Fill(usedCDCHitFraction);
1080  if (usedCDCHitFraction < m_minUsedCDCHitFraction)
1081  return false;
1082  }
1083  } catch (...) {
1084  B2ERROR("Error in checking DAF weights from previous fit to resolve hit ambiguity. Why? Failed fit points in DAF? Skip track to be sure.");
1085  return false;
1086  }
1087 
1088  std::shared_ptr<genfit::GblFitter> gbl(new genfit::GblFitter());
1090  gbl->setTrackSegmentController(new GblMultipleScatteringController);
1091 
1092  MeasurementAdder factory("", "", "", "", "");
1093 
1094  // We need the store arrays
1100 
1101  // Create the genfit::MeasurementFactory
1102  genfit::MeasurementFactory<genfit::AbsMeasurement> genfitMeasurementFactory;
1103 
1104  // Add producer for alignable RecoHits to factory
1105  if (pxdHits.isOptional()) {
1108  genfitMeasurementFactory.addProducer(Const::PXD, PXDProducer);
1109  }
1110 
1111  if (svdHits.isOptional()) {
1114  genfitMeasurementFactory.addProducer(Const::SVD, SVDProducer);
1115  }
1116 
1117  if (cdcHits.isOptional()) {
1120  genfitMeasurementFactory.addProducer(Const::CDC, CDCProducer);
1121  }
1122 
1123  if (bklmHits.isOptional()) {
1126  genfitMeasurementFactory.addProducer(Const::BKLM, BKLMProducer);
1127  }
1128 
1129  if (eklmHits.isOptional()) {
1132  genfitMeasurementFactory.addProducer(Const::EKLM, EKLMProducer);
1133  }
1134 
1135 
1136  // Create the measurement creators
1137  std::vector<std::shared_ptr<PXDBaseMeasurementCreator>> pxdMeasurementCreators = { std::shared_ptr<PXDBaseMeasurementCreator>(new PXDCoordinateMeasurementCreator(genfitMeasurementFactory)) };
1138  std::vector<std::shared_ptr<SVDBaseMeasurementCreator>> svdMeasurementCreators = { std::shared_ptr<SVDBaseMeasurementCreator>(new SVDCoordinateMeasurementCreator(genfitMeasurementFactory)) };
1139  // TODO: Create a new MeasurementCreator based on SVDBaseMeasurementCreator (or on SVDCoordinateMeasurementCreator), which does the combination on the fly.
1140 
1141  std::vector<std::shared_ptr<CDCBaseMeasurementCreator>> cdcMeasurementCreators = { std::shared_ptr<CDCBaseMeasurementCreator>(new CDCCoordinateMeasurementCreator(genfitMeasurementFactory)) };
1142  std::vector<std::shared_ptr<BKLMBaseMeasurementCreator>> bklmMeasurementCreators = { std::shared_ptr<BKLMBaseMeasurementCreator>(new BKLMCoordinateMeasurementCreator(genfitMeasurementFactory)) };
1143  std::vector<std::shared_ptr<EKLMBaseMeasurementCreator>> eklmMeasurementCreators = { std::shared_ptr<EKLMBaseMeasurementCreator>(new EKLMCoordinateMeasurementCreator(genfitMeasurementFactory)) };
1144 
1145  // TODO: Or put it in here and leave the svdMeasurementCreators empty.
1146  std::vector<std::shared_ptr<BaseMeasurementCreator>> additionalMeasurementCreators = {};
1147  factory.resetMeasurementCreators(pxdMeasurementCreators, svdMeasurementCreators, cdcMeasurementCreators, bklmMeasurementCreators,
1148  eklmMeasurementCreators, additionalMeasurementCreators);
1149  factory.addMeasurements(recoTrack);
1150 
1151  auto& gfTrack = RecoTrackGenfitAccess::getGenfitTrack(recoTrack);
1152 
1153  int currentPdgCode = TrackFitter::createCorrectPDGCodeForChargedStable(Const::muon, recoTrack);
1154  if (particle)
1155  currentPdgCode = particle->getPDGCode();
1156 
1157  genfit::AbsTrackRep* trackRep = RecoTrackGenfitAccess::createOrReturnRKTrackRep(recoTrack, currentPdgCode);
1158  gfTrack.setCardinalRep(gfTrack.getIdForRep(trackRep));
1159 
1160  if (particle) {
1161  B2Vector3D vertexPos = particle->getVertex();
1162  B2Vector3D vertexMom = particle->getMomentum();
1163  gfTrack.setStateSeed(vertexPos, vertexMom);
1164 
1165  genfit::StateOnPlane vertexSOP(gfTrack.getCardinalRep());
1166  B2Vector3D vertexRPhiDir(vertexPos[0], vertexPos[1], 0);
1167  B2Vector3D vertexZDir(0, 0, vertexPos[2]);
1168  //FIXME: This causes problem to current GBL version in genfit -> needs update of GBL to re-enable
1169  // genfit::SharedPlanePtr vertexPlane(new genfit::DetPlane(vertexPos, vertexRPhiDir, vertexZDir));
1170  //This works instead fine:
1171  genfit::SharedPlanePtr vertexPlane(new genfit::DetPlane(vertexPos, vertexMom));
1172 
1173  vertexSOP.setPlane(vertexPlane);
1174  vertexSOP.setPosMom(vertexPos, vertexMom);
1175  TMatrixDSym vertexCov(5);
1176  vertexCov.UnitMatrix();
1177  // By using negative covariance no measurement is added to GBL. But this first point
1178  // is then used as additional point in trajectory at the assumed point of its fitted vertex
1179  vertexCov *= -1.;
1180  genfit::MeasuredStateOnPlane mop(vertexSOP, vertexCov);
1181  genfit::FullMeasurement* vertex = new genfit::FullMeasurement(mop, Const::IR);
1182  gfTrack.insertMeasurement(vertex, 0);
1183  }
1184 
1185  try {
1186  for (unsigned int i = 0; i < gfTrack.getNumPoints() - 1; ++i) {
1187  //if (gfTrack.getPointWithMeasurement(i)->getNumRawMeasurements() != 1)
1188  // continue;
1189  genfit::PlanarMeasurement* planarMeas1 = dynamic_cast<genfit::PlanarMeasurement*>(gfTrack.getPointWithMeasurement(
1190  i)->getRawMeasurement(0));
1191  genfit::PlanarMeasurement* planarMeas2 = dynamic_cast<genfit::PlanarMeasurement*>(gfTrack.getPointWithMeasurement(
1192  i + 1)->getRawMeasurement(0));
1193 
1194  if (planarMeas1 != NULL && planarMeas2 != NULL &&
1195  planarMeas1->getDetId() == planarMeas2->getDetId() &&
1196  planarMeas1->getPlaneId() != -1 && // -1 is default plane id
1197  planarMeas1->getPlaneId() == planarMeas2->getPlaneId()) {
1198  Belle2::AlignableSVDRecoHit* hit1 = dynamic_cast<Belle2::AlignableSVDRecoHit*>(planarMeas1);
1199  Belle2::AlignableSVDRecoHit* hit2 = dynamic_cast<Belle2::AlignableSVDRecoHit*>(planarMeas2);
1200  if (hit1 && hit2) {
1201  Belle2::AlignableSVDRecoHit* hitU(NULL);
1202  Belle2::AlignableSVDRecoHit* hitV(NULL);
1203  // We have to decide U/V now (else AlignableSVDRecoHit2D could throw FATAL)
1204  if (hit1->isU() && !hit2->isU()) {
1205  hitU = hit1;
1206  hitV = hit2;
1207  } else if (!hit1->isU() && hit2->isU()) {
1208  hitU = hit2;
1209  hitV = hit1;
1210  } else {
1211  continue;
1212  }
1214  // insert measurement before point i (increases number of currect point to i+1)
1215  gfTrack.insertMeasurement(hit, i);
1216  // now delete current point (at its original place, we have the new 2D recohit)
1217  gfTrack.deletePoint(i + 1);
1218  gfTrack.deletePoint(i + 1);
1219  }
1220  }
1221  }
1222  } catch (std::exception& e) {
1223  B2ERROR(e.what());
1224  B2ERROR("SVD Cluster combination failed. This is symptomatic of pruned tracks. MillepedeCollector cannot process pruned tracks.");
1225  return false;
1226  }
1227 
1228  try {
1229  gbl->processTrackWithRep(&gfTrack, gfTrack.getCardinalRep(), true);
1230  } catch (genfit::Exception& e) {
1231  B2ERROR(e.what());
1232  return false;
1233  } catch (...) {
1234  B2ERROR("GBL fit failed.");
1235  return false;
1236  }
1237 
1238  return true;
1239 }
1240 
1241 std::vector< genfit::Track* > MillepedeCollectorModule::getParticlesTracks(std::vector<Particle*> particles, bool addVertexPoint)
1242 {
1243  std::vector< genfit::Track* > tracks;
1244  for (auto particle : particles) {
1245  auto belle2Track = particle->getTrack();
1246  if (!belle2Track) {
1247  B2WARNING("No Belle2::Track for particle (particle->X");
1248  continue;
1249  }
1250 // auto trackFitResult = belle2Track->getTrackFitResult(Const::chargedStableSet.find(abs(particle->getPDGCode())));
1251 // if (!trackFitResult) {
1252 // B2INFO("No track fit result for track");
1253 // continue;
1254 // }
1255 // auto recoTrack = trackFitResult->getRelatedFrom<RecoTrack>();
1256  auto recoTrack = belle2Track->getRelatedTo<RecoTrack>();
1257 
1258  if (!recoTrack) {
1259  B2WARNING("No related RecoTrack for Belle2::Track (particle->Track->X)");
1260  continue;
1261  }
1262 
1263  // If any track fails, fail completely
1264  if (!fitRecoTrack(*recoTrack, (addVertexPoint) ? particle : nullptr))
1265  return {};
1266 
1267  auto& track = RecoTrackGenfitAccess::getGenfitTrack(*recoTrack);
1268 
1269  if (!track.hasFitStatus()) {
1270  B2WARNING("Track has no fit status");
1271  continue;
1272  }
1273  genfit::GblFitStatus* fs = dynamic_cast<genfit::GblFitStatus*>(track.getFitStatus());
1274  if (!fs) {
1275  B2WARNING("Track FitStatus is not GblFitStatus.");
1276  continue;
1277  }
1278  if (!fs->isFittedWithReferenceTrack()) {
1279  B2WARNING("Track is not fitted with reference track.");
1280  continue;
1281  }
1282 
1283  tracks.push_back(&track);
1284  }
1285 
1286  return tracks;
1287 }
1288 
1290  double motherMass)
1291 {
1292  std::vector<TMatrixD> result;
1293 
1294  double px = mother.getPx();
1295  double py = mother.getPy();
1296  double pz = mother.getPz();
1297  double pt = sqrt(px * px + py * py);
1298  double p = mother.getMomentumMagnitude();
1299  double M = motherMass;
1300  double m = mother.getDaughter(0)->getPDGMass();
1301 
1302  if (mother.getNDaughters() != 2
1303  || m != mother.getDaughter(1)->getPDGMass()) B2FATAL("Only two same-mass daughters (V0->f+f- decays) allowed.");
1304 
1305  // Rotation matrix from mother reference system to lab system
1306  TMatrixD mother2lab(3, 3);
1307  mother2lab(0, 0) = px * pz / pt / p; mother2lab(0, 1) = - py / pt; mother2lab(0, 2) = px / p;
1308  mother2lab(1, 0) = py * pz / pt / p; mother2lab(1, 1) = px / pt; mother2lab(1, 2) = py / p;
1309  mother2lab(2, 0) = - pt / p; mother2lab(2, 1) = 0; mother2lab(2, 2) = pz / p;
1310  ROOT::Math::Rotation3D lab2mother;
1311  lab2mother.SetRotationMatrix(mother2lab); lab2mother.Invert();
1312 
1313  // Need to rotate and boost daughters' momenta to know which goes forward (+sign in decay model)
1314  // and to get the angles theta, phi of the decaying daughter system in mothers' reference frame
1315  RestFrame boostedFrame(&mother);
1316  ROOT::Math::PxPyPzEVector fourVector1 = mother.getDaughter(0)->get4Vector();
1317  ROOT::Math::PxPyPzEVector fourVector2 = mother.getDaughter(1)->get4Vector();
1318 
1319  auto mom1 = lab2mother * boostedFrame.getMomentum(fourVector1).Vect();
1320  auto mom2 = lab2mother * boostedFrame.getMomentum(fourVector2).Vect();
1321  // One momentum has opposite direction (otherwise should be same in CMS of mother), but which?
1322  double sign = 1.;
1323  auto avgMom = 0.5 * (mom1 - mom2);
1324  if (avgMom.Z() < 0.) {
1325  avgMom *= -1.;
1326  // switch meaning of plus/minus trajectories
1327  sign = -1.;
1328  }
1329 
1330  double theta = atan2(avgMom.rho(), avgMom.Z());
1331  double phi = atan2(avgMom.Y(), avgMom.X());
1332  if (phi < 0.) phi += 2. * TMath::Pi();
1333 
1334  double alpha = M / 2. / m;
1335  double c1 = m * sqrt(alpha * alpha - 1.);
1336  double c2 = 0.5 * sqrt((alpha * alpha - 1.) / alpha / alpha * (p * p + M * M));
1337 
1338  double p3 = p * p * p;
1339  double pt3 = pt * pt * pt;
1340 
1341 
1342  for (auto& track : getParticlesTracks(mother.getDaughters())) {
1343 
1344 
1345  TMatrixD R = mother2lab;
1346  B2Vector3D P(sign * c1 * sin(theta) * cos(phi),
1347  sign * c1 * sin(theta) * sin(phi),
1348  p / 2. + sign * c2 * cos(theta));
1349 
1350  TMatrixD dRdpx(3, 3);
1351  dRdpx(0, 0) = - pz * (pow(px, 4.) - pow(py, 4.) - py * py * pz * pz) / pt3 / p3;
1352  dRdpx(0, 1) = px * py / pt3;
1353  dRdpx(0, 2) = (py * py + pz * pz) / p3;
1354 
1355  dRdpx(1, 0) = - px * py * pz * (2. * px * px + 2. * py * py + pz * pz) / pt3 / p3;
1356  dRdpx(1, 1) = - py * py / pt3;
1357  dRdpx(1, 2) = px * py / p3;
1358 
1359  dRdpx(2, 0) = - px * pz * pz / pt / p3;
1360  dRdpx(2, 1) = 0.;
1361  dRdpx(2, 2) = - px * pz / p3;
1362 
1363  TMatrixD dRdpy(3, 3);
1364  dRdpy(0, 0) = - px * py * pz * (2. * px * px + 2. * py * py + pz * pz) / pt3 / p3;
1365  dRdpy(0, 1) = - px * px / pt3;
1366  dRdpy(0, 2) = px * pz / p3;
1367 
1368  dRdpy(1, 0) = - pz * (- pow(px, 4.) - px * px * pz * pz + pow(py, 4.)) / pt3 / p3;
1369  dRdpy(1, 1) = px * py / pt3;
1370  dRdpy(1, 2) = (px * px + pz * pz) / p3;
1371 
1372  dRdpy(2, 0) = - py * pz * pz / pt / p3;
1373  dRdpy(2, 1) = 0.;
1374  dRdpy(2, 2) = - py * pz / p3;
1375 
1376  TMatrixD dRdpz(3, 3);
1377  dRdpz(0, 0) = px * pt / p3;
1378  dRdpz(0, 1) = 0.;
1379  dRdpz(0, 2) = - px * pz / p3;
1380 
1381  dRdpz(1, 0) = py * pt / p3;
1382  dRdpz(1, 1) = 0.;
1383  dRdpz(1, 2) = py * pz / p3;
1384 
1385  dRdpz(2, 0) = pz * pt / p3;
1386  dRdpz(2, 1) = 0.;
1387  dRdpz(2, 2) = (px * px + py * py) / p3;
1388 
1389  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.) *
1390  (M * M + p * p) / M / M);
1391 
1392  B2Vector3D dpdpx = dRdpx * P + R * K * px * B2Vector3D(0., 0., 1.);
1393  B2Vector3D dpdpy = dRdpy * P + R * K * py * B2Vector3D(0., 0., 1.);
1394  B2Vector3D dpdpz = dRdpz * P + R * K * pz * B2Vector3D(0., 0., 1.);
1395 
1396  B2Vector3D dpdtheta = R * B2Vector3D(sign * c1 * cos(theta) * cos(phi),
1397  sign * c1 * cos(theta) * sin(phi),
1398  sign * c2 * (- sin(theta)));
1399 
1400 
1401  B2Vector3D dpdphi = R * B2Vector3D(sign * c1 * sin(theta) * (- sin(phi)),
1402  sign * c1 * sin(theta) * cos(phi),
1403  0.);
1404 
1405  double dc1dM = m * M / (2. * sqrt(M * M - 4. * m * m));
1406  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)));
1407 
1408  B2Vector3D dpdM = R * B2Vector3D(sign * sin(theta) * cos(phi) * dc1dM,
1409  sign * sin(theta) * sin(phi) * dc1dM,
1410  sign * cos(theta) * dc2dM);
1411 
1412  TMatrixD dpdz(3, 6);
1413  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);
1414  dpdz(0, 5) = dpdM(0);
1415  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);
1416  dpdz(1, 5) = dpdM(1);
1417  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);
1418  dpdz(2, 5) = dpdM(2);
1419 
1420  TMatrixD dqdv = getGlobalToLocalTransform(track->getFittedState()).GetSub(0, 4, 0, 2);
1421  TMatrixD dqdp = getGlobalToLocalTransform(track->getFittedState()).GetSub(0, 4, 3, 5);
1422  TMatrixD dfdvz(5, 9);
1423  dfdvz.SetSub(0, 0, dqdv);
1424  dfdvz.SetSub(0, 3, dqdp * dpdz);
1425 
1426  result.push_back(dfdvz);
1427 
1428  // switch sign for second trajectory
1429  sign *= -1.;
1430  }
1431 
1432  return {result[0], result[1]};
1433 }
1434 
1436 {
1437  auto state = msop;
1438  const B2Vector3D& U(state.getPlane()->getU());
1439  const B2Vector3D& V(state.getPlane()->getV());
1440  const B2Vector3D& O(state.getPlane()->getO());
1441  const B2Vector3D& W(state.getPlane()->getNormal());
1442 
1443  const double* state5 = state.getState().GetMatrixArray();
1444 
1445  double spu = 1.;
1446 
1447  const TVectorD& auxInfo = state.getAuxInfo();
1448  if (auxInfo.GetNrows() == 2
1449  || auxInfo.GetNrows() == 1) // backwards compatibility with old RKTrackRep
1450  spu = state.getAuxInfo()(0);
1451 
1452  TVectorD state7(7);
1453 
1454  state7[0] = O.X() + state5[3] * U.X() + state5[4] * V.X(); // x
1455  state7[1] = O.Y() + state5[3] * U.Y() + state5[4] * V.Y(); // y
1456  state7[2] = O.Z() + state5[3] * U.Z() + state5[4] * V.Z(); // z
1457 
1458  state7[3] = spu * (W.X() + state5[1] * U.X() + state5[2] * V.X()); // a_x
1459  state7[4] = spu * (W.Y() + state5[1] * U.Y() + state5[2] * V.Y()); // a_y
1460  state7[5] = spu * (W.Z() + state5[1] * U.Z() + state5[2] * V.Z()); // a_z
1461 
1462  // normalize dir
1463  double norm = 1. / sqrt(state7[3] * state7[3] + state7[4] * state7[4] + state7[5] * state7[5]);
1464  for (unsigned int i = 3; i < 6; ++i) state7[i] *= norm;
1465 
1466  state7[6] = state5[0]; // q/p
1467 
1468  const double AtU = state7[3] * U.X() + state7[4] * U.Y() + state7[5] * U.Z();
1469  const double AtV = state7[3] * V.X() + state7[4] * V.Y() + state7[5] * V.Z();
1470  const double AtW = state7[3] * W.X() + state7[4] * W.Y() + state7[5] * W.Z();
1471 
1472  // J_Mp matrix is d(q/p,u',v',u,v) / d(x,y,z,px,py,pz) (in is 6x6)
1473 
1474  const double qop = state7[6];
1475  const double p = state.getCharge() / qop; // momentum
1476 
1477  TMatrixD J_Mp_6x5(6, 5);
1478  J_Mp_6x5.Zero();
1479 
1480  //d(u)/d(x,y,z)
1481  J_Mp_6x5(0, 3) = U.X(); // [0][3]
1482  J_Mp_6x5(1, 3) = U.Y(); // [1][3]
1483  J_Mp_6x5(2, 3) = U.Z(); // [2][3]
1484  //d(v)/d(x,y,z)
1485  J_Mp_6x5(0, 4) = V.X(); // [0][4]
1486  J_Mp_6x5(1, 4) = V.Y(); // [1][4]
1487  J_Mp_6x5(2, 4) = V.Z(); // [2][4]
1488 
1489  // d(q/p)/d(px,py,pz)
1490  double fact = (-1.) * qop / p;
1491  J_Mp_6x5(3, 0) = fact * state7[3]; // [3][0]
1492  J_Mp_6x5(4, 0) = fact * state7[4]; // [4][0]
1493  J_Mp_6x5(5, 0) = fact * state7[5]; // [5][0]
1494  // d(u')/d(px,py,pz)
1495  fact = 1. / (p * AtW * AtW);
1496  J_Mp_6x5(3, 1) = fact * (U.X() * AtW - W.X() * AtU); // [3][1]
1497  J_Mp_6x5(4, 1) = fact * (U.Y() * AtW - W.Y() * AtU); // [4][1]
1498  J_Mp_6x5(5, 1) = fact * (U.Z() * AtW - W.Z() * AtU); // [5][1]
1499  // d(v')/d(px,py,pz)
1500  J_Mp_6x5(3, 2) = fact * (V.X() * AtW - W.X() * AtV); // [3][2]
1501  J_Mp_6x5(4, 2) = fact * (V.Y() * AtW - W.Y() * AtV); // [4][2]
1502  J_Mp_6x5(5, 2) = fact * (V.Z() * AtW - W.Z() * AtV); // [5][2]
1503 
1504  return J_Mp_6x5.T();
1505 }
1506 
1508 {
1509  auto state = msop;
1510  // get vectors and aux variables
1511  const B2Vector3D& U(state.getPlane()->getU());
1512  const B2Vector3D& V(state.getPlane()->getV());
1513  const B2Vector3D& W(state.getPlane()->getNormal());
1514 
1515  const TVectorD& state5(state.getState());
1516  double spu = 1.;
1517 
1518  const TVectorD& auxInfo = state.getAuxInfo();
1519  if (auxInfo.GetNrows() == 2
1520  || auxInfo.GetNrows() == 1) // backwards compatibility with old RKTrackRep
1521  spu = state.getAuxInfo()(0);
1522 
1523  TVectorD pTilde(3);
1524  pTilde[0] = spu * (W.X() + state5(1) * U.X() + state5(2) * V.X()); // a_x
1525  pTilde[1] = spu * (W.Y() + state5(1) * U.Y() + state5(2) * V.Y()); // a_y
1526  pTilde[2] = spu * (W.Z() + state5(1) * U.Z() + state5(2) * V.Z()); // a_z
1527 
1528  const double pTildeMag = sqrt(pTilde[0] * pTilde[0] + pTilde[1] * pTilde[1] + pTilde[2] * pTilde[2]);
1529  const double pTildeMag2 = pTildeMag * pTildeMag;
1530 
1531  const double utpTildeOverpTildeMag2 = (U.X() * pTilde[0] + U.Y() * pTilde[1] + U.Z() * pTilde[2]) / pTildeMag2;
1532  const double vtpTildeOverpTildeMag2 = (V.X() * pTilde[0] + V.Y() * pTilde[1] + V.Z() * pTilde[2]) / pTildeMag2;
1533 
1534  //J_pM matrix is d(x,y,z,px,py,pz) / d(q/p,u',v',u,v) (out is 6x6)
1535 
1536  const double qop = state5(0);
1537  const double p = state.getCharge() / qop; // momentum
1538 
1539  TMatrixD J_pM_5x6(5, 6);
1540  J_pM_5x6.Zero();
1541 
1542  // d(px,py,pz)/d(q/p)
1543  double fact = -1. * p / (pTildeMag * qop);
1544  J_pM_5x6(0, 3) = fact * pTilde[0]; // [0][3]
1545  J_pM_5x6(0, 4) = fact * pTilde[1]; // [0][4]
1546  J_pM_5x6(0, 5) = fact * pTilde[2]; // [0][5]
1547  // d(px,py,pz)/d(u')
1548  fact = p * spu / pTildeMag;
1549  J_pM_5x6(1, 3) = fact * (U.X() - pTilde[0] * utpTildeOverpTildeMag2); // [1][3]
1550  J_pM_5x6(1, 4) = fact * (U.Y() - pTilde[1] * utpTildeOverpTildeMag2); // [1][4]
1551  J_pM_5x6(1, 5) = fact * (U.Z() - pTilde[2] * utpTildeOverpTildeMag2); // [1][5]
1552  // d(px,py,pz)/d(v')
1553  J_pM_5x6(2, 3) = fact * (V.X() - pTilde[0] * vtpTildeOverpTildeMag2); // [2][3]
1554  J_pM_5x6(2, 4) = fact * (V.Y() - pTilde[1] * vtpTildeOverpTildeMag2); // [2][4]
1555  J_pM_5x6(2, 5) = fact * (V.Z() - pTilde[2] * vtpTildeOverpTildeMag2); // [2][5]
1556  // d(x,y,z)/d(u)
1557  J_pM_5x6(3, 0) = U.X(); // [3][0]
1558  J_pM_5x6(3, 1) = U.Y(); // [3][1]
1559  J_pM_5x6(3, 2) = U.Z(); // [3][2]
1560  // d(x,y,z)/d(v)
1561  J_pM_5x6(4, 0) = V.X(); // [4][0]
1562  J_pM_5x6(4, 1) = V.Y(); // [4][1]
1563  J_pM_5x6(4, 2) = V.Z(); // [4][2]
1564 
1565  return J_pM_5x6.T();
1566 
1567 }
1568 
1569 tuple<B2Vector3D, TMatrixDSym> MillepedeCollectorModule::getPrimaryVertexAndCov() const
1570 {
1571  DBObjPtr<BeamSpot> beam;
1572  return {beam->getIPPosition(), beam->getSizeCovMatrix()};
1573 }
1574 
1575 void MillepedeCollectorModule::updateMassWidthIfSet(string listName, double& mass, double& width)
1576 {
1577  if (m_customMassConfig.find(listName) != m_customMassConfig.end()) {
1578  auto massWidth = m_customMassConfig.at(listName);
1579  mass = std::get<0>(massWidth);
1580  width = std::get<1>(massWidth);
1581  }
1582 }
1583 
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:651
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.
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 useable 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, ...
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
@ 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
const std::string & getName() const
Returns the name of the module.
Definition: Module.h:187
Class to store reconstructed particles.
Definition: Particle.h:75
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
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
unsigned int getNumberOfCDCHits() const
Return the number of cdc hits.
Definition: RecoTrack.h:427
const genfit::TrackPoint * getCreatedTrackPoint(const RecoHitInformation *recoHitInformation) const
Get a pointer to the TrackPoint that was created from this hit.
Definition: RecoTrack.cc:230
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
const std::string & getStoreArrayNameOfRecoHitInformation() const
Name of the store array of the reco hit informations.
Definition: RecoTrack.h:747
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.
bool isValid() const
Check wether 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:96
bool isValid() const
Check whether the object was created.
Definition: StoreObjPtr.h:111
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 comming 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
GBL trajectory.
Definition: GblTrajectory.h:48
unsigned int fit(double &Chi2, int &Ndf, double &lostWeight, std::string optionList="")
Perform fit of (valid) trajectory.
bool isValid() const
Retrieve validity of trajectory.
Abstract base class for a track representation.
Definition: AbsTrackRep.h:66
Detector plane.
Definition: DetPlane.h:59
Exception class for error handling in GENFIT (provides storage for diagnostic information)
Definition: Exception.h:48
bool isFitted() const
Has the track been fitted?
Definition: FitStatus.h:94
Measurement class implementing a measurement of all track parameters.
FitStatus for use with GblFitter.
Definition: GblFitStatus.h:39
Generic GBL implementation.
Definition: GblFitter.h:53
Collects information needed and produced by a AbsKalmanFitter implementations and is specific to one ...
#StateOnPlane with additional covariance matrix.
void addProducer(int detID, AbsMeasurementProducer< measurement_T > *hitProd)
Register a producer module to the factory.
Template class for a measurement producer module.
Measurement class implementing a planar hit geometry (1 or 2D).
A state with arbitrary dimension defined in a DetPlane.
Definition: StateOnPlane.h:47
Object containing AbsMeasurement and AbsFitterInfo objects.
Definition: TrackPoint.h:46
AbsFitterInfo * getFitterInfo(const AbsTrackRep *rep=nullptr) const
Get fitterInfo for rep. Per default, use cardinal rep.
Definition: TrackPoint.cc:170
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:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
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.
Namespace for the general broken lines package.
std::shared_ptr< genfit::DetPlane > SharedPlanePtr
Shared Pointer to a DetPlane.