![]() |
![]() |
#include <SimTransport/HectorProducer/interface/Hector.h>
Public Member Functions | |
void | add (const HepMC::GenEvent *ev, const edm::EventSetup &es) |
HepMC::GenEvent * | addPartToHepMC (HepMC::GenEvent *event) |
void | clear () |
void | clearApertureFlags () |
void | filterD1 () |
void | filterFP420 () |
void | filterZDC () |
int | getDirect (unsigned int part_n) const |
Hector (const edm::ParameterSet &ps, bool verbosity, bool FP420Transport, bool ZDCTransport) | |
void | print () const |
virtual | ~Hector () |
Private Attributes | |
string | beam1filename |
string | beam2filename |
double | etacut |
double | lengthd1 |
double | lengthfp420 |
double | lengthzdc |
H_BeamLine * | m_beamlineD11 |
H_BeamLine * | m_beamlineD12 |
H_BeamLine * | m_beamlineFP4201 |
H_BeamLine * | m_beamlineFP4202 |
H_BeamLine * | m_beamlineZDC1 |
H_BeamLine * | m_beamlineZDC2 |
std::map< unsigned int, H_BeamParticle * > | m_beamPart |
std::map< unsigned int, int > | m_direct |
std::map< unsigned int, double > | m_eAtTrPoint |
std::map< unsigned int, double > | m_eta |
bool | m_FP420Transport |
std::map< unsigned int, bool > | m_isCharged |
std::map< unsigned int, bool > | m_isStoppedd1 |
std::map< unsigned int, bool > | m_isStoppedfp420 |
std::map< unsigned int, bool > | m_isStoppedzdc |
std::map< unsigned int, int > | m_pdg |
std::map< unsigned int, double > | m_pz |
H_RecRPObject * | m_rp420_b |
H_RecRPObject * | m_rp420_f |
float | m_rpp420_b |
float | m_rpp420_f |
float | m_rppd1 |
float | m_rppzdc |
double | m_sig_e |
double | m_sigmaSTX |
double | m_sigmaSTY |
bool | m_smearAng |
bool | m_smearE |
std::map< unsigned int, double > | m_TxAtTrPoint |
std::map< unsigned int, double > | m_TyAtTrPoint |
bool | m_verbosity |
std::map< unsigned int, double > | m_xAtTrPoint |
std::map< unsigned int, double > | m_yAtTrPoint |
bool | m_ZDCTransport |
edm::ESHandle< ParticleDataTable > | pdt |
Definition at line 34 of file Hector.h.
Hector::Hector | ( | const edm::ParameterSet & | ps, | |
bool | verbosity, | |||
bool | FP420Transport, | |||
bool | ZDCTransport | |||
) |
Definition at line 17 of file Hector.cc.
References b1, b2, beam1filename, beam2filename, GenMuonPlsPt100GeV_cfg::cout, e, lat::endl(), etacut, edm::ParameterSet::getParameter(), lengthd1, lengthfp420, lengthzdc, m_beamlineD11, m_beamlineD12, m_beamlineFP4201, m_beamlineFP4202, m_beamlineZDC1, m_beamlineZDC2, m_FP420Transport, m_rpp420_b, m_rpp420_f, m_rppd1, m_rppzdc, m_sig_e, m_sigmaSTX, m_sigmaSTY, m_smearAng, m_smearE, m_verbosity, m_ZDCTransport, alivecheck_mergeAndRegister::msg, and cms::Exception::what().
00017 : 00018 m_verbosity(verbosity), 00019 m_FP420Transport(FP420Transport), 00020 m_ZDCTransport(ZDCTransport) 00021 { 00022 00023 edm::LogInfo ("Hector") << "===================================================================\n" 00024 << "=== Hector Constructor start =====\n"; 00025 00026 // Create LHC beam line 00027 edm::ParameterSet hector_par = param.getParameter<edm::ParameterSet>("Hector"); 00028 00029 // User definitons 00030 lengthfp420 = hector_par.getParameter<double>("BeamLineLengthFP420" ); 00031 m_rpp420_f = (float) hector_par.getParameter<double>("RP420f"); 00032 m_rpp420_b = (float) hector_par.getParameter<double>("RP420b"); 00033 lengthzdc = hector_par.getParameter<double>("BeamLineLengthZDC" ); 00034 lengthd1 = hector_par.getParameter<double>("BeamLineLengthD1" ); 00035 beam1filename = hector_par.getParameter<string>("Beam1"); 00036 beam2filename = hector_par.getParameter<string>("Beam2"); 00037 m_rppzdc = (float) lengthzdc ; 00038 m_rppd1 = (float) lengthd1 ; 00039 m_smearAng = hector_par.getParameter<bool>("smearAng"); 00040 m_sigmaSTX = hector_par.getParameter<double>("sigmaSTX" ); 00041 m_sigmaSTY = hector_par.getParameter<double>("sigmaSTY" ); 00042 m_smearE = hector_par.getParameter<bool>("smearEnergy"); 00043 m_sig_e = hector_par.getParameter<double>("sigmaEnergy"); 00044 etacut = hector_par.getParameter<double>("EtaCutForHector" ); 00045 00046 edm::LogInfo ("Hector") << "Hector parameters: \n" 00047 << " lengthfp420: " << lengthfp420 << "\n" 00048 << " m_rpp420_f: " << m_rpp420_f << "\n" 00049 << " m_rpp420_b: " << m_rpp420_b << "\n" 00050 << " lengthzdc: " << lengthzdc << "\n" 00051 << " lengthd1: " << lengthd1 << "\n"; 00052 00053 if(m_verbosity) { 00054 cout << "===================================================================" << endl; 00055 cout << "=== Hector: m_FP420Transport=" << m_FP420Transport <<"m_ZDCTransport " << m_ZDCTransport << endl; 00056 cout << "=== lengthfp420: " << lengthfp420 << endl; 00057 cout << "=== m_rpp420_f: " << m_rpp420_f << endl; 00058 cout << "=== m_rpp420_b: " << m_rpp420_b << endl; 00059 cout << "=== lengthzdc: " << lengthzdc << endl; 00060 cout << "=== lengthd1: " << lengthd1 << endl; 00061 cout << "===================================================================" << endl; 00062 } 00063 edm::FileInPath b1(beam1filename.c_str()); 00064 edm::FileInPath b2(beam2filename.c_str()); 00065 00066 // construct beam line for FP420: . 00067 if(m_FP420Transport && lengthfp420>0. ) { 00068 m_beamlineFP4201 = new H_BeamLine( 1, lengthfp420 + 0.1 ); // (direction, length) 00069 m_beamlineFP4202 = new H_BeamLine( -1, lengthfp420 + 0.1 ); // 00070 try { 00071 m_beamlineFP4201->fill( b1.fullPath(), 1, "IP5" ); 00072 m_beamlineFP4202->fill( b2.fullPath(), -1, "IP5" ); 00073 } catch ( const edm::Exception& e ) { 00074 std::string msg = e.what(); 00075 msg += " caught in Hector... \nERROR: Could not locate SimTransport/HectorData data files."; 00076 edm::LogError ("DataNotFound") << msg; 00077 } 00078 m_beamlineFP4201->offsetElements( 120, -0.097 ); 00079 m_beamlineFP4202->offsetElements( 120, +0.097 ); 00080 m_beamlineFP4201->calcMatrix(); 00081 m_beamlineFP4202->calcMatrix(); 00082 } 00083 else{ 00084 cout << "=== Hector: WARNING: lengthfp420= " << lengthfp420 << endl; 00085 } 00086 00087 00088 if (m_ZDCTransport && lengthzdc>0. && lengthd1>0.) { 00089 // construct beam line for ZDC: . 00090 m_beamlineZDC1 = new H_BeamLine( 1, lengthzdc + 0.1 ); // (direction, length) 00091 m_beamlineZDC2 = new H_BeamLine( -1, lengthzdc + 0.1 ); // 00092 try { 00093 m_beamlineZDC1->fill( b1.fullPath(), 1, "IP5" ); 00094 m_beamlineZDC2->fill( b2.fullPath(), -1, "IP5" ); 00095 } catch ( const edm::Exception& e ) { 00096 std::string msg = e.what(); 00097 msg += " caught in Hector... \nERROR: Could not locate SimTransport/HectorData data files."; 00098 edm::LogError ("DataNotFound") << msg; 00099 } 00100 m_beamlineZDC1->offsetElements( 120, -0.097 ); 00101 m_beamlineZDC2->offsetElements( 120, +0.097 ); 00102 m_beamlineZDC1->calcMatrix(); 00103 m_beamlineZDC2->calcMatrix(); 00104 00105 00106 // construct beam line for D1: . 00107 m_beamlineD11 = new H_BeamLine( 1, lengthd1 + 0.1 ); // (direction, length) 00108 m_beamlineD12 = new H_BeamLine( -1, lengthd1 + 0.1 ); // 00109 try { 00110 m_beamlineD11->fill( b1.fullPath(), 1, "IP5" ); 00111 m_beamlineD12->fill( b2.fullPath(), -1, "IP5" ); 00112 } catch ( const edm::Exception& e ) { 00113 std::string msg = e.what(); 00114 msg += " caught in Hector... \nERROR: Could not locate SimTransport/HectorData data files."; 00115 edm::LogError ("DataNotFound") << msg; 00116 } 00117 m_beamlineD11->offsetElements( 120, -0.097 ); 00118 m_beamlineD12->offsetElements( 120, +0.097 ); 00119 m_beamlineD11->calcMatrix(); 00120 m_beamlineD12->calcMatrix(); 00121 } 00122 else{ 00123 cout << "=== Hector: WARNING: lengthzdc= " << lengthzdc << "lengthd1= " << lengthd1 << endl; 00124 } 00125 00126 edm::LogInfo ("Hector") << "===================================================================\n"; 00127 00128 }
Hector::~Hector | ( | ) | [virtual] |
Definition at line 130 of file Hector.cc.
References it, m_beamlineD11, m_beamlineD12, m_beamlineFP4201, m_beamlineFP4202, m_beamlineZDC1, m_beamlineZDC2, and m_beamPart.
00130 { 00131 00132 edm::LogInfo ("Hector") << "===================================================================\n" 00133 << "=== Start delete Hector =====\n"; 00134 for (std::map<unsigned int,H_BeamParticle*>::iterator it = m_beamPart.begin(); it != m_beamPart.end(); ++it ) { 00135 delete (*it).second; 00136 } 00137 // 00138 delete m_beamlineFP4201; 00139 delete m_beamlineFP4202; 00140 delete m_beamlineZDC1; 00141 delete m_beamlineZDC2; 00142 delete m_beamlineD11; 00143 delete m_beamlineD12; 00144 00145 edm::LogInfo ("Hector") << "===================================================================\n"; 00146 }
void Hector::add | ( | const HepMC::GenEvent * | ev, | |
const edm::EventSetup & | es | |||
) |
Adds the stable protons from the event ev to a beamline
Definition at line 177 of file Hector.cc.
References funct::abs(), GenMuonPlsPt100GeV_cfg::cout, lat::endl(), etacut, g, edm::EventSetup::getData(), parsecf::pyparsing::line(), m_beamPart, m_direct, m_eta, m_isCharged, m_pdg, m_pz, m_verbosity, pdt, and funct::sqrt().
Referenced by HectorProducer::produce().
00177 { 00178 00179 H_BeamParticle * h_p; 00180 unsigned int line; 00181 00182 // unsigned int npart = ev->nStable(); 00183 00184 /* for (HepMC::GenEvent::particle_const_iterator eventParticle =evt->particles_begin(); 00185 eventParticle != evt->particles_end(); 00186 eventParticle = find_if(++eventParticle, evt->particles_end(), ( (*eventParticle)->status() == 1 ) ) ) {*/ 00187 // for ( HepMC::GenVertex::particle_iterator eventParticle = (*vitr)->particles_begin(HepMC::descendants); 00188 for (HepMC::GenEvent::particle_const_iterator eventParticle =evt->particles_begin(); 00189 eventParticle != evt->particles_end(); 00190 ++eventParticle ) { 00191 // if ( (*eventParticle)->status() == 1 && !(*eventParticle)->end_vertex()) { 00192 if ( (*eventParticle)->status() == 1 ) { 00193 // if ( abs( (*eventParticle)->pdg_id() ) == 2212 ) { // if it's proton 00194 if ( abs( (*eventParticle)->momentum().eta())>etacut){ 00195 line = (*eventParticle)->barcode(); 00196 if ( m_beamPart.find(line) == m_beamPart.end() ) { 00197 double charge=1.; 00198 m_isCharged[line] = false;// neutrals 00199 HepMC::GenParticle * g = (*eventParticle); 00200 iSetup.getData( pdt ); 00201 const ParticleData * part = pdt->particle( g->pdg_id() ); 00202 if (part){ 00203 charge = part->charge(); 00204 if(m_verbosity) cout <<"charge= " << charge << endl; 00205 } 00206 if(charge !=0) m_isCharged[line] = true;//charged 00207 double mass = (*eventParticle)->generatedMass(); 00208 if(m_verbosity) cout << "=== Hector:add: status= "<< g->status() <<"mass = " << mass <<"charge= " << charge << endl; 00209 00210 // h_p = new H_BeamParticle(); 00211 00212 // H_BeamParticle myparticle(mass,charge); 00213 h_p = new H_BeamParticle(mass,charge); 00214 00215 double px,py,pz,pt; 00216 px = (*eventParticle)->momentum().px(); 00217 py = (*eventParticle)->momentum().py(); 00218 pz = (*eventParticle)->momentum().pz(); 00219 00220 h_p->set4Momentum( px, py, pz, (*eventParticle)->momentum().e() ); 00221 00222 pt = sqrt( (px*px) + (py*py) ); 00223 00224 // Clears H_BeamParticle::positions and sets the initial one 00225 h_p->setPosition( (*eventParticle)->production_vertex()->position().x(), (*eventParticle)->production_vertex()->position().y(), std::atan2( px, pt ), std::atan2( py, pt ), (*eventParticle)->production_vertex()->position().z() ); 00226 00228 00229 00230 m_beamPart[line] = h_p; 00231 m_direct[line] = 0; 00232 m_direct[line] = ( pz > 0 ) ? 1 : -1; 00233 00234 m_eta[line] = (*eventParticle)->momentum().eta(); 00235 m_pdg[line] = (*eventParticle)->pdg_id(); 00236 m_pz[line] = (*eventParticle)->momentum().pz(); 00237 00238 00239 if(m_verbosity) { 00240 cout << "=== Hector:add: pz= "<< pz << endl; 00241 cout << "=== Hector:add: pt= "<<pt << endl; 00242 cout << "=== Hector:add: m_isCharged[line]= "<< m_isCharged[line] << endl; 00243 } 00244 }// if find line 00245 }// if 2212 or eta 8.2 00246 }// if status 00247 }// for loop 00248 00249 }
HepMC::GenEvent * Hector::addPartToHepMC | ( | HepMC::GenEvent * | event | ) |
Return vector of the particle lines (HepMC::GenParticle::barcode()) in a beamline
Definition at line 495 of file Hector.cc.
References funct::cos(), GenMuonPlsPt100GeV_cfg::cout, lat::endl(), relval_parameters_module::energy, it, lengthd1, parsecf::pyparsing::line(), m_beamPart, m_direct, m_eAtTrPoint, m_FP420Transport, m_isStoppedd1, m_isStoppedfp420, m_isStoppedzdc, m_rpp420_b, m_rpp420_f, m_TxAtTrPoint, m_TyAtTrPoint, m_verbosity, m_xAtTrPoint, m_yAtTrPoint, m_ZDCTransport, pi, edm::second(), funct::sin(), funct::sqrt(), and theta.
Referenced by HectorProducer::produce().
00495 { 00496 00497 unsigned int line; 00498 // H_BeamParticle * part; 00499 HepMC::GenParticle * gpart; 00500 long double tx,ty,theta,fi,energy,time = 0; 00501 std::map< unsigned int, H_BeamParticle* >::iterator it; 00502 00503 00504 for (it = m_beamPart.begin(); it != m_beamPart.end(); ++it ) { 00505 line = (*it).first; 00506 if(!m_FP420Transport) m_isStoppedfp420[line] = true; 00507 if(!m_ZDCTransport) {m_isStoppedzdc[line] = false;m_isStoppedd1[line] = true;} 00508 if(m_verbosity) { 00509 cout << "=== Hector:addPartToHepMC: isStoppedfp420=" << (*m_isStoppedfp420.find(line)).second << endl; 00510 cout << "=== Hector:addPartToHepMC: isStoppedzdc=" << (*m_isStoppedzdc.find(line)).second << endl; 00511 cout << "=== Hector:addPartToHepMC: isStoppedd1=" << (*m_isStoppedd1.find(line)).second << endl; 00512 } 00513 if (!((*m_isStoppedfp420.find(line)).second) || (!((*m_isStoppedd1.find(line)).second) && ((*m_isStoppedzdc.find(line)).second))){ 00514 // if ( !((*m_isStoppedfp420.find(line)).second) || !((*m_isStoppedd1.find(line)).second) ) { 00515 gpart = evt->barcode_to_particle( line ); 00516 if ( gpart ) { 00517 tx = (*m_TxAtTrPoint.find(line)).second / 1000000.; 00518 ty = (*m_TyAtTrPoint.find(line)).second / 1000000.; 00519 theta = sqrt((tx*tx) + (ty*ty)); 00520 double ddd = 0.; 00521 if ( !((*m_isStoppedfp420.find(line)).second)) { 00522 if( (*m_direct.find( line )).second >0 ) { 00523 ddd = m_rpp420_f; 00524 } 00525 else if((*m_direct.find( line )).second <0 ) { 00526 ddd = m_rpp420_b; 00527 theta= pi-theta; 00528 } 00529 } 00530 else { 00531 ddd = lengthd1; 00532 if((*m_direct.find( line )).second <0 ) theta= pi-theta; 00533 } 00534 00535 fi = std::atan2(tx,ty); // tx, ty never == 0? 00536 energy = (*m_eAtTrPoint.find(line)).second; 00537 00538 // HepMC::GenEvent::vertex_iterator v_it; 00539 // time = ( *evt->vertices_begin() )->position().t(); // does time important? 00540 //long double time_buf = 0; 00541 //for (v_it = evt->vertices_begin(); v_it != evt->vertices_end(); ++v_it) { // since no operator-- 00542 // time_buf = ( *v_it )->position().t(); 00543 // time = ( time > time_buf ) ? time : time_buf; 00544 // } 00545 00546 // time = ( ddd*1000. - gpart->production_vertex()->position().z()*mm )/c_light; 00547 00548 time = ( ddd*1000. - gpart->production_vertex()->position().z()*mm ); // mm 00549 // decaylength/(p.Beta()*p.Gamma()*c_light); 00550 // time = ( ddd*1000. - gpart->production_vertex()->position().z()*mm )/300.; // nsec 00551 00552 if(ddd != 0.) { 00553 if(m_verbosity) { 00554 std::cout<<"=========Hector:: x= "<< (*(m_xAtTrPoint.find(line))).second*0.001<<std::endl; 00555 std::cout<<"=========Hector:: y= "<< (*(m_yAtTrPoint.find(line))).second*0.001<<std::endl; 00556 std::cout<<"=========Hector:: z= "<< ddd * (*(m_direct.find( line ))).second*1000.<<std::endl; 00557 std::cout<<"=========Hector:: t= "<< time <<std::endl; 00558 } 00559 00560 00561 00562 HepMC::GenVertex * vert = new HepMC::GenVertex( HepMC::FourVector( (*(m_xAtTrPoint.find(line))).second*0.001, 00563 (*(m_yAtTrPoint.find(line))).second*0.001, 00564 ddd * (*(m_direct.find( line ))).second*1000., 00565 time + .001*time ) ); 00566 00567 00568 00569 gpart->set_status( 2 ); 00570 vert->add_particle_in( gpart ); 00571 vert->add_particle_out( new HepMC::GenParticle( HepMC::FourVector(energy*std::sin(theta)*std::sin(fi), 00572 energy*std::sin(theta)*std::cos(fi), 00573 energy*std::cos(theta), 00574 energy ), 00575 gpart->pdg_id(), 1, gpart->flow() ) ); 00576 evt->add_vertex( vert ); 00577 00578 if(m_verbosity) std::cout << "Hector::TRANSPORTED pz= "<< gpart->momentum().pz() << " eta= "<< gpart->momentum().eta() << " status= "<< gpart->status() <<std::endl; 00579 00580 00581 }// ddd 00582 }// if gpart 00583 }// if !isStopped 00584 00585 else { 00586 gpart = evt->barcode_to_particle( line ); 00587 if ( gpart ) { 00588 HepMC::GenVertex * vert= new HepMC::GenVertex(); 00589 gpart->set_status( 2 ); 00590 vert->add_particle_in( gpart ); 00591 vert->add_particle_out( gpart ); 00592 evt->add_vertex( vert ); 00593 if(m_verbosity) std::cout << "Hector::NON-transp. pz= "<< gpart->momentum().pz() << " eta= "<< gpart->momentum().eta() << " status= "<< gpart->status() <<std::endl; 00594 } 00595 } 00596 00597 }//for 00598 // cout << "=== Hector:addPartToHepMC: end " << endl; 00599 00600 // if(m_verbosity){ 00601 // std::cout<<"============================================================================"<<std::endl; 00602 // std::cout<<"=========Hector::addPartToHepMC print:"<<std::endl; 00603 // evt->print(); 00604 //} 00605 00606 00607 00608 return evt; 00609 }
Clears BeamParticle, prepares Hector for a next Aperture check or/and a next event
Definition at line 154 of file Hector.cc.
References it, m_beamPart, m_direct, m_eta, m_isCharged, m_pdg, and m_pz.
Referenced by HectorProducer::produce().
00154 { 00155 for ( std::map<unsigned int,H_BeamParticle*>::iterator it = m_beamPart.begin(); it != m_beamPart.end(); ++it ) { 00156 delete (*it).second; 00157 }; 00158 m_beamPart.clear(); 00159 m_direct.clear(); 00160 m_eta.clear(); 00161 m_pdg.clear(); 00162 m_pz.clear(); 00163 m_isCharged.clear(); 00164 }
void Hector::clearApertureFlags | ( | ) |
Clears ApertureFlags, prepares Hector for a next event
Definition at line 148 of file Hector.cc.
References m_isStoppedd1, m_isStoppedfp420, and m_isStoppedzdc.
Referenced by HectorProducer::produce().
00148 { 00149 m_isStoppedfp420.clear(); 00150 m_isStoppedzdc.clear(); 00151 m_isStoppedd1.clear(); 00152 }
void Hector::filterD1 | ( | ) |
propagate the particles through a beamline to ZDC
Definition at line 402 of file Hector.cc.
References GenMuonPlsPt100GeV_cfg::cout, direction, lat::endl(), it, lengthd1, parsecf::pyparsing::line(), m_beamlineD11, m_beamlineD12, m_beamPart, m_direct, m_eAtTrPoint, m_eta, m_isCharged, m_isStoppedd1, m_isStoppedzdc, m_pdg, m_pz, m_sig_e, m_sigmaSTX, m_sigmaSTY, m_smearAng, m_smearE, m_TxAtTrPoint, m_TyAtTrPoint, m_verbosity, m_xAtTrPoint, m_yAtTrPoint, and edm::second().
Referenced by HectorProducer::produce().
00402 { 00403 unsigned int line; 00404 H_BeamParticle * part; 00405 std::map< unsigned int, H_BeamParticle* >::iterator it; 00406 00407 bool is_stop_d1; 00408 int direction; 00409 00410 float x1_d1; 00411 float y1_d1; 00412 00413 if ( m_beamPart.size() && lengthd1>0.) { 00414 00415 for (it = m_beamPart.begin(); it != m_beamPart.end(); ++it ) { 00416 line = (*it).first; 00417 part = (*it).second; 00418 if(m_verbosity) cout << "=== Hector:filterD1: isStopped=" << (*m_isStoppedzdc.find(line)).second << endl; 00419 if ( ((*m_isStoppedzdc.find(line)).second) || !((*m_isCharged.find( line )).second) ) { 00420 00421 if(m_verbosity) cout << "filterD1:pdg_id=" << (*m_pdg.find( line )).second << " momentum().eta=" << (*m_eta.find( line )).second<< " pz=" << (*m_pz.find( line )).second<< endl; 00422 00423 direction = (*m_direct.find( line )).second; 00424 if(m_verbosity) cout << "filterD1:direction=" << direction << endl; 00425 if (m_smearAng) { 00426 if ( m_sigmaSTX>0. && m_sigmaSTY>0.) { 00427 part->smearAng(m_sigmaSTX,m_sigmaSTY); 00428 } 00429 else { 00430 part->smearAng(); 00431 } 00432 } 00433 if (m_smearE) { 00434 if ( m_sig_e ) { 00435 part->smearE(m_sig_e); 00436 } 00437 else { 00438 part->smearE(); 00439 } 00440 } 00441 if ( direction == 1 ) { 00442 part->computePath( m_beamlineD11, 1 ); 00443 is_stop_d1 = part->stopped( m_beamlineD11 ); 00444 m_isStoppedd1[line] = is_stop_d1; 00445 if(m_verbosity) cout << "=== Hector:filterD1: pos is_stop_d1= "<< is_stop_d1 << endl; 00446 } 00447 else if ( direction == -1 ){ 00448 part->computePath( m_beamlineD12, -1 ); 00449 is_stop_d1 = part->stopped( m_beamlineD12 ); 00450 m_isStoppedd1[line] = is_stop_d1; 00451 if(m_verbosity) cout << "=== Hector:filterD1: neg is_stop_d1= "<< is_stop_d1 << endl; 00452 } 00453 else { 00454 is_stop_d1 = true; 00455 m_isStoppedd1[line] = is_stop_d1; 00456 if(m_verbosity) cout << "=== Hector:filterD1: 0 is_stop_d1= "<< endl; 00457 } 00458 //propagating 00459 if (!is_stop_d1 ) { 00460 part->propagate( lengthd1 ); 00461 x1_d1 = part->getX(); 00462 y1_d1 = part->getY(); 00463 m_xAtTrPoint[line] = x1_d1; 00464 m_yAtTrPoint[line] = y1_d1; 00465 m_TxAtTrPoint[line] = part->getTX(); 00466 m_TyAtTrPoint[line] = part->getTY(); 00467 m_eAtTrPoint[line] = part->getE(); 00468 } 00469 }// if stopzdc 00470 else { 00471 m_isStoppedd1[line] = false;// not stopped in propagating to ZDC and therefore in propagation to D1 too. 00472 if(m_verbosity) cout << "=== Hector:filterD1: isStopped=" << (*m_isStoppedd1.find(line)).second << endl; 00473 } 00474 00475 } // for (it = m_beamPart.begin(); it != m_beamPart.end(); it++ ) 00476 } // if ( m_beamPart.size() ) 00477 00478 }
void Hector::filterFP420 | ( | ) |
propagate the particles through a beamline to FP420
Definition at line 254 of file Hector.cc.
References GenMuonPlsPt100GeV_cfg::cout, direction, lat::endl(), it, lengthfp420, parsecf::pyparsing::line(), m_beamlineFP4201, m_beamlineFP4202, m_beamPart, m_direct, m_eAtTrPoint, m_eta, m_isCharged, m_isStoppedfp420, m_pdg, m_pz, m_rpp420_b, m_rpp420_f, m_sig_e, m_sigmaSTX, m_sigmaSTY, m_smearAng, m_smearE, m_TxAtTrPoint, m_TyAtTrPoint, m_verbosity, m_xAtTrPoint, m_yAtTrPoint, and edm::second().
Referenced by HectorProducer::produce().
00254 { 00255 unsigned int line; 00256 H_BeamParticle * part; 00257 std::map< unsigned int, H_BeamParticle* >::iterator it; 00258 00259 bool is_stop; 00260 int direction; 00261 00262 float x1_420; 00263 float y1_420; 00264 00265 if ( m_beamPart.size() && lengthfp420>0. ) { 00266 00267 for (it = m_beamPart.begin(); it != m_beamPart.end(); ++it ) { 00268 line = (*it).first; 00269 part = (*it).second; 00270 00271 if(m_verbosity) { 00272 cout << "filterFP420:pdg_id=" << (*m_pdg.find( line )).second << " momentum().eta=" << (*m_eta.find( line )).second<< " pz=" << (*m_pz.find( line )).second<< endl; 00273 } 00274 if ( (*m_isCharged.find( line )).second ) { 00275 direction = (*m_direct.find( line )).second; 00276 if (m_smearAng) { 00277 if ( m_sigmaSTX>0. && m_sigmaSTY>0.) { 00278 part->smearAng(m_sigmaSTX,m_sigmaSTY); 00279 } 00280 else { 00281 part->smearAng(); 00282 } 00283 } 00284 if (m_smearE) { 00285 if ( m_sig_e ) { 00286 part->smearE(m_sig_e); 00287 } 00288 else { 00289 part->smearE(); 00290 } 00291 } 00292 if ( direction == 1 ) { 00293 part->computePath( m_beamlineFP4201, 1 ); 00294 is_stop = part->stopped( m_beamlineFP4201 ); 00295 if(m_verbosity) cout << "=== Hector:filterFP420: pos is_stop= "<< is_stop << endl; 00296 } 00297 else if ( direction == -1 ){ 00298 part->computePath( m_beamlineFP4202, -1 ); 00299 is_stop = part->stopped( m_beamlineFP4202 ); 00300 if(m_verbosity) cout << "=== Hector:filterFP420: neg is_stop= "<< is_stop << endl; 00301 } 00302 else { 00303 is_stop = true; 00304 if(m_verbosity) cout << "=== Hector:filterFP420: 0 is_stop= "<< is_stop << endl; 00305 } 00306 00307 //propagating 00308 m_isStoppedfp420[line] = is_stop; 00309 if(m_verbosity) cout << "=== Hector:filterFP420: isStopped=" << (*m_isStoppedfp420.find(line)).second << endl; 00310 00311 if (!is_stop) { 00312 if ( direction == 1 ) part->propagate( m_rpp420_f ); 00313 if ( direction == -1 ) part->propagate( m_rpp420_b ); 00314 x1_420 = part->getX(); 00315 y1_420 = part->getY(); 00316 if(m_verbosity) cout << "=== Hector:filterFP420: x= "<< x1_420 <<" y= " << y1_420 << endl; 00317 00318 m_xAtTrPoint[line] = x1_420; 00319 m_yAtTrPoint[line] = y1_420; 00320 m_TxAtTrPoint[line] = part->getTX(); 00321 m_TyAtTrPoint[line] = part->getTY(); 00322 m_eAtTrPoint[line] = part->getE(); 00323 00324 } 00325 }// if isCharged 00326 else { 00327 m_isStoppedfp420[line] = true;// imply that neutral particles stopped to reach 420m 00328 if(m_verbosity) cout << "=== Hector:filterFP420: isStopped=" << (*m_isStoppedfp420.find(line)).second << endl; 00329 } 00330 00331 } // for (it = m_beamPart.begin(); it != m_beamPart.end(); it++ ) 00332 } // if ( m_beamPart.size() ) 00333 00334 }
void Hector::filterZDC | ( | ) |
propagate the particles through a beamline to ZDC
Definition at line 336 of file Hector.cc.
References GenMuonPlsPt100GeV_cfg::cout, direction, lat::endl(), it, lengthzdc, parsecf::pyparsing::line(), m_beamlineZDC1, m_beamlineZDC2, m_beamPart, m_direct, m_eta, m_isCharged, m_isStoppedfp420, m_isStoppedzdc, m_pdg, m_pz, m_sig_e, m_sigmaSTX, m_sigmaSTY, m_smearAng, m_smearE, m_verbosity, and edm::second().
Referenced by HectorProducer::produce().
00336 { 00337 unsigned int line; 00338 H_BeamParticle * part; 00339 std::map< unsigned int, H_BeamParticle* >::iterator it; 00340 00341 bool is_stop_zdc; 00342 int direction; 00343 00344 if ( m_beamPart.size() && lengthzdc>0. ) { 00345 00346 for (it = m_beamPart.begin(); it != m_beamPart.end(); ++it ) { 00347 line = (*it).first; 00348 part = (*it).second; 00349 if(m_verbosity) cout << "=== Hector:filterZDC: isStopped=" << (*m_isStoppedfp420.find(line)).second << endl; 00350 if ( ((*m_isStoppedfp420.find(line)).second) && ((*m_isCharged.find(line)).second) ) { 00351 00352 if(m_verbosity) cout << "filterZDC:pdg_id=" << (*m_pdg.find( line )).second << " momentum().eta=" << (*m_eta.find( line )).second<< " pz=" << (*m_pz.find( line )).second<< endl; 00353 00354 direction = (*m_direct.find( line )).second; 00355 if(m_verbosity) cout << "filterZDC:direction=" << direction << endl; 00356 if (m_smearAng) { 00357 if ( m_sigmaSTX>0. && m_sigmaSTY>0.) { 00358 part->smearAng(m_sigmaSTX,m_sigmaSTY); 00359 } 00360 else { 00361 part->smearAng(); 00362 } 00363 } 00364 if (m_smearE) { 00365 if ( m_sig_e ) { 00366 part->smearE(m_sig_e); 00367 } 00368 else { 00369 part->smearE(); 00370 } 00371 } 00372 if ( direction == 1 ) { 00373 part->computePath( m_beamlineZDC1, 1 ); 00374 is_stop_zdc = part->stopped( m_beamlineZDC1 ); 00375 m_isStoppedzdc[line] = is_stop_zdc; 00376 if(m_verbosity) cout << "=== Hector:filterZDC: pos is_stop_zdc= "<< is_stop_zdc << endl; 00377 } 00378 else if ( direction == -1 ){ 00379 part->computePath( m_beamlineZDC2, -1 ); 00380 is_stop_zdc = part->stopped( m_beamlineZDC2 ); 00381 m_isStoppedzdc[line] = is_stop_zdc; 00382 if(m_verbosity) cout << "=== Hector:filterZDC: neg is_stop_zdc= "<< is_stop_zdc << endl; 00383 } 00384 else { 00385 m_isStoppedzdc[line] = true; 00386 if(m_verbosity) cout << "=== Hector:filterZDC: 0 is_stop_zdc= "<< endl; 00387 } 00388 }// if stopfp420 charged particles 00389 else if ( ((*m_isCharged.find(line)).second) ){ 00390 m_isStoppedzdc[line] = false;// not stopped in propagating to FP420 and therefore in propagation to ZDC too. 00391 if(m_verbosity) cout << "=== Hector:filterZDC: isStopped=" << (*m_isStoppedzdc.find(line)).second << endl; 00392 } 00393 else { 00394 m_isStoppedzdc[line] = true;// neutrals particles considered as stopped in propagating to ZDC 00395 if(m_verbosity) cout << "=== Hector:filterZDC: isStopped=" << (*m_isStoppedzdc.find(line)).second << endl; 00396 } 00397 00398 } // for (it = m_beamPart.begin(); it != m_beamPart.end(); it++ ) 00399 } // if ( m_beamPart.size() ) 00400 00401 }
Prints properties of all particles in a beamline
Definition at line 488 of file Hector.cc.
References it, and m_beamPart.
00488 { 00489 for (std::map<unsigned int,H_BeamParticle*>::const_iterator it = m_beamPart.begin(); it != m_beamPart.end(); ++it ) { 00490 (*it).second->printProperties(); 00491 }; 00492 }
string Hector::beam1filename [private] |
string Hector::beam2filename [private] |
double Hector::etacut [private] |
double Hector::lengthd1 [private] |
double Hector::lengthfp420 [private] |
double Hector::lengthzdc [private] |
H_BeamLine* Hector::m_beamlineD11 [private] |
H_BeamLine* Hector::m_beamlineD12 [private] |
H_BeamLine* Hector::m_beamlineFP4201 [private] |
H_BeamLine* Hector::m_beamlineFP4202 [private] |
H_BeamLine* Hector::m_beamlineZDC1 [private] |
H_BeamLine* Hector::m_beamlineZDC2 [private] |
std::map<unsigned int, H_BeamParticle*> Hector::m_beamPart [private] |
Definition at line 105 of file Hector.h.
Referenced by add(), addPartToHepMC(), clear(), filterD1(), filterFP420(), filterZDC(), print(), and ~Hector().
std::map<unsigned int, int> Hector::m_direct [private] |
Definition at line 106 of file Hector.h.
Referenced by add(), addPartToHepMC(), clear(), filterD1(), filterFP420(), filterZDC(), and getDirect().
std::map<unsigned int, double> Hector::m_eAtTrPoint [private] |
Definition at line 114 of file Hector.h.
Referenced by addPartToHepMC(), filterD1(), and filterFP420().
std::map<unsigned int, double> Hector::m_eta [private] |
Definition at line 116 of file Hector.h.
Referenced by add(), clear(), filterD1(), filterFP420(), and filterZDC().
bool Hector::m_FP420Transport [private] |
std::map<unsigned int, bool> Hector::m_isCharged [private] |
Definition at line 119 of file Hector.h.
Referenced by add(), clear(), filterD1(), filterFP420(), and filterZDC().
std::map<unsigned int, bool> Hector::m_isStoppedd1 [private] |
Definition at line 109 of file Hector.h.
Referenced by addPartToHepMC(), clearApertureFlags(), and filterD1().
std::map<unsigned int, bool> Hector::m_isStoppedfp420 [private] |
Definition at line 107 of file Hector.h.
Referenced by addPartToHepMC(), clearApertureFlags(), filterFP420(), and filterZDC().
std::map<unsigned int, bool> Hector::m_isStoppedzdc [private] |
Definition at line 108 of file Hector.h.
Referenced by addPartToHepMC(), clearApertureFlags(), filterD1(), and filterZDC().
std::map<unsigned int, int> Hector::m_pdg [private] |
Definition at line 117 of file Hector.h.
Referenced by add(), clear(), filterD1(), filterFP420(), and filterZDC().
std::map<unsigned int, double> Hector::m_pz [private] |
Definition at line 118 of file Hector.h.
Referenced by add(), clear(), filterD1(), filterFP420(), and filterZDC().
H_RecRPObject* Hector::m_rp420_b [private] |
H_RecRPObject* Hector::m_rp420_f [private] |
float Hector::m_rpp420_b [private] |
Definition at line 86 of file Hector.h.
Referenced by addPartToHepMC(), filterFP420(), and Hector().
float Hector::m_rpp420_f [private] |
Definition at line 85 of file Hector.h.
Referenced by addPartToHepMC(), filterFP420(), and Hector().
float Hector::m_rppd1 [private] |
float Hector::m_rppzdc [private] |
double Hector::m_sig_e [private] |
Definition at line 80 of file Hector.h.
Referenced by filterD1(), filterFP420(), filterZDC(), and Hector().
double Hector::m_sigmaSTX [private] |
Definition at line 82 of file Hector.h.
Referenced by filterD1(), filterFP420(), filterZDC(), and Hector().
double Hector::m_sigmaSTY [private] |
Definition at line 83 of file Hector.h.
Referenced by filterD1(), filterFP420(), filterZDC(), and Hector().
bool Hector::m_smearAng [private] |
Definition at line 79 of file Hector.h.
Referenced by filterD1(), filterFP420(), filterZDC(), and Hector().
bool Hector::m_smearE [private] |
Definition at line 81 of file Hector.h.
Referenced by filterD1(), filterFP420(), filterZDC(), and Hector().
std::map<unsigned int, double> Hector::m_TxAtTrPoint [private] |
Definition at line 112 of file Hector.h.
Referenced by addPartToHepMC(), filterD1(), and filterFP420().
std::map<unsigned int, double> Hector::m_TyAtTrPoint [private] |
Definition at line 113 of file Hector.h.
Referenced by addPartToHepMC(), filterD1(), and filterFP420().
bool Hector::m_verbosity [private] |
Definition at line 124 of file Hector.h.
Referenced by add(), addPartToHepMC(), filterD1(), filterFP420(), filterZDC(), and Hector().
std::map<unsigned int, double> Hector::m_xAtTrPoint [private] |
Definition at line 110 of file Hector.h.
Referenced by addPartToHepMC(), filterD1(), and filterFP420().
std::map<unsigned int, double> Hector::m_yAtTrPoint [private] |
Definition at line 111 of file Hector.h.
Referenced by addPartToHepMC(), filterD1(), and filterFP420().
bool Hector::m_ZDCTransport [private] |
edm::ESHandle< ParticleDataTable > Hector::pdt [private] |