CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TFitConstraintM.cc
Go to the documentation of this file.
1 // Classname: TFitConstraintM
2 // Author: Jan E. Sundermann, Verena Klose (TU Dresden)
3 
4 
5 //________________________________________________________________
6 //
7 // TFitConstraintM::
8 // --------------------
9 //
10 // Fit constraint: mass conservation ( m_i - m_j - MassConstraint == 0 )
11 //
12 
13 #include <iostream>
14 #include <iomanip>
17 #include "TLorentzVector.h"
18 #include "TClass.h"
19 
20 
21 //----------------
22 // Constructor --
23 //----------------
26  ,_ParList1(0)
27  ,_ParList2(0)
28  ,_TheMassConstraint(0)
29 {
30 
31 }
32 
33 TFitConstraintM::TFitConstraintM(std::vector<TAbsFitParticle*>* ParList1,
34  std::vector<TAbsFitParticle*>* ParList2, Double_t Mass)
36  ,_ParList1(0)
37  ,_ParList2(0)
38 {
39  // ParList1: Vector containing first list of constrained particles
40  // ( sum[ m_i ] - sum[ m_j ] - MassConstraint == 0 )
41  // ParList2: Vector containing second list of constrained particles
42  // ( sum[ m_i ] - sum[ m_j ] - MassConstraint == 0 )
43 
44  if (ParList1) {
45  _ParList1 = (*ParList1);
46  }
47  if (ParList2) {
48  _ParList2 = (*ParList2);
49  }
50  if (Mass >= 0) {
52  }
53  else if(Mass < 0) {
54  edm::LogWarning ("NegativeMassConstr")
55  << "Mass constraint in TFitConstraintM cannot be set to a negative value, will be set to 0.";
56  _TheMassConstraint = 0.;
57  }
58 }
59 
60 TFitConstraintM::TFitConstraintM(const TString &name, const TString &title,
61  std::vector<TAbsFitParticle*>* ParList1,
62  std::vector<TAbsFitParticle*>* ParList2, Double_t Mass)
63  : TAbsFitConstraint(name, title)
64  ,_ParList1(0)
65  ,_ParList2(0)
66 {
67  // ParList1: Vector containing first list of constrained particles
68  // ( sum[ m_i ] - sum[ m_j ] - MassConstraint == 0 )
69  // ParList2: Vector containing second list of constrained particles
70  // ( sum[ m_i ] - sum[ m_j ] - MassConstraint == 0 )
71 
72  if (ParList1) {
73  _ParList1 = (*ParList1);
74  }
75  if (ParList2) {
76  _ParList2 = (*ParList2);
77  }
78  if (Mass >= 0) {
80  }
81  else if(Mass < 0) {
82  edm::LogWarning ("NegativeMassConstr")
83  << "Mass constraint in TFitConstraintM cannot be set to a negative value, will be set to 0.";
84  _TheMassConstraint = 0.;
85  }
86 }
88  // Add one constrained particle to first list of particles
89  // ( sum[ m_i ] - sum[ m_j ] - MassConstraint == 0 )
90 
91  _ParList1.push_back( particle );
92 
93 }
94 
96  // Add one constrained particle to second list of particles
97  // ( sum[ m_i ] - sum[ m_j ] - MassConstraint == 0 )
98 
99  _ParList2.push_back( particle );
100 
101 }
102 
105  TAbsFitParticle* p9, TAbsFitParticle* p10) {
106  // Add many constrained particle to first list of particles
107  // ( sum[ m_i ] - sum[ m_j ] - MassConstraint == 0 )
108 
109  if (p1) addParticle1( p1 );
110  if (p2) addParticle1( p2 );
111  if (p3) addParticle1( p3 );
112  if (p4) addParticle1( p4 );
113  if (p5) addParticle1( p5 );
114  if (p6) addParticle1( p6 );
115  if (p7) addParticle1( p7 );
116  if (p8) addParticle1( p8 );
117  if (p9) addParticle1( p9 );
118  if (p10) addParticle1( p10 );
119 
120 }
121 
124  TAbsFitParticle* p9, TAbsFitParticle* p10) {
125  // Add many constrained particle to second list of particles
126  // ( sum[ m_i ] - sum[ m_j ] - MassConstraint == 0 )
127 
128  if (p1) addParticle2( p1 );
129  if (p2) addParticle2( p2 );
130  if (p3) addParticle2( p3 );
131  if (p4) addParticle2( p4 );
132  if (p5) addParticle2( p5 );
133  if (p6) addParticle2( p6 );
134  if (p7) addParticle2( p7 );
135  if (p8) addParticle2( p8 );
136  if (p9) addParticle2( p9 );
137  if (p10) addParticle2( p10 );
138 
139 }
140 
141 
142 //--------------
143 // Destructor --
144 //--------------
146 
147 }
148 
149 //--------------
150 // Operations --
151 //--------------
153  // returns derivative df/dP with P=(p,E) and f the constraint (f=0).
154  // The matrix contains one row (df/dp, df/dE).
155 
156  TMatrixD* DerivativeMatrix = new TMatrixD(1,4);
157  (*DerivativeMatrix) *= 0.;
158 
159  // Pf[4] is the 4-Mom (p,E) of the sum of particles on
160  // the list particle is part of
161 
162  Double_t Factor = 0.;
163  TLorentzVector Pf(0., 0., 0., 0.);
164 
165  if( OnList( &_ParList1, particle) ) {
166  UInt_t Npart = _ParList1.size();
167  for (unsigned int i=0; i<Npart; i++) {
168  const TLorentzVector* FourVec = (_ParList1[i])->getCurr4Vec();
169  Pf += (*FourVec);
170  }
171  if( Pf.M() == 0. ) {
172  edm::LogInfo ("KinFitter")
173  << "Division by zero in "
174  << IsA()->GetName() << " (named " << GetName() << ", titled " << GetTitle()
175  << ") will lead to Inf in derivative matrix for particle "
176  << particle->GetName() << ".";
177  }
178  Factor = 1./ Pf.M();
179  } else if (OnList( &_ParList2, particle) ) {
180  UInt_t Npart = _ParList2.size();
181  for (unsigned int i=0; i<Npart; i++) {
182  const TLorentzVector* FourVec = (_ParList2[i])->getCurr4Vec();
183  Pf += (*FourVec);
184  }
185  if( Pf.M() == 0. ) {
186  edm::LogInfo ("KinFitter")
187  << "Division by zero in "
188  << IsA()->GetName() << " (named " << GetName() << ", titled " << GetTitle()
189  << ") will lead to Inf in derivative matrix for particle "
190  << particle->GetName() << ".";
191  }
192  Factor = -1./Pf.M();
193  } else {
194  Factor = 0.;
195  }
196 
197  (*DerivativeMatrix)(0,0) = -Pf[0] ;
198  (*DerivativeMatrix)(0,1) = -Pf[1];
199  (*DerivativeMatrix)(0,2) = -Pf[2];
200  (*DerivativeMatrix)(0,3) = +Pf[3];
201  (*DerivativeMatrix) *= Factor;
202 
203  return DerivativeMatrix;
204 
205 }
206 
207 
209  // Get initial value of constraint (before the fit)
210 
211  Double_t InitValue(0) ;
212  InitValue = CalcMass(&_ParList1,true) - CalcMass(&_ParList2,true) - _TheMassConstraint ;
213  return InitValue;
214 }
215 
217  // Get value of constraint after the fit
218 
219  Double_t CurrentValue(0);
220  CurrentValue= CalcMass(&_ParList1,false) - CalcMass(&_ParList2,false) - _TheMassConstraint;
221  return CurrentValue;
222 }
223 
224 
225 Bool_t TFitConstraintM::OnList(std::vector<TAbsFitParticle*>* List,
226  TAbsFitParticle* particle) {
227  // Checks whether list contains given particle
228 
229  Bool_t ok(false);
230  UInt_t Npart = List->size();
231  for (unsigned int i=0;i<Npart;i++) {
232  ok = (particle == (*List)[i]);
233  if (ok) break;
234  }
235  return ok;
236 }
237 
238 Double_t TFitConstraintM::CalcMass(std::vector<TAbsFitParticle*>* List, Bool_t IniVal) {
239  // Calculates initial/current invariant mass of provided list of particles
240 
241  TLorentzVector P(0., 0., 0., 0.);
242  UInt_t Npart = List->size();
243  for (unsigned int i=0;i<Npart;i++) {
244  const TLorentzVector* FourVec = 0;
245  if (IniVal)
246  FourVec = ((*List)[i])->getIni4Vec();
247  else
248  FourVec = ((*List)[i])->getCurr4Vec();
249  P += (*FourVec);
250  }
251  return P.M();
252 }
253 
255  // Collect information to be used for printout
256 
257  std::stringstream info;
258  info << std::scientific << std::setprecision(6);
259 
260  info << "__________________________" << std::endl
261  << std::endl;
262  info <<"OBJ: " << IsA()->GetName() << "\t" << GetName() << "\t" << GetTitle() << std::endl;
263 
264  info << "initial value: " << getInitValue() << std::endl;
265  info << "current value: " << getCurrentValue() << std::endl;
266  info << "mass: " << _TheMassConstraint << std::endl;
267 
268  return info.str();
269 
270 }
271 
273  // Print constraint contents
274 
275  edm::LogVerbatim("KinFitter") << this->getInfoString();
276 
277 }
int i
Definition: DBlmapReader.cc:9
void addParticles1(TAbsFitParticle *p1, TAbsFitParticle *p2=0, TAbsFitParticle *p3=0, TAbsFitParticle *p4=0, TAbsFitParticle *p5=0, TAbsFitParticle *p6=0, TAbsFitParticle *p7=0, TAbsFitParticle *p8=0, TAbsFitParticle *p9=0, TAbsFitParticle *p10=0)
Double_t CalcMass(std::vector< TAbsFitParticle * > *List, Bool_t IniVal)
virtual ~TFitConstraintM()
#define P
Double_t _TheMassConstraint
Bool_t OnList(std::vector< TAbsFitParticle * > *List, TAbsFitParticle *particle)
double p4[4]
Definition: TauolaWrapper.h:92
void addParticles2(TAbsFitParticle *p1, TAbsFitParticle *p2=0, TAbsFitParticle *p3=0, TAbsFitParticle *p4=0, TAbsFitParticle *p5=0, TAbsFitParticle *p6=0, TAbsFitParticle *p7=0, TAbsFitParticle *p8=0, TAbsFitParticle *p9=0, TAbsFitParticle *p10=0)
double p2[4]
Definition: TauolaWrapper.h:90
unsigned int UInt_t
Definition: FUTypes.h:12
virtual Double_t getInitValue()
virtual TMatrixD * getDerivative(TAbsFitParticle *particle)
std::vector< TAbsFitParticle * > _ParList1
double p1[4]
Definition: TauolaWrapper.h:89
void addParticle1(TAbsFitParticle *particle)
virtual Double_t getCurrentValue()
virtual TString getInfoString()
std::vector< TAbsFitParticle * > _ParList2
virtual void print()
void addParticle2(TAbsFitParticle *particle)
double p3[4]
Definition: TauolaWrapper.h:91