BOSS 7.0.6
BESIII Offline Software System
Loading...
Searching...
No Matches
TRungeFitter.cxx
Go to the documentation of this file.
1//-----------------------------------------------------------------------------
2// $Id: TRungeFitter.cxx,v 1.34 2012/05/28 05:16:29 maoh Exp $
3//-----------------------------------------------------------------------------
4// Filename : TRungeFitter.cc
5// Section : Tracking
6// Owner : Kenji Inami
7// Email : [email protected]
8//-----------------------------------------------------------------------------
9// Description : A class to fit a TTrackBase object to a Runge Kutta track
10// See http://bsunsrv1.kek.jp/~yiwasaki/tracking/
11//-----------------------------------------------------------------------------
12
13#ifdef TRKRECO_DEBUG_DETAIL
14#ifndef TRKRECO_DEBUG
15#define TRKRECO_DEBUG
16#endif
17#endif
18#include <float.h>
19#include <iostream>
21#include "TrkReco/TRunge.h"
22#include "GaudiKernel/Service.h"
23#include "GaudiKernel/ISvcLocator.h"
24#include "GaudiKernel/SvcFactory.h"
25#define HEP_SHORT_NAMES
26//#include "panther/panther.h"
27//#include MDC_H
28#include "MdcTables/MdcTables.h"
29#include "CLHEP/Matrix/Vector.h"
30#include "CLHEP/Matrix/SymMatrix.h"
31#include "CLHEP/Matrix/Matrix.h"
32#include "TrkReco/TMLink.h"
33#include "TrkReco/TTrackBase.h"
36#include "GaudiKernel/StatusCode.h"
37#include "GaudiKernel/IInterface.h"
38#include "GaudiKernel/Kernel.h"
39#include "GaudiKernel/Service.h"
40#include "GaudiKernel/ISvcLocator.h"
41#include "GaudiKernel/SvcFactory.h"
42#include "GaudiKernel/IDataProviderSvc.h"
43#include "GaudiKernel/Bootstrap.h"
44#include "GaudiKernel/MsgStream.h"
45#include "GaudiKernel/SmartDataPtr.h"
46#include "GaudiKernel/IMessageSvc.h"
47#include "MdcTables/MdcTables.h"
48#include "CLHEP/Matrix/Vector.h"
49#include "CLHEP/Matrix/SymMatrix.h"
50#include "CLHEP/Matrix/Matrix.h"
51#include "TrkReco/TMLink.h"
52#include "TrkReco/TTrackBase.h"
53//#include "TrkReco/THistAddItem.h"
54#include "GaudiKernel/ISvcLocator.h"
57#include "GaudiKernel/Service.h"
58#include "GaudiKernel/ISvcLocator.h"
59#include "GaudiKernel/SvcFactory.h"
60#include "AIDA/IHistogram1D.h"
61
62#include "AIDA/IHistogram2D.h"
63#include "G4Geo/MdcG4Geo.h"
64#include "G4Geo/BesG4Geo.h"
65#include "G4Material.hh"
66#include "G4Tubs.hh"
67
68#include <time.h>
69
70using CLHEP::HepVector;
71using CLHEP::HepMatrix;
72using CLHEP::HepSymMatrix;
73/*
74extern "C"
75void
76calcdc_driftdist_(int *,
77 int *,
78 int *,
79 float[3],
80 float[3],
81 float *,
82 float *,
83 float *);
84
85extern "C"
86void
87calcdc_tof2_(int *, float *, float *, float *);
88*/
89
90TRungeFitter::TRungeFitter(const std::string& name)
91 : TMFitter(name),
92 _sag(true),_propagation(1),_tof(false){
93
94 StatusCode scmgn = Gaudi::svcLocator()->service("MagneticFieldSvc",m_pmgnIMF);
95 if(scmgn!=StatusCode::SUCCESS) {
96 std::cout<< __FILE__<<" Unable to open Magnetic field service"<<std::endl;
97 }
98}
99TRungeFitter::TRungeFitter(const std::string& name,
100 bool m_sag,int m_prop,bool m_tof)
101 : TMFitter(name),
102 _sag(m_sag),_propagation(m_prop),_tof(m_tof){
103
104 StatusCode scmgn = Gaudi::svcLocator()->service("MagneticFieldSvc",m_pmgnIMF);
105 if(scmgn!=StatusCode::SUCCESS) {
106 std::cout<< "Unable to open Magnetic field service"<<std::endl;
107 }
108}
109
111}
112
113void TRungeFitter::sag(bool _in){
114 _sag = _in;
115}
117 _propagation = _in;
118}
119void TRungeFitter::tof(bool _in){
120 _tof = _in;
121}
123 return fit(tb,0,-1);
124}
125
126int TRungeFitter::fit(TTrackBase& tb, float t0Offset, int layer) const{
127// Argument layer is used for calibration interface in which doca is calculated without measured hit. liucy
128 IMdcCalibFunSvc* l_mdcCalFunSvc;
129 Gaudi::svcLocator() -> service("MdcCalibFunSvc", l_mdcCalFunSvc);
130 //std::cout<<"TRungeFitter::fit start"<<std::endl;
131 //...Type check...
132 if(tb.objectType() != Runge) return TFitUnavailable;
133 TRunge& t = (TRunge&) tb;
134 //...Already fitted ?...
135 if(t.fitted()&&layer==-1) return TFitAlreadyFitted;
136 //...Count # of hits...
137 const AList<TMLink>& cores = t.cores();
138 unsigned nCores = cores.length();
139 unsigned nStereoCores = NStereoHits(cores);
140 TMLink* last=NULL;
141 int lyr=0;
142 int layerid=0;
143 for(int i=0;i<nCores;i++){
144 layerid=(*cores[i]->hit()).wire()->layerId();
145 if(layerid>=lyr){
146 lyr=layerid;
147 last=cores[i];
148 }
149 }
150
151 //...Check # of hits...
152 // if ((nStereoCores < 2) || (nCores - nStereoCores < 3))
153 // return TFitErrorFewHits;
154
155 //...Move pivot to ORIGIN...
156 const HepPoint3D pivot_bak = t.pivot();
157 t.pivot(ORIGIN);
158
159 //...Setup...
160 Vector a(5),da(5);
161 a = t.a();
162 Vector ainitial(a);
163 Vector dDda(5);
164 Vector dchi2da(5);
165 SymMatrix d2chi2d2a(5,0);
166 const SymMatrix zero5(5,0);
167 double chi2;
168 double chi2Old = 0;
169 for(int i=0;i<t.links().length();i++){
170 chi2Old=chi2Old+t.links()[i]->pull();
171 }
172 chi2Old=DBL_MAX;
173 int err = 0;
174
175 double factor = 0.1;
176 Vector beta(5);
177 SymMatrix alpha(5,0);
178 SymMatrix alpha2(5,0);
179
180 double (*d)[5]=new double[nCores][5];
181 // double (*d2)[5]=new double[nCores][5];
182 Vector ea(5);
183
184 float tof;
185 HepVector3D p;
186 unsigned i;
187
188 double distance;
189 double eDistance;
190
191 // ea... init
192 const double ea_init=0.000001;
193 ea[0]=ea_init; //dr
194 ea[1]=ea_init; //phi0
195 ea[2]=ea_init; //kappa
196 ea[3]=ea_init; //dz
197 ea[4]=ea_init; //tanl
198
199 //std::cout<<"TRF ::"<<a[0]<<","<<a[1]<<","<<a[2]<<","<<a[3]<<","<<a[4]<<std::endl;
200
201 //long int lclock0=clock()/1000;
202 //long int lclock=lclock0;
203
204 Vector def_a(t.a());
205 Vector ta(def_a);
206
207 //...Fitting loop...
208 unsigned nTrial = 0;
209 while(nTrial < 100){
210
211 //...Set up...
212 chi2 = 0;
213 for (unsigned j=0;j<5;j++) dchi2da[j]=0;
214 d2chi2d2a=zero5;
215
216 def_a=t.a();
217
218 //#### loop for shifted helix parameter set ####
219 for(unsigned j=0;j<5;j++){
220 ta=def_a;
221 ta[j]+=ea[j];
222 t.a(ta);
223 //...Loop with hits...
224 i=0;
225 //std::cout<<"TRF:: cores="<<cores.length()<<std::endl;
226 while(TMLink * l = cores[i++]){
227 //...Cal. closest points...
228 t.approach(* l,tof,p,_BesBeamPipeMaterials[0],_sag);
229 const HepPoint3D& onTrack=l->positionOnTrack();
230 const HepPoint3D& onWire=l->positionOnWire();
231 d[i-1][j]=(onTrack-onWire).mag();
232 //std::cout<<"TRF:: i="<<i<<std::endl;
233 }//end of loop with hits
234 //lclock=clock()/1000;
235 //std::cout<<"TRF clock="<<lclock-lclock0<<std::endl;
236 //lclock0=lclock;
237 }
238 /*
239 for(int j=0;j<5;j++){
240 ta=def_a;
241 ta[j]-=ea[j];
242 t.a(ta);
243 //...Loop with hits...
244 float tof_dummy;
245 Vector3 p_dummy;
246 unsigned i=0;
247 while(TMLink* l = cores[i++]){
248 //...Cal. closest points...
249 t.approach(*l,tof_dummy,p_dummy,_sag);
250 const HepPoint3D& onTrack=l->positionOnTrack();
251 const HepPoint3D& onWire=l->positionOnWire();
252 d2[i-1][j]=(onTrack-onWire).mag();
253 }//end of loop with hits
254 }
255 */
256 t.a(def_a);
257
258 //#### original parameter set and calc. chi2 ####
259 //...Loop with hits...
260 i=0;
261 while(TMLink * l = cores[i++]){
262 const TMDCWireHit& h = * l->hit();
263 //...Cal. closest points...
264 if(t.approach(* l,tof,p,_BesBeamPipeMaterials[0],_sag)<0){
265 std::cout<<"TrkReco:TRF:: bad wire"<<std::endl;
266 continue;
267 }
268 int layerId=h.wire()->layerId();
269 const HepPoint3D& onTrack=l->positionOnTrack();
270 const HepPoint3D& onWire=l->positionOnWire();
271 unsigned leftRight = WireHitRight;
272 if (onWire.cross(onTrack).z() < 0.) leftRight = WireHitLeft;
273
274 //...Obtain drift distance and its error...
275 if ((t0Offset == 0.) && (_propagation==0) && (! _tof)) {
276 //...No correction...
277 distance = l->drift(leftRight);
278 eDistance = h.dDrift(leftRight);
279 }else{
280 //...T0 and propagation corrections...
281 int wire = h.wire()->id();
282 int wirelocal=h.wire()->localId();
283 int side = leftRight;
284 if (side==0) side = 0;
285 float tp[3] = {p.x(),p.y(),p.z()};
286 float x[3] = {onWire.x(), onWire.y(), onWire.z()};
287 double tprop = l_mdcCalFunSvc->getTprop(layerId,onWire.z()*10.);
288 double timeWalk = l_mdcCalFunSvc->getTimeWalk(layerId, h.reccdc()->adc);
289 double T0 = l_mdcCalFunSvc->getT0(layerId,wirelocal);
290 double drifttime = h.reccdc()->tdc - tof - tprop - timeWalk -T0;
291 l->setDriftTime(drifttime);
292 float dist;
293 float edist;
294 int prop = _propagation;
295 const double x0 = t.helix().pivot().x();
296 const double y0 = t.helix().pivot().y();
297 Hep3Vector pivot_helix(x0,y0,0);
298 const double dr = t.helix().a()[0];
299 const double phi0 = t.helix().a()[1];
300 const double kappa = t.helix().a()[2];
301 const double dz = t.helix().a()[3];
302 const double tanl = t.helix().a()[4];
303
304 const double Bz = -1000*(m_pmgnIMF->getReferField());
305 const double alpha = 333.564095 / Bz;
306
307 const double cox = x0 + dr*cos(phi0)- alpha*cos(phi0)/kappa;
308 const double coy = y0 + dr*sin(phi0) - alpha*sin(phi0)/kappa;
309 HepPoint3D dir(onTrack.y()-coy,cox-onTrack.x(),0);
310 double pos_phi=onWire.phi();
311 double dir_phi=dir.phi();
312 while(pos_phi>pi){pos_phi-=pi;}
313 while(pos_phi<0){pos_phi+=pi;}
314 while(dir_phi>pi){dir_phi-=pi;}
315 while(dir_phi<0){dir_phi+=pi;}
316 double entrangle=dir_phi-pos_phi;
317 while(entrangle>pi/2)entrangle-=pi;
318 while(entrangle<(-1)*pi/2)entrangle+=pi;
319 dist = l_mdcCalFunSvc->driftTimeToDist(drifttime,layerId, wirelocal, side, entrangle);
320 dist= dist/10;//mm->cmo
321 if(layer==-1||layerId!=layer){
322 edist = l_mdcCalFunSvc->getSigma(layerId, side, dist, entrangle,tanl,onWire.z()*10,h.reccdc()->adc); }
323 else if(layerId==layer){edist=99999999;}
324 edist = edist/10;
325 distance = (double) dist;
326 eDistance = (double) edist;
327
328 }
329 double eDistance2=eDistance*eDistance;
330
331 //...Residual...
332 const double d0 = (onTrack-onWire).mag();
333 //if(d0>2) std::cout<<"TRF:: strange dist. d0="<<d0<<" x="<<distance
334 // <<" ex="<<eDistance<<std::endl;
335 double dDistance = d0 - distance;
336
337 //...dDda...
338 for(int j=0;j<5;j++){
339 dDda[j]=(d[i-1][j]-d0)/ea[j];
340 //if(dDda[j]==0) std::cout<<"TRF:: dDda=0 j="<<j<<" ea="<<ea[j]<<std::endl;
341 }
342 // for(int j=0;j<5;j++) dDda[j]=(d[i-1][j]-d2[i-1][j])/ea[j]/2.;
343 //...Chi2 related...
344 dchi2da += (dDistance / eDistance2) * dDda;
345 d2chi2d2a += vT_times_v(dDda) / eDistance2;
346 double pChi2 = dDistance * dDistance / eDistance2;
347 chi2 += pChi2;
348 //if(!(pChi2<0 || pChi2>=0)){
349 // std::cout<<"TRF:: pChi2="<<pChi2<<" X="<<d0<<" fx="<<distance
350 // <<" ex="<<eDistance<<std::endl;
351 //}
352
353 //...Store results...
354 if(layer==-1){
355 l->update(onTrack, onWire, leftRight, pChi2);
356 l->drift(distance,0);
357 l->drift(distance,1);
358 l->dDrift(eDistance,0);
359 l->dDrift(eDistance,1);
360 }
361
362 else if(layerId==layer){
363 l->distance((onTrack-onWire).mag());
364 }
365 }//end of loop with hits
366
367 //...Check condition...
368 double change = chi2Old - chi2;
369 if (0 <= change && change < 0.01) break;
370 if (change < 0.) {
371 factor *= 100;
372 a += da; //recover
373 if(-0.01 < change){
374 d2chi2d2a=alpha;
375 chi2=chi2Old;
376 break;
377 }
378 }else if(change > 0.){
379 chi2Old = chi2;
380 factor *= 0.1;
381 alpha=d2chi2d2a;
382 beta=dchi2da;
383 }else if(nTrial==0){
384 chi2Old = chi2;
385 factor *= 0.1;
386 alpha=d2chi2d2a;
387 beta=dchi2da;
388 }else{
389 std::cout<<"TrkReco:TRF:: bad chi2 = "<<chi2<<std::endl;
390 err=TFitFailed;
391// break; // protection for nan
392 return err;
393 }
394 alpha2=alpha;
395 for(unsigned j=0;j<5;j++) alpha2[j][j]*=(1+factor);
396 //...Cal. helix parameters for next loop...
397 da = solve(alpha2,beta);
398
399 //lclock=clock()/1000;
400 //std::cout<<"TRF "<<nTrial<<": "
401 // <<"cl="<<lclock-lclock0<<": "
402 // <<a[0]<<","<<a[1]<<","<<a[2]<<","<<a[3]<<","<<a[4]<<" "
403 // <<chi2<<"/"<<nCores<<":"<<factor
404 // <<" :da "<<da[0]<<","<<da[1]<<","<<da[2]<<","<<da[3]<<","<<da[4]<<std::endl;
405 //lclock0=lclock;
406
407 a -= da;
408 t.a(a);
409 //ea = 0.0001*da;
410 for(i=0;i<5;i++){
411 ea[i]=0.0001*abs(da[i]);
412 if(ea[i]>ea_init) ea[i]=ea_init;
413 if(ea[i]<1.0e-10) ea[i]=1.0e-10;
414 }
415 ++nTrial;
416 }
417 //std::cout<<"TRungeFitter:: nTrial="<<nTrial<<std::endl;
418
419 //...Cal. error matrix...
420 SymMatrix Ea(5,0);
421 unsigned dim;
422 dim=5;
423 Ea = d2chi2d2a.inverse(err);
424 if(nCores){
425 t.approach(*last,tof,p,_BesBeamPipeMaterials[0],_sag);}
426 // consider the material effect beam pipe.
427 double y_[6]={0,0,0,0,0,0};
428 t.SetFirst(y_);//obtain the momentum of the first layer
429 int lmass=0;
430 innerwall(t,lmass,y_);// consider the material layer by layer
431 int nsign=a[2]/abs(a[2]);
432 // Update the helix parameter (momentum part.)
433 a[4]=y_[5]/sqrt(y_[4]*y_[4]+y_[3]*y_[3]);
434 a[2]=nsign*1/sqrt(y_[4]*y_[4]+y_[3]*y_[3]);
435 //...Store information...
436 if(! err&&layer==-1){
437 t.a(a);
438 t.Ea(Ea);
439 t._fitted = true;
440 }else if(! err&&layer!=-1){
441 t.a(ainitial);
442// t.Ea(Ea);
443 t._fitted = true;
444 }else{
445 err = TFitFailed;
446 }
447
448 t._ndf = nCores - dim;
449 t._chi2 = chi2;
450
451 //...Recover pivot...
452 t.pivot(pivot_bak);
453
454 delete [] d;
455 // delete [] d2;
456
457 return err;
458}
460 int i(0);
461 double Z(0.),A(0.),Ionization(0.),Density(0.),Radlen(0.);
462
463 G4LogicalVolume *logicalMdc = 0;
464 MdcG4Geo* aMdcG4Geo = new MdcG4Geo();
465 logicalMdc = aMdcG4Geo->GetTopVolume();
466
467 /// mdcgas
468 G4Material* mdcMaterial = logicalMdc->GetMaterial();
469
470 for(i=0; i<mdcMaterial->GetElementVector()->size(); i++){
471 Z += (mdcMaterial->GetElement(i)->GetZ())*
472 (mdcMaterial->GetFractionVector()[i]);
473 A += (mdcMaterial->GetElement(i)->GetA())*
474 (mdcMaterial->GetFractionVector()[i]);
475 }
476 Ionization = mdcMaterial->GetIonisation()->GetMeanExcitationEnergy();
477 Density = mdcMaterial->GetDensity()/(g/cm3);
478 Radlen = mdcMaterial->GetRadlen();
479 RkFitMaterial FitMdcMaterial(Z,A/g/mole,Ionization/eV,Density,Radlen/10.);
480 cout<<"mdcgas: Z: "<<Z<<" A: "<<(A/(g/mole))<<" Ionization: "<<(Ionization/eV)<<" Density: "<<Density<<" Radlen: "<<Radlen<<endl;
481 _BesBeamPipeMaterials.push_back(FitMdcMaterial);
482 //RkFitTrack::mdcGasRadlen_ = Radlen/10.;
483
484 /// inner wall aluminium
485 G4LogicalVolume* innerwallVolume = const_cast<G4LogicalVolume*>(GDMLProcessor::GetInstance()->GetLogicalVolume("logicalMdcSegment2"));
486 G4Material* innerwallMaterial = innerwallVolume->GetMaterial();
487 G4Tubs* innerwallTub = dynamic_cast<G4Tubs*>(innerwallVolume->GetSolid());
488
489 Z = 0.;
490 A = 0.;
491 for(i=0; i<innerwallMaterial->GetElementVector()->size(); i++){
492 Z += (innerwallMaterial->GetElement(i)->GetZ())*
493 (innerwallMaterial->GetFractionVector()[i]);
494 A += (innerwallMaterial->GetElement(i)->GetA())*
495 (innerwallMaterial->GetFractionVector()[i]);
496 }
497
498 Ionization = innerwallMaterial->GetIonisation()->GetMeanExcitationEnergy();
499 Density = innerwallMaterial->GetDensity()/(g/cm3);
500 Radlen = innerwallMaterial->GetRadlen();
501 cout<<"Mdc innerwall, Al: Z: "<<Z<<" A: "<<(A/(g/mole))<<" Ionization: "<<(Ionization/eV)<<" Density: "<<Density<<" Radlen: "<<Radlen<<endl;
502 RkFitMaterial FitInnerwallMaterial(Z,A/g/mole,Ionization/eV,Density,Radlen/10.);
503 _BesBeamPipeMaterials.push_back(FitInnerwallMaterial);
504
505 ///////////////////////////////////////////////////////////////////////////////////////////////////
506 G4LogicalVolume *logicalBes = 0;
507 BesG4Geo* aBesG4Geo = new BesG4Geo();
508 logicalBes = aBesG4Geo->GetTopVolume();
509
510 /// air
511 G4LogicalVolume* logicalAirVolume = const_cast<G4LogicalVolume*>(GDMLProcessor::GetInstance()->GetLogicalVolume("logicalWorld"));
512 G4Material* airMaterial = logicalAirVolume->GetMaterial();
513 Z = 0.;
514 A = 0.;
515 for(i=0; i<airMaterial->GetElementVector()->size(); i++){
516 Z += (airMaterial->GetElement(i)->GetZ())*
517 (airMaterial->GetFractionVector()[i]);
518 A += (airMaterial->GetElement(i)->GetA())*
519 (airMaterial->GetFractionVector()[i]);
520 }
521 Ionization = airMaterial->GetIonisation()->GetMeanExcitationEnergy();
522 Density = airMaterial->GetDensity()/(g/cm3);
523 Radlen = airMaterial->GetRadlen();
524 cout<<"air: Z: "<<Z<<" A: "<<(A/(g/mole))<<" Ionization: "<<(Ionization/eV)<<" Density: "<<Density<<" Radlen: "<<Radlen<<endl;
525 RkFitMaterial FitAirMaterial(Z,A/g/mole,Ionization/eV,Density,Radlen/10.);
526 _BesBeamPipeMaterials.push_back(FitAirMaterial);
527
528 /// outer beryllium pipe
529 G4LogicalVolume* logicalOuterBeVolume = const_cast<G4LogicalVolume*>(GDMLProcessor::GetInstance()->GetLogicalVolume("logicalouterBe"));
530 G4Material* outerBeMaterial = logicalOuterBeVolume->GetMaterial();
531
532 G4Tubs* outerBeTub = dynamic_cast<G4Tubs*>(logicalOuterBeVolume->GetSolid());
533 Z = 0.;
534 A = 0.;
535 for(i=0; i<outerBeMaterial->GetElementVector()->size(); i++){
536 Z += (outerBeMaterial->GetElement(i)->GetZ())*
537 (outerBeMaterial->GetFractionVector()[i]);
538 A += (outerBeMaterial->GetElement(i)->GetA())*
539 (outerBeMaterial->GetFractionVector()[i]);
540 }
541 Ionization = outerBeMaterial->GetIonisation()->GetMeanExcitationEnergy();
542 Density = outerBeMaterial->GetDensity()/(g/cm3);
543 Radlen = outerBeMaterial->GetRadlen();
544 cout<<"outer beryllium: Z: "<<Z<<" A: "<<(A/(g/mole))<<" Ionization: "<<(Ionization/eV)<<" Density: "<<Density<<" Radlen: "<<Radlen<<endl;
545 RkFitMaterial FitOuterBeMaterial(Z,A/g/mole,Ionization/eV,Density,Radlen/10.);
546 _BesBeamPipeMaterials.push_back(FitOuterBeMaterial);
547
548 /// cooling oil
549 G4LogicalVolume* logicalOilLayerVolume = const_cast<G4LogicalVolume*>(GDMLProcessor::GetInstance()->GetLogicalVolume("logicaloilLayer"));
550 G4Material* oilLayerMaterial = logicalOilLayerVolume->GetMaterial();
551 G4Tubs* oilLayerTub = dynamic_cast<G4Tubs*>(logicalOilLayerVolume->GetSolid());
552
553 Z = 0.;
554 A = 0.;
555 for(i=0; i<oilLayerMaterial->GetElementVector()->size(); i++){
556 Z += (oilLayerMaterial->GetElement(i)->GetZ())*
557 (oilLayerMaterial->GetFractionVector()[i]);
558 A += (oilLayerMaterial->GetElement(i)->GetA())*
559 (oilLayerMaterial->GetFractionVector()[i]);
560 }
561 Ionization = oilLayerMaterial->GetIonisation()->GetMeanExcitationEnergy();
562 Density = oilLayerMaterial->GetDensity()/(g/cm3);
563 Radlen = oilLayerMaterial->GetRadlen();
564 cout<<"cooling oil: Z: "<<Z<<" A: "<<(A/(g/mole))<<" Ionization: "<<(Ionization/eV)<<" Density: "<<Density<<" Radlen: "<<Radlen<<endl;
565 RkFitMaterial FitOilLayerMaterial(Z,A/g/mole,Ionization/eV,Density,Radlen/10.);
566 _BesBeamPipeMaterials.push_back(FitOilLayerMaterial);
567
568
569 /// inner beryllium pipe
570 G4LogicalVolume* logicalInnerBeVolume = const_cast<G4LogicalVolume*>(GDMLProcessor::GetInstance()->GetLogicalVolume("logicalinnerBe"));
571
572 G4Material* innerBeMaterial = logicalInnerBeVolume->GetMaterial();
573 G4Tubs* innerBeTub = dynamic_cast<G4Tubs*>(logicalInnerBeVolume->GetSolid());
574 Z = 0.;
575 A = 0.;
576 for(i=0; i<innerBeMaterial->GetElementVector()->size(); i++){
577 Z += (innerBeMaterial->GetElement(i)->GetZ())*
578 (innerBeMaterial->GetFractionVector()[i]);
579 A += (innerBeMaterial->GetElement(i)->GetA())*
580 (innerBeMaterial->GetFractionVector()[i]);
581 }
582 Ionization = innerBeMaterial->GetIonisation()->GetMeanExcitationEnergy();
583 Density = innerBeMaterial->GetDensity()/(g/cm3);
584 Radlen = innerBeMaterial->GetRadlen();
585 cout<<"inner beryllium: Z: "<<Z<<" A: "<<(A/(g/mole))<<" Ionization: "<<(Ionization/eV)<<" Density: "<<Density<<" Radlen: "<<Radlen<<endl;
586 RkFitMaterial FitInnerBeMaterial(Z,A/g/mole,Ionization/eV,Density,Radlen/10.);
587 _BesBeamPipeMaterials.push_back(FitInnerBeMaterial);
588
589
590 /// gold
591 G4LogicalVolume* logicalGoldLayerVolume = const_cast<G4LogicalVolume*>(GDMLProcessor::GetInstance()->GetLogicalVolume("logicalgoldLayer"));
592 G4Material* goldLayerMaterial = logicalGoldLayerVolume->GetMaterial();
593 G4Tubs* goldLayerTub = dynamic_cast<G4Tubs*>(logicalGoldLayerVolume->GetSolid());
594
595 Z = 0.;
596 A = 0.;
597 for(i=0; i<goldLayerMaterial->GetElementVector()->size(); i++){
598 Z += (goldLayerMaterial->GetElement(i)->GetZ())*
599 (goldLayerMaterial->GetFractionVector()[i]);
600 A += (goldLayerMaterial->GetElement(i)->GetA())*
601 (goldLayerMaterial->GetFractionVector()[i]);
602 }
603 Ionization = goldLayerMaterial->GetIonisation()->GetMeanExcitationEnergy();
604 Density = goldLayerMaterial->GetDensity()/(g/cm3);
605 Radlen = goldLayerMaterial->GetRadlen();
606 cout<<"gold layer: Z: "<<Z<<" A: "<<(A/(g/mole))<<" Ionization: "<<(Ionization/eV)<<" Density: "<<Density<<" Radlen: "<<Radlen<<endl;
607 RkFitMaterial FitGoldLayerMaterial(Z,A/g/mole,Ionization/eV,Density,Radlen/10.);
608 _BesBeamPipeMaterials.push_back(FitGoldLayerMaterial);
609
610
611 /// now construct the cylinders
612 double radius, thick, length , z0;
613
614 /// innerwall of inner drift chamber
615 radius = innerwallTub->GetInnerRadius()/(cm);
616 thick = innerwallTub->GetOuterRadius()/(cm) - innerwallTub->GetInnerRadius()/(cm);
617 length = 2.0*innerwallTub->GetZHalfLength()/(cm);
618 z0 = 0.0;
619 cout<<"innerwall: "<<" radius: "<<radius<<" thick:"<<thick<<" length: "<<length<<endl;
620 RkFitCylinder innerwallCylinder(&_BesBeamPipeMaterials[1], radius, thick, length , z0);
621 _BesBeamPipeWalls.push_back(innerwallCylinder);
622
623 /// outer air, be attention the calculation of the radius and thick of the air cylinder is special
624 radius = outerBeTub->GetOuterRadius()/(cm);
625 thick = innerwallTub->GetInnerRadius()/(cm) - outerBeTub->GetOuterRadius()/(cm);
626 length = 2.0*innerwallTub->GetZHalfLength()/(cm);
627 z0 = 0.0;
628 cout<<"outer air: "<<" radius: "<<radius<<" thick:"<<thick<<" length: "<<length<<endl;
629 RkFitCylinder outerAirCylinder(&_BesBeamPipeMaterials[2], radius, thick, length , z0);
630 _BesBeamPipeWalls.push_back(outerAirCylinder);
631
632 /// outer Beryllium layer
633 radius = outerBeTub->GetInnerRadius()/(cm);
634 thick = outerBeTub->GetOuterRadius()/(cm) - outerBeTub->GetInnerRadius()/(cm);
635 length = 2.0*outerBeTub->GetZHalfLength()/(cm);
636 z0 = 0.0;
637 cout<<"outer Be: "<<" radius: "<<radius<<" thick:"<<thick<<" length: "<<length<<endl;
638 RkFitCylinder outerBeCylinder(&_BesBeamPipeMaterials[3], radius, thick, length , z0);
639 _BesBeamPipeWalls.push_back(outerBeCylinder);
640
641 /// oil layer
642 radius = oilLayerTub->GetInnerRadius()/(cm);
643 thick = oilLayerTub->GetOuterRadius()/(cm) - oilLayerTub->GetInnerRadius()/(cm);
644 length = 2.0*oilLayerTub->GetZHalfLength()/(cm);
645 z0 = 0.0;
646 cout<<"oil layer: "<<" radius: "<<radius<<" thick:"<<thick<<" length: "<<length<<endl;
647 RkFitCylinder oilLayerCylinder(&_BesBeamPipeMaterials[4], radius, thick, length , z0);
648 _BesBeamPipeWalls.push_back(oilLayerCylinder);
649
650 /// inner Beryllium layer
651 radius = innerBeTub->GetInnerRadius()/(cm);
652 thick = innerBeTub->GetOuterRadius()/(cm) - innerBeTub->GetInnerRadius()/(cm);
653 length = 2.0*innerBeTub->GetZHalfLength()/(cm);
654 z0 = 0.0;
655 cout<<"inner Be: "<<" radius: "<<radius<<" thick:"<<thick<<" length: "<<length<<endl;
656 RkFitCylinder innerBeCylinder(&_BesBeamPipeMaterials[5], radius, thick, length , z0);
657 _BesBeamPipeWalls.push_back(innerBeCylinder);
658
659 /// gold layer
660 radius = goldLayerTub->GetInnerRadius()/(cm);
661 thick = goldLayerTub->GetOuterRadius()/(cm) - goldLayerTub->GetInnerRadius()/(cm);
662 length = 2.0*goldLayerTub->GetZHalfLength()/(cm);
663 z0 = 0.0;
664 cout<<"gold layer: "<<" radius: "<<radius<<" thick:"<<thick<<" length: "<<length<<endl;
665 RkFitCylinder goldLayerCylinder(&_BesBeamPipeMaterials[6], radius, thick, length , z0);
666 _BesBeamPipeWalls.push_back(goldLayerCylinder);
667 delete aMdcG4Geo;
668 delete aBesG4Geo;
669}//end of setBesFromGdml
670
671void TRungeFitter::innerwall(TRunge& track, int l_mass,double y[6])const{
672 for(int i = 0; i < _BesBeamPipeWalls.size(); i++) {
673 _BesBeamPipeWalls[i].updateTrack(track,y);
674 }
675}
double sin(const BesAngle a)
Definition: BesAngle.h:210
double cos(const BesAngle a)
Definition: BesAngle.h:213
const double alpha
#define DBL_MAX
Definition: KalFitAlg.h:13
NTuple::Array< double > m_tof
Definition: MdcHistItem.h:105
const HepPoint3D ORIGIN
Constants.
Definition: TMDCUtil.cxx:47
#define WireHitLeft
Definition: TMDCWireHit.h:21
#define WireHitRight
Definition: TMDCWireHit.h:22
#define TFitAlreadyFitted
Definition: TMFitter.h:28
#define TFitFailed
Definition: TMFitter.h:30
#define TFitUnavailable
Definition: TMFitter.h:31
#define Runge
Definition: TRunge.h:19
TTree * t
Definition: binning.cxx:23
virtual double getReferField()=0
virtual double driftTimeToDist(double drifttime, int layid, int cellid, int lr, double entrance=0.0) const =0
virtual double getSigma(int layid, int lr, double dist, double entrance=0.0, double tanlam=0.0, double z=0.0, double Q=1000.0) const =0
virtual double getTprop(int lay, double z) const =0
virtual double getT0(int layid, int cellid) const =0
virtual double getTimeWalk(int layid, double Q) const =0
Cylinder is an Element whose shape is a cylinder.
Definition: RkFitCylinder.h:21
G4LogicalVolume * GetTopVolume()
Get the top(world) volume;.
struct MdcRec_wirhit * reccdc(void) const
returns a pointer to RECMDC_WIRHIT.
Definition: TMDCWireHit.h:224
float dDrift(unsigned) const
returns drift distance error.
Definition: TMDCWireHit.h:243
const TMDCWire *const wire(void) const
returns a pointer to a TMDCWire.
Definition: TMDCWireHit.h:218
unsigned id(void) const
returns id.
Definition: TMDCWire.h:207
unsigned localId(void) const
returns local id in a wire layer.
Definition: TMDCWire.h:213
unsigned layerId(void) const
returns layer id.
Definition: TMDCWire.h:219
A class to fit a TTrackBase object.
Definition: TMFitter.h:34
void setBesFromGdml(void)
TRungeFitter(const std::string &name)
Constructor.
void innerwall(TRunge &trk, int l_mass, double y[6]) const
std::vector< RkFitMaterial > _BesBeamPipeMaterials
Definition: TRungeFitter.h:58
virtual int fit(TTrackBase &) const
void tof(bool)
virtual ~TRungeFitter()
Destructor.
std::vector< RkFitCylinder > _BesBeamPipeWalls
Definition: TRungeFitter.h:57
void propagation(int)
void sag(bool)
A class to represent a track in tracking.
Definition: TRunge.h:65
A virtual class for a track class in tracking.
Definition: TTrackBase.h:46
virtual unsigned objectType(void) const
returns object type.
Definition: TTrackBase.h:268
double y[1000]
double x[1000]
const float pi
Definition: vector3.h:133