BOSS 7.0.7
BESIII Offline Software System
Loading...
Searching...
No Matches
BesMdcSD.cc
Go to the documentation of this file.
1//#include "DedxPar.hh"
2#include "BesMdcSD.hh"
3#include "G4HCofThisEvent.hh"
4#include "G4Step.hh"
5#include "G4ThreeVector.hh"
6#include "G4SDManager.hh"
7#include "G4UnitsTable.hh"
8#include "G4ios.hh"
9#include "G4RunManager.hh"
10#include "ReadBoostRoot.hh"
11#include "G4Svc/IG4Svc.h"
12#include "G4Svc/G4Svc.h"
16#include "GaudiKernel/DataSvc.h"
17#include "TFile.h"
18#include "TH1F.h"
19#include "TH2D.h"
20
21#include "GaudiKernel/Bootstrap.h"
22#include "GaudiKernel/IService.h"
23#include "GaudiKernel/Service.h"
24#include "GaudiKernel/SmartDataPtr.h"
25
26#include <iostream>
27
28
29BesMdcSD::BesMdcSD(G4String name)
31{
32 collectionName.insert("BesMdcHitsCollection");
33 collectionName.insert("BesMdcTruthCollection");
34
35 mdcGeoPointer=BesMdcGeoParameter::GetGeo();
36 mdcCalPointer=new BesMdcCalTransfer;
37
38 IMdcGeomSvc* ISvc;
39 StatusCode sc=Gaudi::svcLocator()->service("MdcGeomSvc", ISvc);
40 if (!sc.isSuccess())
41 std::cout<<"BesMdcSD::Could not open MdcGeomSvc"<<std::endl;
42 mdcGeomSvc=dynamic_cast<MdcGeomSvc *>(ISvc);
43
44 IG4Svc* tmpSvc;
45 sc=Gaudi::svcLocator()->service("G4Svc", tmpSvc);
46 if (!sc.isSuccess())
47 G4cout <<" MdcSD::Error,could not open G4Svc"<<G4endl;
48 m_G4Svc=dynamic_cast<G4Svc *>(tmpSvc);
49
50 if(m_G4Svc->GetMdcDedxFlag()==1){
51 G4cout <<" MdcSD: Use sampled dedx instead of Geant4 value"<<G4endl;
53 }
54
55 ////dedx sim
56
57 //get DedxSimData
58 std::string dedxTDSPath = "/Calib/DedxSim";
59 IDataProviderSvc* pCalibDataSvc;
60 sc = Gaudi::svcLocator()->service("CalibDataSvc", pCalibDataSvc, true);
61 if (!sc.isSuccess())
62 std::cout << "BesMdcSD::Could not open CalibDataSvc" << std::endl;
63 m_calibDataSvc = dynamic_cast<CalibDataSvc*>(pCalibDataSvc);
64 if (!sc.isSuccess())
65 {
66 std::cout << "Could not get CalibDataSvc"
67 << std::endl;
68 }
69 SmartDataPtr<CalibData::DedxSimData> pDedxSimData(m_calibDataSvc, dedxTDSPath);
70 m_version = pDedxSimData->getVersion();
71 m_numDedxHists = pDedxSimData->gethistNo();
72 m_numBg = pDedxSimData->getRangeNo();
73 if (m_version == 0) m_numTheta = 10;
74 else m_numTheta = pDedxSimData->getThetaNo();
75 m_dedx_hists = new TH1F[m_numDedxHists];
76 for (G4int i = 0; i < m_numBg; i++)
77 {
78 m_bgRange.push_back(pDedxSimData->getRange(i));
79 }
80 for (G4int i = 0; i < m_numDedxHists; i++)
81 {
82 m_dedx_hists[i] = pDedxSimData->getHist(i);
83 }
84
85 //get CalibCurSvc
86 IDedxCurSvc* tmp_dedxCurSvc;
87 sc = Gaudi::svcLocator()->service("DedxCurSvc", tmp_dedxCurSvc, true);
88 if (!sc.isSuccess())
89 {
90 std::cout << "Could not get DedxCurSvc"
91 << std::endl;
92 }
93 m_pDedxCurSvc = dynamic_cast<DedxCurSvc*>(tmp_dedxCurSvc);
94
95 if(m_G4Svc->MdcRootFlag())
96 {
97 m_tupleMdc = m_G4Svc->GetTupleMdc();
98 sc = m_tupleMdc->addItem("betaGamma",m_betaGamma);
99 sc = m_tupleMdc->addItem("fitval",m_fitval);
100 sc = m_tupleMdc->addItem("dedx",m_dedx);
101 sc = m_tupleMdc->addItem("de",m_de);
102 //sc = m_tupleMdc->addItem("length",m_length);
103 //sc = m_tupleMdc->addItem("random",m_random);
104 sc = m_tupleMdc->addItem("charge", m_charge);
105 sc = m_tupleMdc->addItem("costheta", m_costheta);
106 }
107}
108
110 delete []m_dedx_hists;
111}
112
113void BesMdcSD::Initialize(G4HCofThisEvent* HCE)
114{
115 hitsCollection = new BesMdcHitsCollection
116 (SensitiveDetectorName,collectionName[0]);
117 static G4int HCID = -1;
118 if(HCID<0)
119 { HCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]); }
120 HCE->AddHitsCollection( HCID, hitsCollection );
121 G4int i,j;
122 for(i=0; i<43;i++){
123 for(j=0;j<288;j++){
124 hitPointer[i][j]=-1;
125 truthPointer[i][j]=-1;
126 }
127 }
128}
129
130//for MC Truth
131void BesMdcSD::BeginOfTruthEvent(const G4Event* evt)
132{
133 truthCollection = new BesMdcHitsCollection
134 (SensitiveDetectorName,collectionName[1]);
135 // G4cout<<" begin event "<<evt->GetEventID()<<G4endl;
136}
137
138void BesMdcSD::EndOfTruthEvent(const G4Event* evt)
139{ static G4int HLID=-1;
140 if(HLID<0)
141 {
142 HLID = G4SDManager::GetSDMpointer()->
143 GetCollectionID(collectionName[1]);
144 }
145 G4HCofThisEvent* HCE = evt->GetHCofThisEvent();
146 HCE->AddHitsCollection(HLID,truthCollection);
147}
148
149G4bool BesMdcSD::ProcessHits(G4Step* aStep,G4TouchableHistory*)
150{
151 G4Track* gTrack = aStep->GetTrack() ;
152
153 G4double globalT=gTrack->GetGlobalTime();//Time since the event in which the track belongs is created
154 if(isnan(globalT)){
155 G4cout<<"MdcSD:error, globalT is nan "<<G4endl;
156 return false;
157 }
158 if(globalT > 2000)return false; //MDC T window is 2 microsecond
159
160 //skip neutral tracks
161 G4int charge = gTrack->GetDefinition()->GetPDGCharge();
162 if (charge == 0) return false;
163
164 //skip no energy deposit tracks
165 G4double stepLength=aStep->GetStepLength();
166 if(stepLength==0){
167 // G4cout<<"step length equal 0!!"<<G4endl;
168 return false;
169 }
170
171 G4double edep = aStep->GetTotalEnergyDeposit() / stepLength;
172 if(edep==0.) return false;
173
174 // get position of the track at the beginning and at the end of step
175 G4StepPoint* prePoint = aStep->GetPreStepPoint() ;
176 G4StepPoint* postPoint = aStep->GetPostStepPoint() ;
177
178 //get position coordinate
179 G4ThreeVector pointIn = prePoint->GetPosition();
180 G4ThreeVector pointOut = postPoint->GetPosition();
181
182 // get physical volumes
183 const G4VTouchable *touchable = prePoint->GetTouchable();
184 G4VPhysicalVolume *volume = touchable->GetVolume(0);
185
186 G4double driftD = 0;
187 G4double driftT = 0;
188 G4double edepTemp = 0;
189 G4double lengthTemp = 0;
190 G4int cellId=0;
191 G4int layerId = touchable->GetVolume(1)->GetCopyNo();
192 if(volume->IsReplicated()){
193 cellId = touchable->GetReplicaNumber();
194 }else{
195 cellId=touchable->GetVolume(0)->GetCopyNo();
196 }
197 if(layerId==18&&(cellId==27||cellId==42))return false; // no sense wire
198
199 if(ReadBoostRoot::GetMdc() == 2) { //gdml
200 //layerId 0-35 -> CopyNo 0-35 in gdml
201 //layerId 36-42 -> CopyNo (36,37),(38,39),...(48,49)
202 if(layerId >= 36) layerId = 36 + (layerId-36)/2;
203 }
204
205 G4double halfLayerLength=mdcGeomSvc->Layer(layerId)->Length()/2.;
206 if(((fabs(pointIn.z())-halfLayerLength)>-7.)
207 &&((fabs(pointOut.z())-halfLayerLength)>-7.))return false;//Out sensitive area
208
209 G4int trackID= gTrack->GetTrackID(); //G4 track ID of current track.
210 G4int truthID, g4TrackID;
211 GetCurrentTrackIndex(truthID, g4TrackID); //ID of current primary track.
212
213 G4double theta=gTrack->GetMomentumDirection().theta();
214
215 G4ThreeVector hitPosition=0;
216 G4double transferT=0;
217 driftD = Distance(layerId,cellId,pointIn,pointOut,hitPosition,transferT);
218
219 G4double posPhi, wirePhi;
220 posPhi=hitPosition.phi();//from -pi to pi
221 if(posPhi<0)posPhi += 2*pi;
222 wirePhi=mdcGeoPointer->SignalWire(layerId,cellId).Phi(hitPosition.z());//from 0 to 2pi
223
224 G4int posFlag;
225 if(posPhi<=wirePhi){
226 posFlag = 0;
227 }else{
228 posFlag = 1;
229 }
230 // if x axis is between pos and wire, phi will has a jump of one of them.
231 if(posPhi < 1. && wirePhi > 5.)posFlag = 1;
232 if(posPhi > 5. && wirePhi < 1.)posFlag = 0;
233
234 G4ThreeVector hitLine=pointOut-pointIn;
235 G4double enterAngle=hitLine.phi()-wirePhi;
236 while(enterAngle<-pi/2.)enterAngle+=pi;
237 while(enterAngle>pi/2.)enterAngle-=pi;
238
239 if(m_G4Svc->GetMdcDedxFlag()==1){
240 G4double betaGamma=aStep->GetPreStepPoint()->GetBeta() * aStep->GetPreStepPoint()->GetGamma();
241 if(betaGamma<0.01)return false;//too low momentum
242 //if (betaGamma < 10.0) betaGamma = 10.0;
243
244 G4double eCount=dedxSample(betaGamma, charge, theta);
245 edep=eCount;
246 }
247
248 BesMdcHit* newHit = new BesMdcHit();
249 newHit->SetTrackID(truthID);
250 //newHit->SetTrkID(trackID);
251 newHit->SetLayerNo(layerId);
252 newHit->SetCellNo(cellId);
253 newHit->SetEdep(edep);
254 newHit->SetPos(hitPosition);
255 newHit->SetDriftD(driftD);
256 newHit->SetTheta(theta);
257 newHit->SetPosFlag(posFlag);
258 newHit->SetEnterAngle(enterAngle);
259
260 //Transfer hit pointer to BesMdcCalTransfer
261 mdcCalPointer->SetHitPointer(newHit);
262
263 driftT=mdcCalPointer->D2T(driftD);
264 globalT+=transferT;
265 driftT+=globalT;
266
267 newHit->SetDriftT (driftT);
268 newHit->SetGlobalT(globalT);
269
270 if (hitPointer[layerId][cellId] == -1) {
271 hitsCollection->insert(newHit);
272 G4int NbHits = hitsCollection->entries();
273 hitPointer[layerId][cellId]=NbHits-1;
274 }
275 else
276 {
277 G4int pointer=hitPointer[layerId][cellId];
278 if (g4TrackID == trackID) {
279 G4double preDriftT=(*hitsCollection)[pointer]->GetDriftT();
280 }
281
282 G4double preDriftT = (*hitsCollection)[pointer]->GetDriftT();
283 if (driftT < preDriftT) {
284 (*hitsCollection)[pointer]->SetTrackID(truthID);
285 //(*hitsCollection)[pointer]->SetTrkID(trackID);
286 (*hitsCollection)[pointer]->SetDriftD(driftD);
287 (*hitsCollection)[pointer]->SetDriftT(driftT);
288 (*hitsCollection)[pointer]->SetPos(hitPosition);
289 (*hitsCollection)[pointer]->SetGlobalT(globalT);
290 (*hitsCollection)[pointer]->SetTheta(theta);
291 (*hitsCollection)[pointer]->SetPosFlag(posFlag);
292 (*hitsCollection)[pointer]->SetEnterAngle(enterAngle);
293 }
294
295 delete newHit;
296 }
297
298 //for mc truth
299 if(truthCollection){
300 if(g4TrackID==trackID){ //This track is the primary track & will appear in MC truth
301 G4int pointer=truthPointer[layerId][cellId];
302 if(pointer==-1){
303 BesMdcHit* truthHit = new BesMdcHit();
304 truthHit->SetTrackID (truthID);
305 truthHit->SetLayerNo(layerId);
306 truthHit->SetCellNo(cellId);
307 truthHit->SetEdep (edep);
308 truthHit->SetPos (hitPosition);
309 truthHit->SetDriftD (driftD);
310 truthHit->SetDriftT (driftT);
311 truthHit->SetGlobalT(globalT);
312 truthHit->SetTheta(theta);
313 truthHit->SetPosFlag(posFlag);
314 truthHit->SetEnterAngle(enterAngle);
315
316 truthCollection->insert(truthHit);
317 G4int NbHits = truthCollection->entries();
318 truthPointer[layerId][cellId]=NbHits-1;
319 }
320 else {
321 if(truthID == (*truthCollection)[pointer]->GetTrackID()){
322 G4double preDriftT=(*truthCollection)[pointer]->GetDriftT();
323 if(driftT<preDriftT){
324 (*truthCollection)[pointer]->SetDriftD(driftD);
325 (*truthCollection)[pointer]->SetDriftT(driftT);
326 (*truthCollection)[pointer]->SetPos(hitPosition);
327 (*truthCollection)[pointer]->SetGlobalT(globalT);
328 (*truthCollection)[pointer]->SetTheta(theta);
329 (*truthCollection)[pointer]->SetPosFlag(posFlag);
330 (*truthCollection)[pointer]->SetEnterAngle(enterAngle);
331 }
332 edepTemp=(*truthCollection)[pointer]->GetEdep();
333 (*truthCollection)[pointer]->SetEdep(edepTemp+edep);
334 } else {
335 BesMdcHit* truthHit = new BesMdcHit();
336 truthHit->SetTrackID (truthID);
337 truthHit->SetLayerNo(layerId);
338 truthHit->SetCellNo(cellId);
339 truthHit->SetEdep(edep);
340 truthHit->SetPos(hitPosition);
341 truthHit->SetDriftD (driftD);
342 truthHit->SetDriftT (driftT);
343 truthHit->SetGlobalT(globalT);
344 truthHit->SetTheta(theta);
345 truthHit->SetPosFlag(posFlag);
346 truthHit->SetEnterAngle(enterAngle);
347
348 truthCollection->insert(truthHit);
349 G4int NbHits = truthCollection->entries();
350 truthPointer[layerId][cellId]=NbHits-1;
351 }
352 }
353 }
354 }
355
356 //newHit->Print();
357// newHit->Draw();
358
359 return true;
360}
361
362void BesMdcSD::EndOfEvent(G4HCofThisEvent*)
363{
364 if (verboseLevel>0) {
365 hitsCollection->PrintAllHits();
366 /*
367 G4int NbHits = hitsCollection->entries();
368 G4cout << "\n-------->Hits Collection: in this event they are " << NbHits
369 << " hits in the MDC chambers: " << G4endl;
370 for (G4int i=0;i<NbHits;i++) (*hitsCollection)[i]->Print();
371 */
372 }
373}
374
375G4double BesMdcSD::Distance(G4int layerId, G4int cellId, G4ThreeVector pointIn, G4ThreeVector pointOut,G4ThreeVector& hitPosition,G4double& transferT)
376{
377 //For two lines r=r1+t1.v1 & r=r2+t2.v2
378 //the closest approach is d=|(r2-r1).(v1 x v2)|/|v1 x v2|
379 //the point where closest approach are
380 //t1=(v1 x v2).[(r2-r1) x v2]/[(v1 x v2).(v1 x v2)]
381 //t2=(v1 x v2).[(r2-r1) x v1]/[(v1 x v2).(v1 x v2)]
382 //if v1 x v2=0 means two lines are parallel
383 //d=|(r2-r1) x v1|/|v1|
384
385 G4double t1,distance,dInOut,dHitIn,dHitOut;
386 //Get wirepoint @ endplate
387 G4ThreeVector east=mdcGeomSvc->Wire(layerId,cellId)->Backward();
388 G4ThreeVector west=mdcGeomSvc->Wire(layerId,cellId)->Forward();
389 G4ThreeVector wireLine=east-west;
390 G4ThreeVector hitLine=pointOut-pointIn;
391
392 G4ThreeVector hitXwire=hitLine.cross(wireLine);
393 G4ThreeVector wire2hit=east-pointOut;
394 //Hitposition is the position on hit line where closest approach
395 //of two lines, but it may out the area from pointIn to pointOut
396 if(hitXwire.mag()==0){
397 distance=wireLine.cross(wire2hit).mag()/wireLine.mag();
398 hitPosition=pointIn;
399 }else{
400 t1=hitXwire.dot(wire2hit.cross(wireLine))/hitXwire.mag2();
401 hitPosition=pointOut+t1*hitLine;
402
403 dInOut=(pointOut-pointIn).mag();
404 dHitIn=(hitPosition-pointIn).mag();
405 dHitOut=(hitPosition-pointOut).mag();
406 if(dHitIn<=dInOut && dHitOut<=dInOut){ //Between point in & out
407 distance=fabs(wire2hit.dot(hitXwire)/hitXwire.mag());
408 }else if(dHitOut>dHitIn){ // out pointIn
409 distance=wireLine.cross(pointIn-east).mag()/wireLine.mag();
410 hitPosition=pointIn;
411 }else{ // out pointOut
412 distance=wireLine.cross(pointOut-east).mag()/wireLine.mag();
413 hitPosition=pointOut;
414 }
415 }
416
417 //Calculate signal transferT on wire
418 G4double halfLayerLength=mdcGeomSvc->Layer(layerId)->Length()/2.;
419 G4double halfWireLength=wireLine.mag()/2.;
420 G4double transferZ=0;
421 if(layerId%2==0){
422 transferZ=halfLayerLength+hitPosition.z(); //West readout
423 }else{
424 transferZ=halfLayerLength-hitPosition.z(); //East readout
425 }
426 if(layerId<8){
427 transferT=transferZ*halfWireLength/halfLayerLength/220;
428 }else{
429 transferT=transferZ*halfWireLength/halfLayerLength/240;
430 }
431
432 return distance;
433
434}
435
437{
438 dEdE_mylanfunc = new TF1("dEdE_mylanfunc",
439 "[3]*TMath::Exp([2]*((x[0]-[0])/[1]+TMath::Exp(-1*((x[0]-[0])/[1]))))",
440 0,
441 7500);
442 //dEdE_mylanfunc->SetParameters(2009.35,559.776,-1.0932,6327.38);
443 dEdE_mylanfunc->SetParNames("MPV","Sigma","constant1","constant2");
444}
445
446G4double BesMdcSD::dedxSample(G4double betagamma, G4double charge, G4double theta)
447{
448 G4double x = betagamma;
449 G4double fitval = GetValDedxCurve(x, charge);
450 if(fitval <= 0)return 0;
451
452 G4double random1, random2, dedx1, dedx2, de;
453 G4double standard1, standard2, beta_temp1, beta_temp2;
454 G4double dedx = -1;
455
456 G4int range_idx, bg_idx, angle_idx, charge_idx, hist_idx;
457 range_idx = GetBetagammaIndex(betagamma);
458 angle_idx = GetAngleIndex(theta);
459 charge_idx = GetChargeIndex(charge);
460
461 if (range_idx == -1)
462 {
463 while (dedx <= 0)
464 {
465 bg_idx = 0;
466 hist_idx = bg_idx * 2 * m_numTheta + angle_idx * 2 + charge_idx;
467 random1 = m_dedx_hists[hist_idx].GetRandom();
468 beta_temp1 = (m_bgRange[0] + m_bgRange[1]) / 2;
469 standard1 = GetValDedxCurve(beta_temp1, charge);
470 dedx = random1 + fitval - standard1;
471 }
472 }
473 else if (range_idx == m_numBg - 1)
474 {
475 while (dedx <= 0)
476 {
477 bg_idx = (G4int)(range_idx / 2);
478 hist_idx = bg_idx * 2 * m_numTheta + angle_idx * 2 + charge_idx;
479 random1 = m_dedx_hists[hist_idx].GetRandom();
480 beta_temp1 = (m_bgRange[m_numBg - 2] +
481 m_bgRange[m_numBg - 1]) / 2;
482 standard1 = GetValDedxCurve(beta_temp1, charge);
483 dedx = random1 + fitval - standard1;
484 }
485 }
486 else
487 {
488 //case 1: Given betagamma fall in one histograph range
489 if (range_idx % 2 == 0)
490 {
491 while (dedx <= 0)
492 {
493 bg_idx = (G4int)(range_idx / 2);
494 hist_idx = bg_idx * 2 * m_numTheta + angle_idx * 2 + charge_idx;
495 random1 = m_dedx_hists[hist_idx].GetRandom();
496 beta_temp1 = (m_bgRange[range_idx] +
497 m_bgRange[range_idx + 1]) / 2;
498 standard1 = GetValDedxCurve(beta_temp1, charge);
499 dedx1 = random1 + fitval - standard1;
500 dedx = dedx1;
501 }
502 }
503 //case 2: Given betagamma fall in interval between
504 // two histographs.
505 else
506 {
507 while (dedx <= 0)
508 {
509 //standard1
510 beta_temp1 = (m_bgRange[range_idx - 1] +
511 m_bgRange[range_idx]) / 2;
512 standard1 = GetValDedxCurve(beta_temp1, charge);
513
514 //stardard2
515 beta_temp2 = (m_bgRange[range_idx + 1] +
516 m_bgRange[range_idx + 2]) / 2;
517 standard2 = GetValDedxCurve(beta_temp2, charge);
518
519 //random1
520 bg_idx = (G4int)(range_idx / 2);
521 hist_idx = bg_idx * 2 * m_numTheta + angle_idx * 2 + charge_idx;
522 random1 = m_dedx_hists[hist_idx].GetRandom();
523
524 //random2
525 bg_idx++;
526 hist_idx = bg_idx * 2 * m_numTheta + angle_idx * 2 + charge_idx;
527 random2 = m_dedx_hists[hist_idx].GetRandom();
528
529 //combine dedx1 and dedx2
530 dedx1 = random1 + fitval - standard1;
531 dedx2 = random2 + fitval - standard2;
532 dedx = (dedx2 * (x - m_bgRange[range_idx]) +
533 dedx1 * (m_bgRange[range_idx + 1] - x)) /
534 (m_bgRange[range_idx + 1] - m_bgRange[range_idx]);
535 }
536 }
537 }
538
539
540 de= dedx;// * y/10./1.5;// stepLength unit is mm in Geant4, cm in Rec.& Cal. software
541 // dedx counts *1.5 in dedx Cal.
542
543 if(m_G4Svc->MdcRootFlag())
544 {
545 m_betaGamma= x;
546 m_fitval= fitval;
547 //m_trackId = trackId;
548 //m_layer = layerId;
549 //m_wire = cellId;
550 //m_random= random;
551 m_dedx= dedx;
552 m_de= de;
553 //m_length=y;
554 m_charge = charge;
555 m_costheta = cos(theta);
556 m_tupleMdc->write();
557 }
558 return de;
559}
560
561/*-----------------------------------------------------
562Func: GetValDedxCurve
563Pre: A betagamma is given.
564Post: Return dE/dx value from betagamma~dE/dx curve.
565-----------------------------------------------------*/
566G4double BesMdcSD::GetValDedxCurve(G4double x, G4double charge)
567{
568 G4int dedxflag = -1;
569 G4int size = -1;
570 G4double A = 0.;
571 G4double B = 0.;
572 G4double C = 0.;
573 std::vector<G4double> par;
574 G4double val;
575 G4int index = -1;
576
577 par.clear();
578 dedxflag = m_pDedxCurSvc->getCurve(0);
579 size = m_pDedxCurSvc->getCurveSize();
580 for (G4int i = 1; i < size; i++) {
581 par.push_back(m_pDedxCurSvc->getCurve(i));
582 }
583
584 if (x < 4.5)
585 A = 1;
586 else if(x < 10)
587 B = 1;
588 else
589 C = 1;
590
591 G4double partA = par[0] * pow(sqrt(x * x + 1), par[2]) / pow(x, par[2]) *
592 (par[1] - par[5] * log(pow(1 / x, par[3]))) -
593 par[4] + exp(par[6] + par[7] * x);
594 G4double partB = par[8] * pow(x, 3) + par[9] * pow(x, 2) + par[10] * x + par[11];
595 G4double partC = - par[12] * log(par[15] + pow(1 / x, par[13])) + par[14];
596
597 val = 550 * (A * partA + B * partB + C * partC);
598 return val;
599
600}
601
602/*-----------------------------------------------------
603Func: GetBetagammaIndex
604Pre : A betagamma of a track is given.
605Post: Return index of the betagamma range.
606-----------------------------------------------------*/
607G4int BesMdcSD::GetBetagammaIndex(G4double bg)
608{
609 if (bg < m_bgRange[0]) return -1;
610 for (G4int i = 0; i < m_numBg - 1; i++)
611 {
612 if (bg > m_bgRange[i] && bg < m_bgRange[i + 1])
613 {
614 return i;
615 }
616 }
617 if (bg > m_bgRange[m_numBg - 1])
618 return (m_numBg - 1);
619 return -1;
620}
621
622/*-----------------------------------------------------
623Func: GetAngleIndex
624Pre : A theta of a track is given.
625Post: version 0 - Return index of the angle (|cos(theta)| in 10 bins),
626 version 1 - Return index of the angle ( cos(theta) in m_numTheta bins).
627-----------------------------------------------------*/
628G4int BesMdcSD::GetAngleIndex(G4double theta)
629{
630 if (m_version == 0) {
631 if (fabs(cos(theta)) >= 0.93) return 9;
632 return (G4int)(fabs(cos(theta)) * 10 / 0.93);
633 }
634 else {
635 G4double cos_min = -0.93;
636 G4double cos_max = 0.93;
637 G4int nbin = m_numTheta;
638 G4double cos_step = (cos_max - cos_min) / nbin;
639
640 if (cos(theta) < cos_min) return 0;
641 else if (cos_min < cos(theta) && cos(theta) < cos_max) {
642 return (G4int)((cos(theta) - cos_min) / cos_step);
643 }
644 else return nbin - 1;
645 }
646}
647
648/*-----------------------------------------------------
649Func: GetChargeIndex
650Pre : A charge of a track is given.
651Post: Return index of charge (pos->0 ~ neg->1).
652-----------------------------------------------------*/
653G4int BesMdcSD::GetChargeIndex(G4int charge)
654{
655 if (charge > 0) return 0;
656 if (charge == 0) return -1; // warning: -1 is illegal, for charged tracks are expected.
657 if (charge < 0) return 1;
658 return -1;
659}
660
double cos(const BesAngle a)
Definition: BesAngle.h:213
G4THitsCollection< BesMdcHit > BesMdcHitsCollection
Definition: BesMdcHit.hh:78
EvtComplex exp(const EvtComplex &c)
Definition: EvtComplex.hh:252
***************************************************************************************Pseudo Class RRes *****************************************************************************************Parameters and physical constants **Maarten sept ************************************************************************DOUBLE PRECISION xsmu **************************************************************************PARTICLE DATA all others are from PDG *Only resonances with known widths into electron pairs are sept ************************************************************************C Declarations C
Definition: RRes.h:29
double D2T(double driftDNew)
void SetHitPointer(BesMdcHit *hit)
static BesMdcGeoParameter * GetGeo(void)
BesMdcWire SignalWire(int, int)
void SetEdep(G4double de)
Definition: BesMdcHit.hh:41
void SetDriftT(G4double time)
Definition: BesMdcHit.hh:44
void SetEnterAngle(G4double angle)
Definition: BesMdcHit.hh:47
void SetCellNo(G4int cell)
Definition: BesMdcHit.hh:40
void SetPos(G4ThreeVector xyz)
Definition: BesMdcHit.hh:42
void SetTrackID(G4int track)
Definition: BesMdcHit.hh:38
void SetLayerNo(G4int layer)
Definition: BesMdcHit.hh:39
void SetTheta(G4double angle)
Definition: BesMdcHit.hh:46
void SetDriftD(G4double distance)
Definition: BesMdcHit.hh:43
void SetGlobalT(G4double time)
Definition: BesMdcHit.hh:45
void SetPosFlag(G4int flag)
Definition: BesMdcHit.hh:48
void dedxFuncInti(void)
Definition: BesMdcSD.cc:436
void BeginOfTruthEvent(const G4Event *)
Definition: BesMdcSD.cc:131
void EndOfTruthEvent(const G4Event *)
Definition: BesMdcSD.cc:138
G4double Distance(G4int, G4int, G4ThreeVector, G4ThreeVector, G4ThreeVector &, G4double &)
Definition: BesMdcSD.cc:375
G4bool ProcessHits(G4Step *, G4TouchableHistory *)
Definition: BesMdcSD.cc:149
void EndOfEvent(G4HCofThisEvent *)
Definition: BesMdcSD.cc:362
BesMdcSD(G4String)
Definition: BesMdcSD.cc:29
void Initialize(G4HCofThisEvent *)
Definition: BesMdcSD.cc:113
~BesMdcSD()
Definition: BesMdcSD.cc:109
double Phi(void) const
void GetCurrentTrackIndex(G4int &trackIndex, G4int &g4TrackId) const
Definition: G4Svc.h:32
bool MdcRootFlag()
Definition: G4Svc.h:117
NTuple::Tuple * GetTupleMdc()
Definition: G4Svc.h:96
int GetMdcDedxFlag()
Definition: G4Svc.h:94
virtual const int getCurveSize()=0
virtual const double getCurve(int i)=0
Definition: IG4Svc.h:30
double Length(void) const
Definition: MdcGeoLayer.h:161
HepPoint3D Forward(void) const
Definition: MdcGeoWire.h:129
HepPoint3D Backward(void) const
Definition: MdcGeoWire.h:128
const MdcGeoWire *const Wire(unsigned id)
Definition: MdcGeomSvc.cxx:768
const MdcGeoLayer *const Layer(unsigned id)
Definition: MdcGeomSvc.cxx:784
static G4int GetMdc()
double x[1000]
float charge
float bg
const float pi
Definition: vector3.h:133