#include <MollerAnalysis.hh>
Public Member Functions | |
MollerAnalysis () | |
MollerAnalysis (MollerDetectorConstruction *) | |
MollerAnalysis (const MollerAnalysis &right) | |
const MollerAnalysis & | operator= (const MollerAnalysis &right) |
~MollerAnalysis () | |
void | BeginOfRunAction (const G4Run *) |
void | EndOfRunAction (const G4Run *) |
void | BeginOfEventAction (const G4Event *) |
void | EndOfEventAction (const G4Event *) |
void | UserSteppingAction (const G4Step *) |
void | SetRootFileName (const G4String &) |
void | SetMomentum0 (G4double E0, G4double px, G4double py, G4double pz) |
void | SetMomentum1 (G4double E1, G4double px, G4double py, G4double pz) |
void | SetMomentum2 (G4double E2, G4double px, G4double py, G4double pz) |
void | SetDiffXS (G4double diffXS_in) |
void | SetTotXS (G4double totXS_in) |
void | SetRate (G4double rate_in) |
void | SetProcess (G4int i) |
void | AddData (MollerDetectorHit *aHit) |
Private Attributes | |
MollerDetectorConstruction * | myDetector |
G4int | NUM_DETS |
MollerAnalysisMessenger * | analysisMessenger |
G4String | rootfileName |
G4int | ev_num |
G4int | index |
G4int | partType |
G4int | splitProcess |
G4int | creatorProcess |
TNtuple * | ntup |
Float_t | fntup [39] |
std::vector< G4int > | hitsCollID |
G4int | verboseLevel |
G4double | generator |
G4double | event |
G4double | kineE0 |
G4double | kineE1 |
G4double | kineE2 |
G4double | px0 |
G4double | py0 |
G4double | pz0 |
G4double | px1 |
G4double | py1 |
G4double | pz1 |
G4double | px2 |
G4double | py2 |
G4double | pz2 |
G4double | theta0 |
G4double | theta1 |
G4double | theta2 |
G4double | x0 |
G4double | y0 |
G4double | z0 |
G4double | rate |
G4double | totXS |
G4double | diffXS |
Static Private Attributes | |
static const G4bool | doDIPS = 0 |
static const G4bool | doCOLS = 0 |
static const G4bool | doQUADS = 0 |
Definition at line 41 of file MollerAnalysis.hh.
MollerAnalysis::MollerAnalysis | ( | ) |
Definition at line 42 of file MollerAnalysis.cc.
References analysisMessenger, creatorProcess, ev_num, index, kineE0, kineE1, kineE2, partType, px0, px1, px2, py0, py1, py2, pz0, pz1, pz2, splitProcess, theta0, theta1, theta2, verboseLevel, x0, y0, and z0.
00043 { 00044 if (gSystem) gSystem->ProcessEvents(); 00045 00046 ev_num=0; 00047 partType=-1; 00048 splitProcess=-1; 00049 creatorProcess=-1; 00050 index=-1; 00051 00052 // for (int i=0; i<NUM_DETS; i++) { 00053 // hitsCollID[i] = -1; 00054 // } 00055 00056 verboseLevel=0; 00057 00058 kineE0=0; kineE1=0; kineE2=0; 00059 00060 px0=0; py0=0; pz0=0; 00061 px1=0; py1=0; pz1=0; 00062 px2=0; py2=0; pz2=0; 00063 00064 theta0=0; theta1=0; theta2=0; 00065 00066 x0=0; y0=0; z0=0; 00067 00068 analysisMessenger = new MollerAnalysisMessenger(this); 00069 //G4cout << "Root file name is "<<rootfileName << G4endl; 00070 }
MollerAnalysis::MollerAnalysis | ( | MollerDetectorConstruction * | detector | ) |
Definition at line 71 of file MollerAnalysis.cc.
References analysisMessenger, creatorProcess, ev_num, MollerDetectorConstruction::getNumColl(), index, kineE0, kineE1, kineE2, myDetector, NUM_DETS, partType, px0, px1, px2, py0, py1, py2, pz0, pz1, pz2, splitProcess, theta0, theta1, theta2, verboseLevel, x0, y0, and z0.
00072 :myDetector(detector) 00073 { 00074 if (gSystem) gSystem->ProcessEvents(); 00075 00076 G4int NUM_COLS = myDetector->getNumColl(); 00077 00078 NUM_DETS = NUM_COLS; 00079 G4cout<<G4endl<<G4endl<<"NUM_COLS = "<<NUM_COLS<<G4endl<<G4endl; 00080 00081 ev_num=0; 00082 partType=-1; 00083 splitProcess=-1; 00084 creatorProcess=-1; 00085 index=-1; 00086 00087 // for (int i=0; i<NUM_DETS; i++) { 00088 // hitsCollID[i] = -1; 00089 // } 00090 00091 verboseLevel=0; 00092 00093 kineE0=0; kineE1=0; kineE2=0; 00094 00095 px0=0; py0=0; pz0=0; 00096 px1=0; py1=0; pz1=0; 00097 px2=0; py2=0; pz2=0; 00098 00099 theta0=0; theta1=0; theta2=0; 00100 00101 x0=0; y0=0; z0=0; 00102 00103 analysisMessenger = new MollerAnalysisMessenger(this); 00104 //G4cout << "Root file name is "<<rootfileName << G4endl; 00105 }
MollerAnalysis::MollerAnalysis | ( | const MollerAnalysis & | right | ) |
Definition at line 107 of file MollerAnalysis.cc.
References creatorProcess, ev_num, index, kineE0, kineE1, kineE2, partType, px0, px1, px2, py0, py1, py2, pz0, pz1, pz2, splitProcess, theta0, theta1, theta2, verboseLevel, x0, y0, and z0.
00108 : RootAnalysis() 00109 { 00110 if (gSystem) gSystem->ProcessEvents(); 00111 00112 ev_num=right.ev_num; 00113 partType=right.partType; 00114 splitProcess=right.splitProcess; 00115 creatorProcess=right.creatorProcess; 00116 index=right.index; 00117 00118 00119 // for (int i=0; i<NUM_DETS; i++) { 00120 // hitsCollID[i]=right.hitsCollID[i]; 00121 // } 00122 00123 verboseLevel=right.verboseLevel; 00124 00125 kineE0=right.kineE0; 00126 kineE1=right.kineE1; 00127 kineE2=right.kineE2; 00128 00129 px0=right.px0; 00130 py0=right.py0; 00131 pz0=right.pz0; 00132 00133 px1=right.px1; 00134 py1=right.py1; 00135 pz1=right.pz1; 00136 00137 px2=right.px2; 00138 py2=right.py2; 00139 pz2=right.pz2; 00140 00141 theta0=right.theta0; 00142 theta1=right.theta1; 00143 theta2=right.theta2; 00144 00145 x0=right.x0; 00146 y0=right.y0; 00147 z0=right.z0; 00148 00149 }
MollerAnalysis::~MollerAnalysis | ( | ) |
Definition at line 195 of file MollerAnalysis.cc.
References analysisMessenger.
00196 { 00197 delete analysisMessenger; 00198 }
void MollerAnalysis::AddData | ( | MollerDetectorHit * | aHit | ) | [inline, virtual] |
Reimplemented from RootAnalysis.
Definition at line 102 of file MollerAnalysis.hh.
References creatorProcess, diffXS, ev_num, event, fntup, MollerDetectorHit::GetIon(), MollerDetectorHit::GetKineticEnergy(), MollerDetectorHit::GetMomentum(), MollerDetectorHit::GetTrackID(), MollerDetectorHit::GetType(), MollerDetectorHit::GetVertexPos(), MollerDetectorHit::GetVolume(), MollerDetectorHit::GetWorldPos(), index, kineE0, kineE1, kineE2, ntup, px0, px1, px2, py0, py1, py2, pz0, pz1, pz2, rate, splitProcess, theta0, theta1, theta2, and totXS.
Referenced by EndOfEventAction().
00102 { 00103 00104 fntup[0] = aHit->GetIon(); 00105 fntup[1] = aHit->GetWorldPos().x(); 00106 fntup[2] = aHit->GetWorldPos().y(); 00107 fntup[3] = aHit->GetWorldPos().z(); 00108 //fntup[4] = (Float_t)x0; 00109 //fntup[5] = (Float_t)y0; 00110 //fntup[6] = (Float_t)z0; 00111 fntup[4] = aHit->GetVertexPos().x(); 00112 fntup[5] = aHit->GetVertexPos().y(); 00113 fntup[6] = aHit->GetVertexPos().z(); 00114 fntup[7] = aHit->GetKineticEnergy(); 00115 fntup[8] = aHit->GetMomentum().x(); 00116 fntup[9] = aHit->GetMomentum().y(); 00117 fntup[10] = aHit->GetMomentum().z(); 00118 00119 fntup[11] = (Float_t)kineE0; 00120 fntup[12] = (Float_t)px0; 00121 fntup[13] = (Float_t)py0; 00122 fntup[14] = (Float_t)pz0; 00123 fntup[15] = (Float_t)kineE1; 00124 fntup[16] = (Float_t)px1; 00125 fntup[17] = (Float_t)py1; 00126 fntup[18] = (Float_t)pz1; 00127 00128 fntup[19] = (Float_t)kineE2; 00129 fntup[20] = (Float_t)px2; 00130 fntup[21] = (Float_t)py2; 00131 fntup[22] = (Float_t)pz2; 00132 00133 fntup[23] = aHit->GetType(); 00134 fntup[24] = aHit->GetVolume(); 00135 fntup[25] = (Float_t)theta0; 00136 fntup[26] = (Float_t)theta1; 00137 fntup[27] = (Float_t)theta2; 00138 fntup[28] = (Float_t)ev_num; 00139 fntup[29] = (Float_t)splitProcess; 00140 fntup[30] = (Float_t)event; 00141 fntup[31] = (Float_t)creatorProcess; 00142 00143 if (!index&& 00144 (aHit->GetIon()||splitProcess==0)&& 00145 kineE0>0&& 00146 aHit->GetType()==0 00147 ) { 00148 00149 if (aHit->GetIon()) { 00150 fntup[32] = 2; 00151 fntup[33] = kineE2; 00152 fntup[34] = theta2; 00153 } 00154 else { 00155 fntup[32] = 1; 00156 fntup[33] = kineE1; 00157 fntup[34] = theta1; 00158 } 00159 } 00160 else if (aHit->GetTrackID()==1){ 00161 fntup[32] = 1; 00162 fntup[33] = kineE1; 00163 fntup[34] = theta1; 00164 } 00165 else if (aHit->GetTrackID()==2){ 00166 fntup[32] = 2; 00167 fntup[33] = kineE2; 00168 fntup[34] = theta2; 00169 } 00170 else { 00171 fntup[32] = 0; 00172 fntup[33] = 0; 00173 fntup[34] = 0; 00174 } 00175 00176 fntup[35] = aHit->GetTrackID(); 00177 00178 fntup[36] = diffXS; 00179 fntup[37] = totXS; 00180 fntup[38] = rate; 00181 00182 ntup->Fill(fntup); 00183 }
void MollerAnalysis::BeginOfEventAction | ( | const G4Event * | ) | [virtual] |
Reimplemented from RootAnalysis.
Definition at line 256 of file MollerAnalysis.cc.
References ev_num, gRootAnalysis, hitsCollID, NUM_DETS, and RootAnalysis::SetProcess().
00257 { 00258 ev_num++; 00259 if (gSystem) gSystem->ProcessEvents(); 00260 00261 G4SDManager * SDman = G4SDManager::GetSDMpointer(); 00262 00263 /* if(hitsCollID[0]<0 || hitsCollID[1]<0 || 00264 hitsCollID[2]<0 || hitsCollID[3]<0 || 00265 hitsCollID[4]<0 || hitsCollID[5]<0 || 00266 hitsCollID[6]<0 || hitsCollID[7]<0 || 00267 hitsCollID[8]<0 || hitsCollID[9]<0 || 00268 hitsCollID[10]<0 || hitsCollID[11]<0 || 00269 hitsCollID[12]<0 || hitsCollID[13]<0 || 00270 hitsCollID[14]<0 || hitsCollID[15]<0 || hitsCollID[16]<0 00271 ) 00272 */ 00273 00274 if((int) hitsCollID.size()==0) 00275 { 00276 G4cout<<G4endl<<G4endl<<"Setting hitsCollID"<<G4endl<<G4endl; 00277 00278 for (int i=0; i<NUM_DETS; i++) { 00279 G4String colNam = "detector"; 00280 colNam+=G4UIcommand::ConvertToString(i+1); 00281 colNam+="/hitsColl"; 00282 00283 G4cout << colNam << G4endl; 00284 00285 hitsCollID.push_back (SDman->GetCollectionID(colNam)); 00286 } 00287 00288 G4cout<<G4endl<<"Finished setting hitsCollID"<<G4endl<<G4endl; 00289 00290 } 00291 00292 00293 gRootAnalysis->SetProcess(1); 00294 /* These are you lines Dustin! Only these three lines. */ 00295 //gRootAnalysis->SetMomentum0(); 00296 //gRootAnalysis->SetMomentum1(); 00297 //gRootAnalysis->SetMomentum2(); 00298 /* **************** */ 00299 }
void MollerAnalysis::BeginOfRunAction | ( | const G4Run * | ) | [virtual] |
Reimplemented from RootAnalysis.
Definition at line 205 of file MollerAnalysis.cc.
References RootAnalysis::hfile, ntup, and rootfileName.
00206 { 00207 if (gSystem) gSystem->ProcessEvents(); 00208 00209 G4cout << "Run Action has been entered\n" << G4endl; 00210 00211 G4cout << "Root file name is "<<rootfileName << G4endl; 00212 00213 00214 // This works: 00215 //G4String filenamechar=("~/scratch/ROOTfiles/"+rootfileName+".root"); 00216 G4String filenamechar=(rootfileName+".root"); 00217 00218 // This works: 00219 /* 00220 G4String filenamechar; 00221 filenamechar="~/scratch/ROOTfiles/"; 00222 filenamechar += rootfileName; 00223 filenamechar += ".root"; 00224 G4cout<<filenamechar<<G4endl; 00225 */ 00226 00227 //char filenamechar[200]; 00228 /* 00229 time_t rawtime; 00230 tm * ptm; 00231 time ( &rawtime ); 00232 ptm = gmtime ( &rawtime ); 00233 snprintf(filenamechar,200,"~/scratch/ROOTfiles/coll_test0_%02i%02i%02i_%02i%02i%02i.root", 00234 ptm->tm_year-100,ptm->tm_mon+1,ptm->tm_mday,(ptm->tm_hour+19)%24,ptm->tm_min,ptm->tm_sec); 00235 // snprintf(filenamechar,200,"~/scratch/ROOTfiles/moller_%02i%02i%02i_%02i%02i%02i.root", 00236 // ptm->tm_year-100,ptm->tm_mon+1,ptm->tm_mday,(ptm->tm_hour+19)%24,ptm->tm_min,ptm->tm_sec); 00237 printf(filenamechar); 00238 */ 00239 //snprintf(filenamechar,200,"~/scratch/ROOTfiles/%s.root",rootfileName); 00240 00241 hfile = new TFile(filenamechar,"RECREATE","Root/G4 analysis"); 00242 //hfile = new TFile(rootfileName,"RECREATE","Root/G4 analysis"); 00243 00244 00245 string myntuple_string; 00246 myntuple_string="ion:x:y:z:x0:y0:z0:kineE:px:py:pz:kineE0:px0:py0:pz0:kineE1:px1:py1:pz1:kineE2:px2:py2:pz2:type:volume:theta0:theta1:theta2:ev_num:process:event:creator:hit:kineE_org:theta_org:track:diffXS:totXS:rate"; 00247 00248 char *string_ntup = new char[myntuple_string.length()+1]; 00249 strncpy (string_ntup,myntuple_string.c_str(),myntuple_string.length()+1); 00250 00251 ntup = new TNtuple("geant","G4 Results",(char *)string_ntup); 00252 00253 G4cout << "Run Action has been completed\n" << G4endl; 00254 }
void MollerAnalysis::EndOfEventAction | ( | const G4Event * | anEvent | ) | [virtual] |
Reimplemented from RootAnalysis.
Definition at line 483 of file MollerAnalysis.cc.
References AddData(), MollerDetectorHit::GetTrackID(), MollerDetectorHit::GetType(), gRootAnalysis, hitsCollID, NUM_DETS, RootAnalysis::SetMomentum0(), RootAnalysis::SetMomentum1(), RootAnalysis::SetMomentum2(), and RootAnalysis::SetProcess().
00484 { 00485 00486 if (gSystem) gSystem->ProcessEvents(); 00487 00488 G4HCofThisEvent * HCE = anEvent->GetHCofThisEvent(); 00489 00490 event = anEvent->GetEventID(); 00491 00492 std::vector<MollerDetectorHitsCollection*> THC; 00493 00494 00495 if(HCE) { 00496 for (int i=0; i<NUM_DETS; i++) { 00497 THC.push_back ((MollerDetectorHitsCollection*)(HCE->GetHC(hitsCollID[i]))); 00498 } 00499 } 00500 G4int checkTrack1=0; 00501 G4int checkTrack2=0; 00502 // G4cout<<"Checking for hits..."<<G4endl; 00503 for (int i=0; i<NUM_DETS; i++) { 00504 // somehow use: aHit->GetTrackID()==1 and aHit->GetTrackID()==2 with an indicator of first entry for that event/track? 00505 if(THC[i]){ 00506 int n_hit = THC[i]->entries(); 00507 /* // Just do first hit; rate ~65GHz - missing second electron? 00508 if (n_hit>0) { 00509 MollerDetectorHit* aHit = (*THC[i])[0]; 00510 //Decide here what kind of data we want to keep 00511 if (aHit->GetType()==0) AddData(aHit); 00512 if (aHit->GetType()==4) AddData(aHit); 00513 } 00514 } 00515 }*/ 00516 // G4cout<<"Done checking for hits..."<<G4endl; 00517 checkTrack1=0; //need to reset for each detector... 00518 checkTrack2=0; 00519 if (n_hit>0) { 00520 00521 for(int i1=0;i1<n_hit;i1++) { 00522 MollerDetectorHit* aHit = (*THC[i])[i1]; 00523 //aHit->Print(); 00524 //aHit->Draw(); 00525 // G4cout<<"Ev "<<event<<", det: "<<i<<", hit: "<<i1<<", track: "<<aHit->GetTrackID()<<", check: "<<checkTrack1<<", "<<checkTrack2<<", "<<G4endl; 00526 //Try to get only the first hit, but for both primary electrons... 00527 //Really need to define an ntuple variable for the n_hit somehow... 00528 //and/or do this checking in the analysis scripts... 00529 00530 if ((aHit->GetTrackID()==1)&&(checkTrack1==0)) { 00531 if (aHit->GetType()==0) AddData(aHit); 00532 if (aHit->GetType()==4) AddData(aHit); 00533 checkTrack1++; 00534 } else if ((aHit->GetTrackID()==2)&&(checkTrack2==0)){ 00535 if (aHit->GetType()==0) AddData(aHit); 00536 if (aHit->GetType()==4) AddData(aHit); 00537 checkTrack2++; 00538 } 00539 00540 } 00541 00542 } 00543 } 00544 } 00545 00546 gRootAnalysis->SetMomentum0(0,0,0,0); 00547 gRootAnalysis->SetMomentum1(0,0,0,0); 00548 gRootAnalysis->SetMomentum2(0,0,0,0); 00549 gRootAnalysis->SetProcess(0); 00550 00551 }
void MollerAnalysis::EndOfRunAction | ( | const G4Run * | ) | [virtual] |
Reimplemented from RootAnalysis.
Definition at line 555 of file MollerAnalysis.cc.
References RootAnalysis::hfile.
const MollerAnalysis & MollerAnalysis::operator= | ( | const MollerAnalysis & | right | ) |
Definition at line 151 of file MollerAnalysis.cc.
References creatorProcess, ev_num, index, kineE0, kineE1, kineE2, partType, px0, px1, px2, py0, py1, py2, pz0, pz1, pz2, splitProcess, theta0, theta1, theta2, verboseLevel, x0, y0, and z0.
00152 { 00153 if (gSystem) gSystem->ProcessEvents(); 00154 00155 ev_num=right.ev_num; 00156 partType=right.partType; 00157 splitProcess=right.splitProcess; 00158 creatorProcess=right.creatorProcess; 00159 index=right.index; 00160 00161 // for (int i=0; i<NUM_DETS; i++) { 00162 // hitsCollID[i]=right.hitsCollID[i]; 00163 // } 00164 00165 verboseLevel=right.verboseLevel; 00166 00167 kineE0=right.kineE0; 00168 kineE1=right.kineE1; 00169 kineE2=right.kineE2; 00170 00171 px0=right.px0; 00172 py0=right.py0; 00173 pz0=right.pz0; 00174 00175 px1=right.px1; 00176 py1=right.py1; 00177 pz1=right.pz1; 00178 00179 px2=right.px2; 00180 py2=right.py2; 00181 pz2=right.pz2; 00182 00183 theta0=right.theta0; 00184 theta1=right.theta1; 00185 theta2=right.theta2; 00186 00187 x0=right.x0; 00188 y0=right.y0; 00189 z0=right.z0; 00190 00191 return *this; 00192 00193 }
void MollerAnalysis::SetDiffXS | ( | G4double | diffXS_in | ) | [inline, virtual] |
Reimplemented from RootAnalysis.
Definition at line 86 of file MollerAnalysis.hh.
References diffXS.
00086 { 00087 diffXS = diffXS_in; 00088 }
void MollerAnalysis::SetMomentum0 | ( | G4double | E0, | |
G4double | px, | |||
G4double | py, | |||
G4double | pz | |||
) | [inline, virtual] |
void MollerAnalysis::SetMomentum1 | ( | G4double | E1, | |
G4double | px, | |||
G4double | py, | |||
G4double | pz | |||
) | [inline, virtual] |
void MollerAnalysis::SetMomentum2 | ( | G4double | E2, | |
G4double | px, | |||
G4double | py, | |||
G4double | pz | |||
) | [inline, virtual] |
void MollerAnalysis::SetProcess | ( | G4int | i | ) | [inline, virtual] |
Reimplemented from RootAnalysis.
Definition at line 98 of file MollerAnalysis.hh.
References index.
00098 { 00099 index=i; 00100 }
void MollerAnalysis::SetRate | ( | G4double | rate_in | ) | [inline, virtual] |
Reimplemented from RootAnalysis.
Definition at line 94 of file MollerAnalysis.hh.
References rate.
00094 { 00095 rate = rate_in; 00096 }
void MollerAnalysis::SetRootFileName | ( | const G4String & | nam | ) |
Definition at line 199 of file MollerAnalysis.cc.
References rootfileName.
Referenced by MollerAnalysisMessenger::SetNewValue().
00200 { 00201 rootfileName = nam; 00202 }
void MollerAnalysis::SetTotXS | ( | G4double | totXS_in | ) | [inline, virtual] |
Reimplemented from RootAnalysis.
Definition at line 90 of file MollerAnalysis.hh.
References totXS.
00090 { 00091 totXS = totXS_in; 00092 }
void MollerAnalysis::UserSteppingAction | ( | const G4Step * | aStep | ) | [virtual] |
Get data that will go into ROOT tree
Reimplemented from RootAnalysis.
Definition at line 302 of file MollerAnalysis.cc.
References creatorProcess, kineE0, kineE1, kineE2, partType, px0, px1, px2, py0, py1, py2, pz0, pz1, pz2, splitProcess, theta0, theta1, theta2, x0, y0, and z0.
00303 { 00304 /** Get data that will go into ROOT tree **/ 00305 00306 G4Track* fTrack = aStep->GetTrack(); 00307 00308 00309 if (fTrack->GetTrackStatus()!=fAlive) { 00310 // printf("track is not 'alive.'\n"); 00311 return; 00312 } 00313 00314 G4int tracking = 1; 00315 00316 G4StepPoint* thePrePoint = aStep->GetPreStepPoint(); 00317 G4StepPoint* thePostPoint = aStep->GetPostStepPoint(); 00318 G4VPhysicalVolume* thePrePV = thePrePoint->GetPhysicalVolume(); 00319 G4Material* material = fTrack->GetMaterial(); 00320 00321 const G4VProcess* originProcess = fTrack->GetCreatorProcess(); 00322 const G4VProcess* stepProcess = thePostPoint->GetProcessDefinedStep(); 00323 00324 G4String process; 00325 if (originProcess!=0) process = originProcess->GetProcessName(); 00326 G4String process1; 00327 if (stepProcess!=0) process1 = stepProcess->GetProcessName(); 00328 00329 00330 G4String volume = fTrack->GetVolume()->GetName(); 00331 G4String particle = fTrack->GetDefinition()->GetParticleName(); 00332 00333 00334 if (process=="eIoni") creatorProcess=0; 00335 else if (process=="conv") creatorProcess=1; 00336 else if (process=="Decay") creatorProcess=2; 00337 else if (process=="annihil") creatorProcess=3; 00338 else if (process=="eBrem") creatorProcess=4; 00339 else if (process=="phot") creatorProcess=5; 00340 else if (process=="compt") creatorProcess=6; 00341 else if (process=="Transportation") creatorProcess=7; 00342 else if (process=="msc") creatorProcess=8; 00343 else creatorProcess=9; 00344 00345 //G4bool ion = process=="eIoni"; 00346 00347 if ((particle == "e-") 00348 && (fTrack->GetVertexKineticEnergy()==11*GeV) 00349 && ((thePrePoint->GetKineticEnergy()- 00350 thePostPoint->GetKineticEnergy())>1.00*GeV) 00351 && (process1=="eIoni")) 00352 { 00353 if (kineE0==0 || kineE1==0) { 00354 kineE0 = thePrePoint->GetKineticEnergy(); 00355 kineE1 = thePostPoint->GetKineticEnergy(); 00356 px0 = thePrePoint->GetMomentumDirection().x(); 00357 py0 = thePrePoint->GetMomentumDirection().y(); 00358 pz0 = thePrePoint->GetMomentumDirection().z(); 00359 px1 = thePostPoint->GetMomentumDirection().x(); 00360 py1 = thePostPoint->GetMomentumDirection().y(); 00361 pz1 = thePostPoint->GetMomentumDirection().z(); 00362 theta0 = sqrt(px0*px0 + py0*py0); 00363 theta1 = sqrt(px1*px1 + py1*py1); 00364 x0 = thePostPoint->GetPosition().x(); 00365 y0 = thePostPoint->GetPosition().y(); 00366 z0 = thePostPoint->GetPosition().z(); 00367 } 00368 00369 if (process1=="eIoni") splitProcess=0; 00370 else if (process1=="conv") splitProcess=1; 00371 else if (process1=="Decay") splitProcess=2; 00372 else if (process1=="annihil") splitProcess=3; 00373 else if (process1=="eBrem") splitProcess=4; 00374 else if (process1=="photo") splitProcess=5; 00375 else if (process1=="compton") splitProcess=6; 00376 else if (process1=="Transportation") splitProcess=7; 00377 else splitProcess=8; 00378 } 00379 00380 if (tracking) { 00381 00382 G4ThreeVector worldPos = fTrack->GetPosition(); 00383 G4ThreeVector worldMomentum = fTrack->GetMomentumDirection(); 00384 G4ThreeVector vertexMomentum = fTrack->GetVertexMomentumDirection(); 00385 00386 G4double scat_ang = sqrt(vertexMomentum.x()*vertexMomentum.x() + 00387 vertexMomentum.y()*vertexMomentum.y()); 00388 00389 /* if (scat_ang<0.002 || scat_ang>0.020) { 00390 return; 00391 }*/ 00392 00393 G4double kineE = fTrack->GetKineticEnergy(); 00394 G4double kineE0 = fTrack->GetVertexKineticEnergy(); 00395 00396 G4bool ion = (process=="eIoni"); 00397 00398 if (particle=="e-") partType=0; 00399 else if (particle=="e+") partType=1; 00400 else if (particle=="proton") partType=2; 00401 else if (particle=="antiproton") partType=3; 00402 else if (particle=="gamma") partType=4; 00403 else partType=5; 00404 /* 00405 if (partType!=4) { 00406 00407 fntup[0] = (Float_t)ion; 00408 fntup[1] = (Float_t)worldPos.x(); 00409 fntup[2] = (Float_t)worldPos.y(); 00410 fntup[3] = (Float_t)worldPos.z(); 00411 fntup[4] = 0; 00412 fntup[5] = 0; 00413 fntup[6] = 0; 00414 fntup[7] = (Float_t)kineE; 00415 fntup[8] = (Float_t)worldMomentum.x(); 00416 fntup[9] = (Float_t)worldMomentum.y(); 00417 fntup[10] = (Float_t)worldMomentum.z(); 00418 00419 fntup[11] = (Float_t)kineE0; 00420 fntup[12] = (Float_t)px0; 00421 fntup[13] = (Float_t)py0; 00422 fntup[14] = (Float_t)pz0; 00423 fntup[15] = (Float_t)kineE1; 00424 fntup[16] = (Float_t)px1; 00425 fntup[17] = (Float_t)py1; 00426 fntup[18] = (Float_t)pz1; 00427 00428 fntup[19] = (Float_t)fTrack->GetVertexKineticEnergy(); 00429 fntup[20] = (Float_t)fTrack->GetVertexMomentumDirection().x(); 00430 fntup[21] = (Float_t)fTrack->GetVertexMomentumDirection().y(); 00431 fntup[22] = (Float_t)fTrack->GetVertexMomentumDirection().z(); 00432 00433 fntup[23] = (Float_t)partType; 00434 fntup[24] = (Float_t)atof(volume); 00435 fntup[25] = (Float_t)theta0; 00436 fntup[26] = (Float_t)theta1; 00437 fntup[27] = (Float_t)theta2; 00438 fntup[28] = (Float_t)ev_num; 00439 fntup[29] = (Float_t)ion; 00440 fntup[30] = (Float_t)ev_num; 00441 00442 00443 ntup->Fill(fntup); 00444 } 00445 */ 00446 } 00447 00448 /* KILL particle if it hits a detector */ 00449 /* made of tungsten! */ 00450 /* 00451 // if(thePrePV->GetName()=="1"||thePrePV->GetName()=="2"|| 00452 thePrePV->GetName()=="3"||thePrePV->GetName()=="4"|| 00453 thePrePV->GetName()=="5"||thePrePV->GetName()=="6"|| 00454 thePrePV->GetName()=="7"||thePrePV->GetName()=="8"|| 00455 thePrePV->GetName()=="9"||thePrePV->GetName()=="10"|| 00456 thePrePV->GetName()=="11"||thePrePV->GetName()=="12"|| 00457 thePrePV->GetName()=="13"||thePrePV->GetName()=="14"|| 00458 thePrePV->GetName()=="15"||thePrePV->GetName()=="16"||thePrePV->GetName()=="17" 00459 ){ 00460 */ 00461 // G4cout<<"Material = "<<material->GetName(); 00462 if((material->GetName()=="Tungsten")||(material->GetName()=="Pb")||(material->GetName()=="Copper")){ 00463 if (kineE2 == 0 && process=="eIoni") { 00464 kineE2 = fTrack->GetVertexKineticEnergy(); 00465 px2 = fTrack->GetVertexMomentumDirection().x(); 00466 py2 = fTrack->GetVertexMomentumDirection().y(); 00467 pz2 = fTrack->GetVertexMomentumDirection().z(); 00468 theta2 = sqrt(px2*px2 + py2*py2); 00469 x0 = fTrack->GetVertexPosition().x(); 00470 y0 = fTrack->GetVertexPosition().y(); 00471 z0 = fTrack->GetVertexPosition().z(); 00472 00473 } 00474 // G4cout<<" Stopping the track..."<<G4endl; 00475 // fTrack->SetTrackStatus(fStopAndKill); 00476 fTrack->SetTrackStatus(fKillTrackAndSecondaries); 00477 return; 00478 } 00479 // G4cout<<endl; 00480 }
Definition at line 197 of file MollerAnalysis.hh.
Referenced by MollerAnalysis(), and ~MollerAnalysis().
G4int MollerAnalysis::creatorProcess [private] |
Definition at line 201 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), and UserSteppingAction().
G4double MollerAnalysis::diffXS [private] |
Definition at line 236 of file MollerAnalysis.hh.
Referenced by AddData(), and SetDiffXS().
const G4bool MollerAnalysis::doCOLS = 0 [static, private] |
Definition at line 188 of file MollerAnalysis.hh.
const G4bool MollerAnalysis::doDIPS = 0 [static, private] |
Definition at line 187 of file MollerAnalysis.hh.
const G4bool MollerAnalysis::doQUADS = 0 [static, private] |
Definition at line 189 of file MollerAnalysis.hh.
G4int MollerAnalysis::ev_num [private] |
Definition at line 200 of file MollerAnalysis.hh.
Referenced by AddData(), BeginOfEventAction(), MollerAnalysis(), and operator=().
G4double MollerAnalysis::event [private] |
Definition at line 210 of file MollerAnalysis.hh.
Referenced by AddData().
Float_t MollerAnalysis::fntup[39] [private] |
Definition at line 204 of file MollerAnalysis.hh.
Referenced by AddData().
G4double MollerAnalysis::generator [private] |
Definition at line 210 of file MollerAnalysis.hh.
std::vector<G4int> MollerAnalysis::hitsCollID [private] |
Definition at line 206 of file MollerAnalysis.hh.
Referenced by BeginOfEventAction(), and EndOfEventAction().
G4int MollerAnalysis::index [private] |
Definition at line 200 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), and SetProcess().
G4double MollerAnalysis::kineE0 [private] |
Definition at line 212 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), SetMomentum0(), and UserSteppingAction().
G4double MollerAnalysis::kineE1 [private] |
Definition at line 213 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), SetMomentum1(), and UserSteppingAction().
G4double MollerAnalysis::kineE2 [private] |
Definition at line 214 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), SetMomentum2(), and UserSteppingAction().
Definition at line 193 of file MollerAnalysis.hh.
Referenced by MollerAnalysis().
TNtuple* MollerAnalysis::ntup [private] |
Definition at line 203 of file MollerAnalysis.hh.
Referenced by AddData(), and BeginOfRunAction().
G4int MollerAnalysis::NUM_DETS [private] |
Definition at line 195 of file MollerAnalysis.hh.
Referenced by BeginOfEventAction(), EndOfEventAction(), and MollerAnalysis().
G4int MollerAnalysis::partType [private] |
Definition at line 201 of file MollerAnalysis.hh.
Referenced by MollerAnalysis(), operator=(), and UserSteppingAction().
G4double MollerAnalysis::px0 [private] |
Definition at line 216 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), SetMomentum0(), and UserSteppingAction().
G4double MollerAnalysis::px1 [private] |
Definition at line 219 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), SetMomentum1(), and UserSteppingAction().
G4double MollerAnalysis::px2 [private] |
Definition at line 222 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), SetMomentum2(), and UserSteppingAction().
G4double MollerAnalysis::py0 [private] |
Definition at line 217 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), SetMomentum0(), and UserSteppingAction().
G4double MollerAnalysis::py1 [private] |
Definition at line 220 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), SetMomentum1(), and UserSteppingAction().
G4double MollerAnalysis::py2 [private] |
Definition at line 223 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), SetMomentum2(), and UserSteppingAction().
G4double MollerAnalysis::pz0 [private] |
Definition at line 218 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), SetMomentum0(), and UserSteppingAction().
G4double MollerAnalysis::pz1 [private] |
Definition at line 221 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), SetMomentum1(), and UserSteppingAction().
G4double MollerAnalysis::pz2 [private] |
Definition at line 224 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), SetMomentum2(), and UserSteppingAction().
G4double MollerAnalysis::rate [private] |
Definition at line 234 of file MollerAnalysis.hh.
G4String MollerAnalysis::rootfileName [private] |
Definition at line 198 of file MollerAnalysis.hh.
Referenced by BeginOfRunAction(), and SetRootFileName().
G4int MollerAnalysis::splitProcess [private] |
Definition at line 201 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), and UserSteppingAction().
G4double MollerAnalysis::theta0 [private] |
Definition at line 226 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), SetMomentum0(), and UserSteppingAction().
G4double MollerAnalysis::theta1 [private] |
Definition at line 227 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), SetMomentum1(), and UserSteppingAction().
G4double MollerAnalysis::theta2 [private] |
Definition at line 228 of file MollerAnalysis.hh.
Referenced by AddData(), MollerAnalysis(), operator=(), SetMomentum2(), and UserSteppingAction().
G4double MollerAnalysis::totXS [private] |
Definition at line 235 of file MollerAnalysis.hh.
Referenced by AddData(), and SetTotXS().
G4int MollerAnalysis::verboseLevel [private] |
Definition at line 208 of file MollerAnalysis.hh.
Referenced by MollerAnalysis(), and operator=().
G4double MollerAnalysis::x0 [private] |
Definition at line 230 of file MollerAnalysis.hh.
Referenced by MollerAnalysis(), operator=(), and UserSteppingAction().
G4double MollerAnalysis::y0 [private] |
Definition at line 231 of file MollerAnalysis.hh.
Referenced by MollerAnalysis(), operator=(), and UserSteppingAction().
G4double MollerAnalysis::z0 [private] |
Definition at line 232 of file MollerAnalysis.hh.
Referenced by MollerAnalysis(), operator=(), and UserSteppingAction().