| Classes | Job Modules | Data Objects | Services | Algorithms | Tools | Packages | Directories | Tracs |

In This Package:

B8DecayGen.cc

Go to the documentation of this file.
00001 // This is a HepEVT generator of B8 decay (resulting positrons & alphas). 
00002  // Modified from K-40 generator by Raymond Tsang, h0237236@hkusua.hku.hk, 30 Aug, 2006.
00003 
00004 //\author Fengpeng An
00005 //anfp@ihep.ac.cn, 8 Nov., 2009
00006 #include <cstdlib>
00007 #include <cstdio>
00008 #include <cstring>
00009 #include <iostream>
00010 #include <math.h>
00011 
00012 #include <CLHEP/Vector/ThreeVector.h>
00013 #include <CLHEP/Random/Randomize.h>
00014 #include <CLHEP/Units/PhysicalConstants.h>
00015 
00016 #include "TFile.h"
00017 #include "TH1F.h"
00018 #include "TCanvas.h"
00019 #include "TRandom.h"
00020 #include "TMath.h"
00021 
00022 using namespace std;
00023 using namespace CLHEP;
00024 
00025 void ProcessArgs(int argc, char** argv, long* rseed, char** outfilename, 
00026                  unsigned int* nevents );
00027 void Usage();
00028 
00029 
00030 int main(int argc, char** argv) {
00031   long rseed = 0;
00032   char* outFilename = NULL;
00033   unsigned int nEvents = 1000000000; // a billion. Default to something big for piping 
00034   ProcessArgs(argc, argv, &rseed, &outFilename, &nEvents);
00035   Hep3Vector p0; // the beta  momentum vector. Unit GeV/c
00036   Hep3Vector p1; // the alpha 1 momentum vector. Unit GeV/c
00037   
00038   FILE* stream = stdout;
00039   if( outFilename!=NULL ) {
00040     stream = fopen(outFilename, "w");
00041     if( stream==NULL ) {
00042       printf("Please enter a valid filename.\n");
00043       Usage();
00044       exit(0);
00045     }
00046   }
00047   HepRandom::setTheSeed(rseed);
00048    gRandom->SetSeed(rseed);
00049   //start by printing some information to comment lines
00050   fprintf(stream, "# File generated by %s.\n", argv[0]);
00051   fprintf(stream, "# Random seed for generator = %ld.\n", rseed);
00052   
00053 
00054 // The root file used by this generator should be located at 
00055 // $SITEROOT/dybgaudi/Generators/RadioActivity/B8/
00056     std::string specfilename = "B8BetaSpectrum.root";
00057     const char* prefix = getenv("SITEROOT");
00058     if(prefix) {
00059         if(strlen(prefix)>0) {
00060             std::string newname = prefix;
00061             newname += "/dybgaudi/Generators/RadioActivity/B8/";
00062             newname += specfilename;
00063             specfilename = newname;
00064         }
00065     }
00066 
00067 
00068 //cout<<specfilename<<endl;
00069 
00070   TFile  *betaSpec = TFile::Open(specfilename.c_str());
00071 
00072     if(!betaSpec)
00073         cout << "Can not open B8BetaSpectrum.root. "
00074                << "Make sure SITEROOT or B8BetaSpectrum.root is set to it's location."
00075                << "rootfile location is at $SITEROOT/dybgaudi/Generators/RadioActivity/B8/" << endl;
00076     if(!(betaSpec->IsOpen()))
00077         cout << "Can not open B8BetaSpectrum.root "
00078                << "Make sure SITEROOT is set to it's location."
00079                << "rootfile location is at $SITEROOT/dybgaudi/Generators/RadioActivity/B8/" << endl;
00080 
00081 
00082   double cosTheta, theta, azimuth; // angles used
00083   unsigned int i;
00084 
00085   TH1F *h0 = (TH1F*)betaSpec->Get("B8");
00086   
00087   for( i=0 ; i<nEvents ; i++ ) {
00088 
00089       cosTheta = RandFlat::shoot(-1, 1);
00090       azimuth = RandFlat::shoot( 2.0*M_PI );
00091       theta = acos(cosTheta);
00092       p0.setRThetaPhi(1, acos(cosTheta), azimuth);
00093 
00094       cosTheta = RandFlat::shoot(-1, 1);
00095       azimuth = RandFlat::shoot( 2.0*M_PI );
00096       theta = acos(cosTheta);
00097       p1.setRThetaPhi(1, acos(cosTheta), azimuth);
00098 
00099       const double alpha_mass_c2=3727.0*MeV;
00100       double ke = h0->GetRandom();
00101       double ka =1.513*MeV/keV;
00102       double pe = sqrt(ke*ke + 2*ke*electron_mass_c2/keV);
00103       double pa = sqrt(ka*ka + 2*ka*alpha_mass_c2/keV);
00104       p0=p0*pe*keV;
00105       p1=p1*pa*keV;
00106       fprintf(stream, "3\n");
00107       fprintf(stream, "1\t-11 0 0 %e %e %e %e\n", p0.x()/GeV, p0.y()/GeV, p0.z()/GeV, electron_mass_c2/GeV );
00108       fprintf(stream, "1\t1000020040 0 0 %e %e %e %e\n", p1.x()/GeV, p1.y()/GeV, p1.z()/GeV,alpha_mass_c2/GeV );
00109       fprintf(stream, "1\t1000020040 0 0 %e %e %e %e\n", -p1.x()/GeV, -p1.y()/GeV, -p1.z()/GeV,alpha_mass_c2/GeV );
00110 
00111   }
00112 
00113 
00114   return 0;
00115 }
00116 
00117 
00118 
00119 
00120 
00121 void ProcessArgs(int argc, char** argv, long *rseed, 
00122                                  char** outfilename, unsigned int* nevents) {
00123         int i;
00124         for( i=1 ; i<argc ; i++ ) {
00125                 if( strcmp(argv[i], "-seed")==0 ) {
00126                         i++;
00127                         sscanf(argv[i], "%ld", rseed);
00128                                                   }
00129                  else if( strcmp(argv[i], "-o")==0 ) {
00130                         i++;
00131                         *outfilename = new char[strlen(argv[i]) +1];
00132                         strcpy(*outfilename, argv[i]);
00133                                                      } 
00134                 else if( strcmp(argv[i], "-n")==0 ) { 
00135                         i++;
00136                         sscanf(argv[i], "%ud", nevents);
00137                                                     }  
00138                 else {
00139                         Usage();
00140                         exit(0);
00141                      }
00142                                    }
00143 }
00144 
00145 void Usage() {
00146         printf("B8.exe [-seed seed] [-o outputfilename] [-n nevents]\n");
00147 
00148 }
| Classes | Job Modules | Data Objects | Services | Algorithms | Tools | Packages | Directories | Tracs |

Generated on Mon Apr 11 20:06:49 2011 for B8 by doxygen 1.4.7