CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
KinematicConstrainedVertexFitterT.h
Go to the documentation of this file.
1 #ifndef KinematicConstrainedVertexFitterT_H
2 #define KinematicConstrainedVertexFitterT_H
3 
9 
11 
14 
24 template < int nTrk, int nConstraint> class KinematicConstrainedVertexFitterT{
25 
26 public:
27 
31  explicit KinematicConstrainedVertexFitterT(const MagneticField* ifield);
32 
37 
39 
45  void setParameters(const edm::ParameterSet& pSet);
46 
51  RefCountedKinematicTree fit(const std::vector<RefCountedKinematicParticle> &part) {
52  return fit(part, 0, 0);
53  }
54 
58  RefCountedKinematicTree fit(const std::vector<RefCountedKinematicParticle> &part,
60  return fit(part, cs, 0);
61  };
62 
66  RefCountedKinematicTree fit(const std::vector<RefCountedKinematicParticle> &part,
68  GlobalPoint * pt);
69 
70  //return the number of iterations
71  int getNit() const;
72  //return the value of the constraint equation
73  float getCSum() const;
74 
75 private:
76 
77  void defaultParameters();
78 
84 
85  float theMaxDelta; //maximum (delta parameter)^2/(sigma parameter)^2 per iteration for convergence
86  int theMaxStep;
87  float theMaxReducedChiSq; //max of initial (after 2 iterations) chisq/dof value
88  float theMinChiSqImprovement; //minimum required improvement in chisq to avoid fit termination for cases exceeding theMaxReducedChiSq
89 
90 
92  float csum;
93 };
94 
95 
100 
101 
102 
103 template < int nTrk, int nConstraint>
105  field(ifield)
106 {
112  iterations = -1;
113  csum = -1000.0;
114 }
115 
116 template < int nTrk, int nConstraint>
118  field(ifield)
119 {
120  finder = fnd.clone();
125  iterations = -1;
126  csum = -1000.0;
127 }
128 
129 template < int nTrk, int nConstraint>
131 {
132  delete finder;
133  delete vCons;
134  delete updator;
135  delete tBuilder;
136 }
137 
138 template < int nTrk, int nConstraint>
140 {
141  theMaxDelta = pSet.getParameter<double>("maxDelta");
142  theMaxStep = pSet.getParameter<int>("maxNbrOfIterations");
143  theMaxReducedChiSq = pSet.getParameter<double>("maxReducedChiSq");
144  theMinChiSqImprovement = pSet.getParameter<double>("minChiSqImprovement");
145 }
146 
147 template < int nTrk, int nConstraint>
149 {
150  theMaxDelta = 0.01;
151  theMaxStep = 1000;
152  theMaxReducedChiSq = 225.;
153  theMinChiSqImprovement = 50.;
154 
155 }
156 
157 template < int nTrk, int nConstraint>
159 KinematicConstrainedVertexFitterT< nTrk, nConstraint>::fit(const std::vector<RefCountedKinematicParticle> &part,
161  GlobalPoint * pt)
162 {
163  assert( nConstraint==0 || cs!=0);
164  if(part.size()!=nTrk) throw VertexException("KinematicConstrainedVertexFitterT::input states are not nTrk");
165 
166  //sorting out the input particles
167  InputSort iSort;
168  std::pair<std::vector<RefCountedKinematicParticle>, std::vector<FreeTrajectoryState> > input = iSort.sort(part);
169  const std::vector<RefCountedKinematicParticle> & particles = input.first;
170  const std::vector<FreeTrajectoryState> & fStates = input.second;
171 
172  // linearization point:
173  GlobalPoint linPoint = (pt!=0) ? *pt : finder->getLinearizationPoint(fStates);
174 
175  //initial parameters:
176  ROOT::Math::SVector<double,3+7*nTrk> inPar; //3+ 7*ntracks
177  ROOT::Math::SVector<double,3+7*nTrk> finPar; //3+ 7*ntracks
178 
179  ROOT::Math::SMatrix<double, 3+7*nTrk,3+7*nTrk ,ROOT::Math::MatRepSym<double,3+7*nTrk> > inCov;
180 
181  //making initial vector of parameters and initial particle-related covariance
182  int nSt = 0;
183  std::vector<KinematicState> lStates(nTrk);
184  for(std::vector<RefCountedKinematicParticle>::const_iterator i = particles.begin(); i!=particles.end(); i++)
185  {
186  lStates[nSt] = (*i)->stateAtPoint(linPoint);
187  KinematicState const & state = lStates[nSt];
188  if (!state.isValid()) {
189  LogDebug("KinematicConstrainedVertexFitter")
190  << "State is invalid at point: "<<linPoint<<std::endl;
192  }
193  inPar.Place_at(state.kinematicParameters().vector(),3+7*nSt);
194  inCov.Place_at(state.kinematicParametersError().matrix(),3 + 7*nSt,3 + 7*nSt);
195  ++nSt;
196  }
197 
198  //initial vertex error matrix components (huge error method)
199  //and vertex related initial vector components
200  double in_er = 100.;
201  inCov(0,0) = in_er;
202  inCov(1,1) = in_er;
203  inCov(2,2) = in_er;
204 
205  inPar(0) = linPoint.x();
206  inPar(1) = linPoint.y();
207  inPar(2) = linPoint.z();
208 
209  //constraint equations value and number of iterations
210  double eq;
211  int nit = 0;
212  iterations = 0;
213  csum = 0.0;
214 
215  GlobalPoint lPoint = linPoint;
217  ROOT::Math::SMatrix<double, 3+7*nTrk,3+7*nTrk ,ROOT::Math::MatRepSym<double,3+7*nTrk> > refCCov;
218 
219  double chisq = 1e6;
220  bool convergence = false;
221 
222  //iterarions over the updator: each time updated parameters
223  //are taken as new linearization point
224  do{
225  eq = 0.;
226  refCCov = inCov;
227  std::vector<KinematicState> oldStates = lStates;
228  GlobalVector mf = field->inInverseGeV(lPoint);
229  rVtx = updator->update(inPar,refCCov,lStates,lPoint,mf,cs);
230  if (particles.size() != lStates.size() || rVtx == 0) {
231  LogDebug("KinematicConstrainedVertexFitter")
232  << "updator failure\n";
234  }
235 
236  double newchisq = rVtx->chiSquared();
237  if ( nit>2 && newchisq > theMaxReducedChiSq*rVtx->degreesOfFreedom() && (newchisq-chisq) > (-theMinChiSqImprovement) ) {
238  LogDebug("KinematicConstrainedVertexFitter")
239  << "bad chisq and insufficient improvement, bailing\n";
241  }
242  chisq = newchisq;
243 
244 
245 
246  const GlobalPoint &newPoint = rVtx->position();
247 
248  double maxDelta = 0.0;
249 
250  double deltapos[3];
251  deltapos[0] = newPoint.x() - lPoint.x();
252  deltapos[1] = newPoint.y() - lPoint.y();
253  deltapos[2] = newPoint.z() - lPoint.z();
254  for (int i=0; i<3; ++i) {
255  double delta = deltapos[i]*deltapos[i]/rVtx->error().matrix_new()(i,i);
256  if (delta>maxDelta) maxDelta = delta;
257  }
258 
259  for (std::vector<KinematicState>::const_iterator itold = oldStates.begin(), itnew = lStates.begin();
260  itnew!=lStates.end(); ++itold,++itnew) {
261  for (int i=0; i<7; ++i) {
262  double deltapar = itnew->kinematicParameters()(i) - itold->kinematicParameters()(i);
263  double delta = deltapar*deltapar/itnew->kinematicParametersError().matrix()(i,i);
264  if (delta>maxDelta) maxDelta = delta;
265  }
266  }
267 
268  lPoint = newPoint;
269 
270  nit++;
271  convergence = maxDelta<theMaxDelta || (nit==theMaxStep && maxDelta<4.0*theMaxDelta);
272 
273  }while(nit<theMaxStep && !convergence);
274 
275  if (!convergence) {
277  }
278 
279  // std::cout << "new full cov matrix" << std::endl;
280  // std::cout << refCCov << std::endl;
281 
282  iterations = nit;
283  csum = eq;
284 
285  return tBuilder->buildTree<nTrk>(particles, lStates, rVtx, refCCov);
286 
287 }
288 
289 template < int nTrk, int nConstraint>
291  return iterations;
292 }
293 
294 template < int nTrk, int nConstraint>
296  return csum;
297 }
298 
299 #endif
#define LogDebug(id)
dbl * delta
Definition: mlp_gen.cc:36
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
AlgebraicVector7 const & vector() const
The full vector (7 elements)
auto_ptr< ClusterSequence > cs
bool isValid() const
Common base class.
T y() const
Definition: PV3DBase.h:63
std::pair< std::vector< RefCountedKinematicParticle >, std::vector< FreeTrajectoryState > > sort(const std::vector< RefCountedKinematicParticle > &particles) const
Definition: InputSort.cc:6
void setParameters(const edm::ParameterSet &pSet)
KinematicConstrainedVertexFitterT(const MagneticField *ifield)
T z() const
Definition: PV3DBase.h:64
KinematicParametersError const & kinematicParametersError() const
KinematicParameters const & kinematicParameters() const
part
Definition: HCALResponse.h:20
RefCountedKinematicTree fit(const std::vector< RefCountedKinematicParticle > &part, MultiTrackKinematicConstraintT< nTrk, nConstraint > *cs)
char state
Definition: procUtils.cc:75
KinematicConstrainedVertexUpdatorT< nTrk, nConstraint > * updator
virtual LinearizationPointFinder * clone() const =0
RefCountedKinematicTree fit(const std::vector< RefCountedKinematicParticle > &part)
T x() const
Definition: PV3DBase.h:62
AlgebraicSymMatrix77 matrix() const