CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
LagrangeParentParticleFitter Class Reference

#include <LagrangeParentParticleFitter.h>

Inheritance diagram for LagrangeParentParticleFitter:
ParentParticleFitter

Public Member Functions

LagrangeParentParticleFitterclone () const override
 
std::vector< RefCountedKinematicTreefit (const std::vector< RefCountedKinematicTree > &trees, KinematicConstraint *cs) const override
 
 LagrangeParentParticleFitter ()
 
void setParameters (const edm::ParameterSet &pSet)
 
 ~LagrangeParentParticleFitter () override
 
- Public Member Functions inherited from ParentParticleFitter
 ParentParticleFitter ()
 
virtual ~ParentParticleFitter ()
 

Private Member Functions

void defaultParameters ()
 

Private Attributes

float theMaxDiff
 
int theMaxStep
 

Detailed Description

KinematicParticle refit via LMS minimization with Lagrange multipliers method. Tunable for number of iterations and stopping condition. Every new iteration takes the result of previous iteration as a linearization point. This is only needed in case of nonlinear constraints. In the linear case result will be reached after the first step. Stopping condition at the moment: sum of abs. constraint equations values at given point

Definition at line 18 of file LagrangeParentParticleFitter.h.

Constructor & Destructor Documentation

◆ LagrangeParentParticleFitter()

LagrangeParentParticleFitter::LagrangeParentParticleFitter ( )

Definition at line 6 of file LagrangeParentParticleFitter.cc.

References defaultParameters().

Referenced by clone().

◆ ~LagrangeParentParticleFitter()

LagrangeParentParticleFitter::~LagrangeParentParticleFitter ( )
inlineoverride

Definition at line 22 of file LagrangeParentParticleFitter.h.

22 {}

Member Function Documentation

◆ clone()

LagrangeParentParticleFitter* LagrangeParentParticleFitter::clone ( ) const
inlineoverridevirtual

Clone method

Implements ParentParticleFitter.

Definition at line 41 of file LagrangeParentParticleFitter.h.

References LagrangeParentParticleFitter().

◆ defaultParameters()

void LagrangeParentParticleFitter::defaultParameters ( )
private

◆ fit()

std::vector< RefCountedKinematicTree > LagrangeParentParticleFitter::fit ( const std::vector< RefCountedKinematicTree > &  trees,
KinematicConstraint cs 
) const
overridevirtual

Refit method taking KinematicConstraint and vector of trees as imput. Only top particles of corresponding trees are refitted, vertex is not created. Number of trees should be one (single track refit) or greater (multiple track refit). Some constraints may not work with single tracks (back to back for ex.)

Implements ParentParticleFitter.

Definition at line 120 of file LagrangeParentParticleFitter.cc.

References funct::abs(), callgraph::cs, hgcalPerformanceValidation::df, l1ctLayer1_cff::dr, mps_fire::i, dqmiolumiharvest::j, dqmdumpme::k, LogDebug, KinematicParametersError::matrix(), L1TEGammaOffline_cfi::nVertex, InputSort::sort(), theMaxDiff, theMaxStep, dqmdumpme::trees, and KinematicVertexFactory::vertex().

Referenced by trackingPlots.Iteration::modules().

121  {
122  InputSort iSort;
123  std::vector<RefCountedKinematicParticle> prt = iSort.sort(trees);
124  int nStates = prt.size();
125 
126  //making full initial parameters and covariance
127  AlgebraicVector part(7 * nStates, 0);
128  AlgebraicSymMatrix cov(7 * nStates, 0);
129 
130  AlgebraicVector chi_in(nStates, 0);
131  AlgebraicVector ndf_in(nStates, 0);
132  int l_c = 0;
133  for (std::vector<RefCountedKinematicParticle>::const_iterator i = prt.begin(); i != prt.end(); i++) {
134  AlgebraicVector7 lp = (*i)->currentState().kinematicParameters().vector();
135  for (int j = 1; j != 8; j++) {
136  part(7 * l_c + j) = lp(j - 1);
137  }
138  AlgebraicSymMatrix lc = asHepMatrix<7>((*i)->currentState().kinematicParametersError().matrix());
139  cov.sub(7 * l_c + 1, lc);
140  chi_in(l_c + 1) = (*i)->chiSquared();
141  ndf_in(l_c + 1) = (*i)->degreesOfFreedom();
142  l_c++;
143  }
144  //refitted parameters and covariance matrix:
145  //simple and symmetric
146  AlgebraicVector refPar;
147  AlgebraicSymMatrix refCovS;
148 
149  //constraint values, derivatives and deviations:
150  AlgebraicVector vl;
152  AlgebraicVector dev;
153  int nstep = 0;
154  double df = 0.;
155  AlgebraicVector exPoint = part;
156 
157  //this piece of code should later be refactored:
158  // The algorithm is the same as above, but the
159  // number of refitted particles is >1. Smart way of
160  // refactoring should be chosen for it.
161  AlgebraicVector chi;
162  AlgebraicVector ndf;
163  // cout<<"Starting the main loop"<<endl;
164  do {
165  df = 0.;
166  chi = chi_in;
167  ndf = ndf_in;
168  // cout<<"Iterational linearization point: "<<exPoint<<endl;
169  vl = cs->value(exPoint).first;
170  dr = cs->derivative(exPoint).first;
171  dev = cs->deviations(nStates);
172 
173  // cout<<"The value : "<<vl<<endl;
174  // cout<<"The derivative: "<<dr<<endl;
175  // cout<<"deviations: "<<dev<<endl;
176  // cout<<"covariance "<<cov<<endl;
177 
178  //residual between expansion and current parameters
179  //0 at the first step
180  AlgebraicVector delta_alpha = part - exPoint;
181 
182  //parameters needed for refit
183  // v_d = (D * V_alpha & * D.T)^(-1)
184  AlgebraicMatrix drt = dr.T();
185  AlgebraicMatrix v_d = dr * cov * drt;
186  int ifail = 0;
187  v_d.invert(ifail);
188  if (ifail != 0) {
189  LogDebug("KinematicConstrainedVertexFitter") << "Fit failed: unable to invert covariance matrix\n";
190  return std::vector<RefCountedKinematicTree>();
191  }
192 
193  //lagrangian multipliers
194  //lambda = V_d * (D * delta_alpha + d)
195  AlgebraicVector lambda = v_d * (dr * delta_alpha + vl);
196 
197  //refitted parameters
198  refPar = part - (cov * drt * lambda);
199 
200  //refitted covariance: simple and SymMatrix
201  refCovS = cov;
202  AlgebraicMatrix sPart = drt * v_d * dr;
203  AlgebraicMatrix covF = cov * sPart * cov;
204 
205  //total covariance symmatrix
206  AlgebraicSymMatrix sCovF(nStates * 7, 0);
207  for (int i = 1; i < nStates * 7 + 1; ++i) {
208  for (int j = 1; j < nStates * 7 + 1; j++) {
209  if (i <= j)
210  sCovF(i, j) = covF(i, j);
211  }
212  }
213 
214  refCovS -= sCovF;
215 
216  // cout<<"Fitter: refitted covariance "<<refCovS<<endl;
217  for (int i = 1; i < nStates + 1; i++) {
218  for (int j = 1; j < 8; j++) {
219  refCovS((i - 1) + j, (i - 1) + j) += dev(j);
220  }
221  }
222 
223  //chiSquared
224  for (int k = 1; k < nStates + 1; k++) {
225  chi(k) += (lambda.T() * (dr * delta_alpha + vl))(1);
226  ndf(k) += cs->numberOfEquations();
227  }
228  //new expansionPoint
229  exPoint = refPar;
230  AlgebraicVector vlp = cs->value(exPoint).first;
231  for (int i = 1; i < vl.num_row(); i++) {
232  df += std::abs(vlp(i));
233  }
234  nstep++;
235  } while (df > theMaxDiff && nstep < theMaxStep);
236  //here math and iterative part is finished, starting an output production
237  //creating an output KinematicParticle and putting it on its place in the tree
238 
239  //vector of refitted particles and trees
240  std::vector<RefCountedKinematicParticle> refPart;
241  std::vector<RefCountedKinematicTree> refTrees = trees;
242 
243  int j = 1;
244  std::vector<RefCountedKinematicTree>::const_iterator tr = refTrees.begin();
245  for (std::vector<RefCountedKinematicParticle>::const_iterator i = prt.begin(); i != prt.end(); i++) {
246  AlgebraicVector7 lRefPar;
247  for (int k = 1; k < 8; k++) {
248  lRefPar(k - 1) = refPar((j - 1) * 7 + k);
249  }
250  AlgebraicSymMatrix77 lRefCovS = asSMatrix<7>(refCovS.sub((j - 1) * 7 + 1, (j - 1) * 7 + 7));
251 
252  //new refitted parameters and covariance
253  KinematicParameters param(lRefPar);
254  KinematicParametersError er(lRefCovS);
255  KinematicState kState(param, er, (*i)->initialState().particleCharge(), (**i).magneticField());
256  RefCountedKinematicParticle refParticle = (*i)->refittedParticle(kState, chi(j), ndf(j), cs->clone());
257 
258  //replacing particle itself
259  (*tr)->findParticle(*i);
260  RefCountedKinematicVertex inVertex = (*tr)->currentDecayVertex();
261  (*tr)->replaceCurrentParticle(refParticle);
262 
263  //replacing the vertex with its refitted version
264  GlobalPoint nvPos(param.position());
265  AlgebraicSymMatrix nvMatrix = asHepMatrix<7>(er.matrix()).sub(1, 3);
266  GlobalError nvError(asSMatrix<3>(nvMatrix));
267  VertexState vState(nvPos, nvError, 1.0);
268  KinematicVertexFactory vFactory;
269  RefCountedKinematicVertex nVertex = vFactory.vertex(vState, inVertex, chi(j), ndf(j));
270  (*tr)->replaceCurrentVertex(nVertex);
271  tr++;
272  j++;
273  }
274  return refTrees;
275 }
std::pair< std::vector< RefCountedKinematicParticle >, std::vector< FreeTrajectoryState > > sort(const std::vector< RefCountedKinematicParticle > &particles) const
Definition: InputSort.cc:5
ROOT::Math::SVector< double, 7 > AlgebraicVector7
Definition: Matrices.h:8
CLHEP::HepMatrix AlgebraicMatrix
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
CLHEP::HepVector AlgebraicVector
part
Definition: HCALResponse.h:20
ROOT::Math::SMatrix< double, 7, 7, ROOT::Math::MatRepSym< double, 7 > > AlgebraicSymMatrix77
Definition: Matrices.h:9
CLHEP::HepSymMatrix AlgebraicSymMatrix
static RefCountedKinematicVertex vertex(const VertexState &state, float totalChiSq, float degreesOfFr)
#define LogDebug(id)

◆ setParameters()

void LagrangeParentParticleFitter::setParameters ( const edm::ParameterSet pSet)

Configuration through PSet: number of iterations(maxDistance) and stopping condition (maxNbrOfIterations)

Definition at line 277 of file LagrangeParentParticleFitter.cc.

References edm::ParameterSet::getParameter(), theMaxDiff, and theMaxStep.

277  {
278  theMaxDiff = pSet.getParameter<double>("maxDistance");
279  theMaxStep = pSet.getParameter<int>("maxNbrOfIterations");
280  ;
281 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307

Member Data Documentation

◆ theMaxDiff

float LagrangeParentParticleFitter::theMaxDiff
private

Definition at line 46 of file LagrangeParentParticleFitter.h.

Referenced by defaultParameters(), fit(), and setParameters().

◆ theMaxStep

int LagrangeParentParticleFitter::theMaxStep
private

Definition at line 47 of file LagrangeParentParticleFitter.h.

Referenced by defaultParameters(), fit(), and setParameters().