ATLAS Offline Software
Loading...
Searching...
No Matches
TrkVKalVrtCoreBase.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include <algorithm>
6#include <cmath>
7#include <memory>
8#include <iostream>
9
15
16namespace Trk {
17
19 const basePropagator* baseP, const addrPropagator addrP,
20 IVKalState* istate):
21 vk_objMagFld(baseFld),
22 vk_funcMagFld(addrFld),
23 vk_objProp(baseP),
24 vk_funcProp(addrP),
25 vk_istate(istate)
26 {}
27
48
49 VKTrack::VKTrack(long int iniId,
50 const double Perigee[],
51 const double Covariance[],
52 VKVertex* vk,
53 double m)
54 : Id(iniId)
55 , Charge(1)
56 , Chi2(0)
57 , m_mass(m)
58 {
59 if (Perigee[4] < 0)
60 Charge = -1;
61 for (int i = 0; i < 3; i++) {
62 fitP[i] = 0;
63 cnstP[i] = 0.;
64 iniP[i] = 0.;}
65 std::copy(Perigee, Perigee+5, Perig);
66 std::copy(Perigee, Perigee+5, refPerig);
67 for(int i=0; i<5; i++) { rmnd[i]=0;}
68 std::copy(Covariance, Covariance+15, refCovar);
69 m_originVertex = vk;
70 }
71
72 void VKTrack::setCurrent(const double Perigee[], const double Weight[])
73 {
74 std::copy(Perigee, Perigee+5, Perig);
75 std::copy(Weight, Weight+15, WgtM);
76 std::copy(Weight, Weight+15, WgtM_save);
77 }
78 void VKTrack::setReference(const double Perigee[], const double Covariance[])
79 {
80 std::copy(Perigee, Perigee+5, refPerig);
81 std::copy(Covariance, Covariance+15, refCovar);
82 }
84 {
85 std::copy(WgtM_save, WgtM_save+15, WgtM);
86 }
87
88// VKTrack::VKTrack(const VKTrack & src ) // copy operator
89// { Charge =src.Charge;
90// m_mass =src.m_mass;
91// for(int i=0; i<5; i++) {Perig[i]=src.Perig[i];refPerig[i]=src.refPerig[i];}
92// for(int i=0; i<15; i++) {refCovar[i]=src.refCovar[i];WgtM[i]=src.WgtM[i];}
93// for(int i=0; i<3; i++) {fitP[i]=src.fitP[i];}
94// }
95
96 std::ostream & operator << ( std::ostream& out, const VKTrack& track )
97 {
98 //out.setf( std::ios::defaultfloat ); out.precision(5); out << std::endl;
99 out.precision(5); out << std::defaultfloat;
100 out << " Track par: Iteration <-> Ref" << std::endl;
101 out << " * a_0 : "<< track.a0() <<" "<< track.r_a0() << std::endl;
102 out << " * z_0 : "<< track.z() <<" "<< track.r_z() << std::endl;
103 out << " * theta : "<< track.theta()<<" "<< track.r_theta()<< std::endl;
104 out << " * phi : "<< track.phi() <<" "<< track.r_phi() << std::endl;
105 out << " * q/R : "<< track.invR() <<" "<< track.r_invR() << std::endl;
106 out << " * Charge: "<< track.Charge <<" Mass: "<<track.m_mass<< std::endl;
107 out << "-> with ref ErrorMatrix: " << std::endl;out.precision(5);
108 out << track.refCovar[0] <<std::endl;
109 out << track.refCovar[1] <<", "<<track.refCovar[2] <<std::endl;
110 out << track.refCovar[3] <<", "<<track.refCovar[4] <<", "<<track.refCovar[5] <<std::endl;
111 out << track.refCovar[6] <<", "<<track.refCovar[7] <<", "<<track.refCovar[8] <<", "
112 << track.refCovar[9]<<std::endl;
113 out << track.refCovar[10]<<", "<<track.refCovar[11]<<", "<<track.refCovar[12]<<", "
114 << track.refCovar[13]<<", "<<track.refCovar[14]<<std::endl;
115 out << "-> and iteration WeightMatrix: " << std::endl;
116 out << track.WgtM[0] <<std::endl;
117 out << track.WgtM[1] <<", "<<track.WgtM[2] <<std::endl;
118 out << track.WgtM[3] <<", "<<track.WgtM[4] <<", "<<track.WgtM[5] <<std::endl;
119 out << track.WgtM[6] <<", "<<track.WgtM[7] <<", "<<track.WgtM[8] <<", "
120 << track.WgtM[9]<<std::endl;
121 out << track.WgtM[10]<<", "<<track.WgtM[11]<<", "<<track.WgtM[12]<<", "
122 << track.WgtM[13]<<", "<<track.WgtM[14]<<std::endl;
123 return out;
124 }
125
126
127
128 // cppcheck-suppress uninitMemberVar; large ader array not initialized
130 { vk_fitterControl = std::make_unique<VKalVrtControl>(FitControl); }
131
132 // cppcheck-suppress uninitMemberVar; large ader array not initialized
136 {
137 for(int i=0; i<3; i++) apriorV[i]=refV[i]=iniV[i]=fitV[i]=cnstV[i]=fitMom[i]=0.;
138 for(int i=0; i<6; i++) apriorVWGT[i]=0.;
139 for(int i=0; i<21; i++) fitCovXYZMom[i]=savedVrtMomCov[i]=0.;
142 Chi2=0.;
143 truncatedStep=false;
144 existFullCov=0;
145 }
146
147 VKVertex::~VKVertex() = default;
148
149 void VKVertex::setRefV(double v[3]) noexcept { std::copy(v, v+3, refV);}
150 void VKVertex::setRefIterV(double v[]) noexcept { std::copy(v, v+3, refIterV); }
151 void VKVertex::setIniV(double v[3]) noexcept { std::copy(v, v+3, iniV); }
152 void VKVertex::setFitV(double v[3]) noexcept { std::copy(v, v+3, fitV); }
153 void VKVertex::setCnstV(double v[3]) noexcept { std::copy(v, v+3, cnstV);}
154
155
156 VKVertex::VKVertex(const VKVertex & src): //copy constructor
157 Chi2(src.Chi2), // vertex Chi2
158 useApriorVertex(src.useApriorVertex), //for a priory vertex position knowledge usage
159 passNearVertex(src.passNearVertex), // needed for "passing near vertex" constraint
160 passWithTrkCov(src.passWithTrkCov), // Vertex, CovVertex, Charge and derivatives
161 FVC(src.FVC),
162 TrackList(0),
164 {
165 for( int i=0; i<6; i++) {
166 fitVcov[i] =src.fitVcov[i]; // range[0:5]
167 apriorVWGT[i] =src.apriorVWGT[i];
168 if(i>=3) continue;
169 fitV[i] =src.fitV[i]; // range[0:2]
170 iniV[i] =src.iniV[i];
171 cnstV[i] =src.cnstV[i];
172 refIterV[i]=src.refIterV[i];
173 refV[i] =src.refV[i];
174 apriorV[i] =src.apriorV[i];
175 fitMom[i] = src.fitMom[i];
176 }
177 for( int i=0; i<21; i++){
178 savedVrtMomCov[i]=src.savedVrtMomCov[i];
179 fitCovXYZMom[i]=src.fitCovXYZMom[i];
180 }
181 nextCascadeVrt = nullptr;
182 truncatedStep = src.truncatedStep;
183 existFullCov = src.existFullCov;
184
185 int NTrack=src.TrackList.size();
186 int FULLSIZE=3*vkalNTrkM+3;
187 int SIZE=3*NTrack+3;
188 for(int i=0; i<SIZE; i++) for(int j=0; j<SIZE; j++) { ader[i*FULLSIZE+j]=src.ader[i*FULLSIZE+j]; }
189 //----- Creation of track and constraint copies
190 for( int i=0; i<NTrack; i++) TrackList.emplace_back( new VKTrack(*(src.TrackList[i])) );
191 ConstraintList.reserve(src.ConstraintList.size());
192 for( int ic=0; ic<(int)src.ConstraintList.size(); ic++){
193 ConstraintList.emplace_back(src.ConstraintList[ic]->clone());
194 }
195 vk_fitterControl = std::make_unique<VKalVrtControl>(*src.vk_fitterControl);
196 }
197
198 VKVertex& VKVertex::operator= (const VKVertex & src) //Assignment operator
199 {
200 if (this!=&src){
201 vk_fitterControl.reset();
202 vk_fitterControl=std::make_unique<VKalVrtControl>(*(src.vk_fitterControl));
203 Chi2=src.Chi2; // vertex Chi2
204 useApriorVertex=src.useApriorVertex; //for a priory vertex position knowledge usage
205 passNearVertex=src.passNearVertex; // needed for "passing near vertex" constraint
206 passWithTrkCov=src.passWithTrkCov; // Vertex, CovVertex, Charge and derivatives
207 FVC=src.FVC;
208 for( int i=0; i<6; i++) {
209 fitVcov[i] =src.fitVcov[i]; // range[0:5]
210 apriorVWGT[i] =src.apriorVWGT[i];
211 if(i>=3) continue;
212 fitV[i] =src.fitV[i]; // range[0:2]
213 iniV[i] =src.iniV[i];
214 cnstV[i] =src.cnstV[i];
215 refIterV[i]=src.refIterV[i];
216 refV[i] =src.refV[i];
217 apriorV[i] =src.apriorV[i];
218 fitMom[i] = src.fitMom[i];
219 }
220 for( int i=0; i<21; i++){
221 savedVrtMomCov[i]=src.savedVrtMomCov[i];
222 fitCovXYZMom[i]=src.fitCovXYZMom[i];
223 }
224 nextCascadeVrt = nullptr;
225 truncatedStep = src.truncatedStep;
226 existFullCov = src.existFullCov;
227
228 int NTrack=src.TrackList.size();
229 int FULLSIZE=3*vkalNTrkM+3;
230 int SIZE=3*NTrack+3;
231 for(int i=0; i<SIZE; i++) for(int j=0; j<SIZE; j++) { ader[i*FULLSIZE+j]=src.ader[i*FULLSIZE+j]; }
232 //----- Creation of track and constraint copies
233 TrackList.clear();
234 tmpArr.clear();
235 ConstraintList.clear();
236 for( int i=0; i<NTrack; i++) TrackList.emplace_back( new VKTrack(*(src.TrackList[i])) );
237 ConstraintList.reserve(src.ConstraintList.size());
238 for( int ic=0; ic<(int)src.ConstraintList.size(); ic++){
239 ConstraintList.emplace_back(src.ConstraintList[ic]->clone());
240 }
241
242 std::copy (std::begin(src.T), std::end(src.T), std::begin(T));
243 std::copy (std::begin(src.wa), std::end(src.wa), std::begin(wa));
244 std::copy (std::begin(src.dxyz0), std::end(src.dxyz0), std::begin(dxyz0));
245 }
246 return *this;
247 }
248
249
251 {
252 if (Iter<3) Iter=3;
253 if (Iter>100) Iter=100;
254 vk_forcft.IterationNumber = Iter;
255 }
256
258 {
259 if (Prec<1.e-5) Prec=1.e-5;
260 if (Prec>1.e-1) Prec=1.e-1;
261 vk_forcft.IterationPrecision = Prec;
262 }
263
265 {
266 if (Scale<0.01) Scale=0.01;
267 if (Scale>100.) Scale=100.;
268 vk_forcft.RobustScale = Scale;
269 }
270
272 {
273 if (Rob<0) Rob=0;
274 if (Rob>7) Rob=7;
275 vk_forcft.irob = Rob;
276 }
277 void VKalVrtControl::setMassCnstData(int Ntrk, double Mass){ // Define global mass constraint. It must be first
278 double sumM(0.);
279 for(int it=0; it<Ntrk; it++) sumM += vk_forcft.wm[it]; //sum of particle masses
280 if(sumM<Mass) {
281 vk_forcft.wmfit[0]=Mass;
282 for(int it=0; it<Ntrk; it++) vk_forcft.indtrkmc[0][it]=1; //Set participating particles
283 vk_forcft.nmcnst=1;
284 }
285 vk_forcft.useMassCnst = 1;
286 }
287 void VKalVrtControl::setMassCnstData(int NtrkTot, const std::vector<int> & Index, double Mass){
288 double sumM(0.);
289 int Ntrk=std::min((int)Index.size(),NtrkTot);
290 for(int it=0; it<Ntrk; it++) sumM += vk_forcft.wm[Index[it]]; //sum of particle masses
291 if(sumM<Mass) {
292 vk_forcft.wmfit[vk_forcft.nmcnst]=Mass;
293 for(int it=0; it<Ntrk; it++) vk_forcft.indtrkmc[vk_forcft.nmcnst][Index[it]-1]=1; //Set participating particles
294 vk_forcft.nmcnst++;
295 }
296 vk_forcft.useMassCnst = 1;
297 }
298
299 void VKalVrtControl::setUsePhiCnst() { vk_forcft.usePhiCnst = 1;}
300 void VKalVrtControl::setUsePlaneCnst(double a, double b, double c, double d) {
301 if(a+b+c+d == 0.){ vk_forcft.usePlaneCnst = 0;
302 }else{ vk_forcft.usePlaneCnst = 1; }
303 vk_forcft.Ap = a; vk_forcft.Bp = b; vk_forcft.Cp = c; vk_forcft.Dp = d;
304 }
305 void VKalVrtControl::setUseRadiusCnst(double R, double RefP[2]) {
306 if(R == 0.){ vk_forcft.useRadiusCnst = 0;
307 }else{ vk_forcft.useRadiusCnst = 1; }
308 vk_forcft.RC = R;
309 vk_forcft.radiusRefP[0]=RefP[0];
310 vk_forcft.radiusRefP[1]=RefP[1];
311 }
312 void VKalVrtControl::setUseThetaCnst() { vk_forcft.useThetaCnst = 1;}
313 void VKalVrtControl::setUseAprioriVrt(){ vk_forcft.useAprioriVrt = 1;}
314 void VKalVrtControl::setUsePointingCnst(int iType = 1 ) { vk_forcft.usePointingCnst = iType<2 ? 1 : 2 ;}
315 void VKalVrtControl::setUsePassNear(int iType = 1 ) { vk_forcft.usePassNear = iType<2 ? 1 : 2 ;}
320 void VKalVrtControl::renewFullCovariance(double * newarray) {
321 m_fullCovariance.reset(newarray);
322 }
323
324} /* End of namespace */
#define vkalNTrkM
Definition CommonPars.h:22
IndexedConstituentUserInfo::Index Index
static Double_t a
double WgtM_save[15]
double refCovar[15]
void setReference(const double[], const double[])
void setCurrent(const double[], const double[])
VKVertex * m_originVertex
VKTrack(long int, const double[], const double[], VKVertex *, double)
double savedVrtMomCov[21]
VKVertex(const VKalVrtControl &)
void setRefIterV(double v[]) noexcept
void setFitV(double v[3]) noexcept
void setCnstV(double v[3]) noexcept
std::vector< VKVertex * > includedVrt
VKVertex & operator=(const VKVertex &src)
void setIniV(double v[3]) noexcept
std::vector< std::unique_ptr< VKTrack > > TrackList
double ader[(3 *vkalNTrkM+3) *(3 *vkalNTrkM+3)]
std::vector< std::unique_ptr< VKConstraintBase > > ConstraintList
double fitCovXYZMom[21]
void setRefV(double v[3]) noexcept
std::vector< std::unique_ptr< TWRK > > tmpArr
VKVertex * nextCascadeVrt
std::unique_ptr< VKalVrtControl > vk_fitterControl
const basePropagator * vk_objProp
const addrPropagator vk_funcProp
VKalVrtControlBase(baseMagFld *, const addrMagHandler, const basePropagator *, const addrPropagator, IVKalState *istate=nullptr)
const addrMagHandler vk_funcMagFld
void renewFullCovariance(double *)
void renewCascadeEvent(CascadeEvent *)
void setIterationPrec(double Prec)
void setMassCnstData(int Ntrk, double Mass)
VKalVrtControl(const VKalVrtControlBase &)
void setRobustScale(double Scale)
CascadeEvent * m_cascadeEvent
void setUsePlaneCnst(double a, double b, double c, double d)
void setUseRadiusCnst(double R, double RefP[2])
std::unique_ptr< double[]> m_fullCovariance
void Scale(TH1 *h, double d=1)
std::string base
Definition hcg.cxx:81
Ensure that the ATLAS eigen extensions are properly loaded.
void(* addrMagHandler)(double, double, double, double &, double &, double &)
Definition VKalVrtBMag.h:23
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
MsgStream & operator<<(MsgStream &sl, const AlignModule &alignModule)
overload of << operator for MsgStream for debug output
@ v
Definition ParamDefs.h:78
void(* addrPropagator)(long int, long int, double *, double *, double *, double *, double *, double *)
Definition Propagator.h:23