// // ******************************************************************** // * DISCLAIMER * // * * // * The following disclaimer summarizes all the specific disclaimers * // * of contributors to this software. The specific disclaimers,which * // * govern, are listed with their locations in: * // * http://cern.ch/geant4/license * // * * // * Neither the authors of this software system, nor their employing * // * institutes,nor the agencies providing financial support for this * // * work make any representation or warranty, express or implied, * // * regarding this software system or assume any liability for its * // * use. * // * * // * This code implementation is the intellectual property of the * // * GEANT4 collaboration. * // * By copying, distributing or modifying the Program (or any work * // * based on the Program) you indicate your acceptance of this * // * statement, and all its terms. * // ******************************************************************** // // #include "KM3PrimaryGeneratorAction.hh" #include "globals.hh" #include "G4Event.hh" #include "G4HEPEvtInterface.hh" #include "G4PrimaryParticle.hh" #include "G4PrimaryVertex.hh" #include "Randomize.hh" #include "G4ParticleDefinition.hh" #include "G4ParticleTypes.hh" #include "G4SystemOfUnits.hh" KM3PrimaryGeneratorAction::KM3PrimaryGeneratorAction() {} KM3PrimaryGeneratorAction::~KM3PrimaryGeneratorAction() {} void KM3PrimaryGeneratorAction::Initialize() { nevents = evtAANET->GetNumberOfEvents(); useHEPEvt = true; } // primary particle generation. Single, multiple (many vertexes) particles and // neutrino interaction events (single vertex) are supported // that covers almost everything, except exotic particles (monopoles etc) void KM3PrimaryGeneratorAction::GeneratePrimaries(G4Event *anEvent) { static G4int ievent = 0; G4int idneu; // type of neutrino interacting (PDG Code) G4int idtarget; // type of target if neutrino interaction (PDG Code) G4double xneu, yneu, zneu; // neutrino vertex (cm) if neutrino interaction G4double pxneu, pyneu, pzneu; // neutrino momentum (GeV/c) if neutrino interaction G4double xx0, yy0, zz0; // Vertex of injected or produced particles (in cm) G4double pxx0, pyy0, pzz0; // Momentum of injected or produced particles (in GeV/c) G4double t0; // initial time of injected particles(ns) ievent++; event_action->Initialize(); evtAANET->ReadEvent(); bool isneutrino = evtAANET->GetNeutrinoInfo(idneu, idtarget, xneu, yneu, zneu, pxneu, pyneu, pzneu, t0); evtAANET->GeneratePrimaryVertex(anEvent); G4int numberofVertices = anEvent->GetNumberOfPrimaryVertex(); G4int ipart = 0; numberofParticles = 0; for (G4int iv = 0; iv < numberofVertices; iv++) { numberofParticles += anEvent->GetPrimaryVertex(iv)->GetNumberOfParticle(); xx0 = anEvent->GetPrimaryVertex(iv)->GetX0(); yy0 = anEvent->GetPrimaryVertex(iv)->GetY0(); zz0 = anEvent->GetPrimaryVertex(iv)->GetZ0(); t0 = anEvent->GetPrimaryVertex(iv)->GetT0(); if(verbosity >2){ G4cout<<"Genereting vertex at ("<GetPrimaryVertex(iv)->Print(); } for (G4int ip = 0; ip < anEvent->GetPrimaryVertex(iv)->GetNumberOfParticle(); ip++) { idbeam = anEvent->GetPrimaryVertex(iv)->GetPrimary(ip)->GetPDGcode(); pxx0 = anEvent->GetPrimaryVertex(iv)->GetPrimary(ip)->GetPx() / GeV; pyy0 = anEvent->GetPrimaryVertex(iv)->GetPrimary(ip)->GetPy() / GeV; pzz0 = anEvent->GetPrimaryVertex(iv)->GetPrimary(ip)->GetPz() / GeV; if(verbosity >2){ G4cout<<"Particle "<AddPrimaryNumber(ipart + 1); } ipart++; } } if(verbosity>1 && isneutrino) G4cout << "Generating one neutrino event: evnum= " << ievent << " neutrinotype= " << idneu << " target= " << idtarget << " energy= " << sqrt(pxneu * pxneu + pyneu * pyneu + pzneu * pzneu) << " GeV" << G4endl; myTracking->numofInitialParticles = numberofParticles; if(verbosity>1) G4cout << "Number of particles in the event: "<