// // ******************************************************************** // * License and Disclaimer * // * * // * The Geant4 software is copyright of the Copyright Holders of * // * the Geant4 Collaboration. It is provided under the terms and * // * conditions of the Geant4 Software License, included in the file * // * LICENSE and available at http://cern.ch/geant4/license . These * // * include a list of copyright holders. * // * * // * 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. Please see the license in the file LICENSE and URL above * // * for the full disclaimer and the limitation of liability. * // * * // * This code implementation is the result of the scientific and * // * technical work of the GEANT4 collaboration. * // * By using, copying, modifying or distributing the software (or * // * any work based on the software) you agree to acknowledge its * // * use in resulting scientific publications, and indicate your * // * acceptance of all terms of the Geant4 Software license. * // ******************************************************************** // // // $Id: G4SmoothTrajectory.cc 69003 2013-04-15 09:25:23Z gcosmo $ // // // --------------------------------------------------------------- // // G4SmoothTrajectory.cc // // Contact: // Questions and comments to this code should be sent to // Katsuya Amako (e-mail: Katsuya.Amako@kek.jp) // Makoto Asai (e-mail: asai@kekvax.kek.jp) // Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp) // // --------------------------------------------------------------- #include "G4SmoothTrajectory.hh" #include "G4SmoothTrajectoryPoint.hh" #include "G4ParticleTable.hh" #include "G4AttDefStore.hh" #include "G4AttDef.hh" #include "G4AttValue.hh" #include "G4UIcommand.hh" #include "G4UnitsTable.hh" //#define G4ATTDEBUG #ifdef G4ATTDEBUG #include "G4AttCheck.hh" #endif G4ThreadLocal G4Allocator *aSmoothTrajectoryAllocator = 0; G4SmoothTrajectory::G4SmoothTrajectory() : positionRecord(0), fTrackID(0), fParentID(0), PDGEncoding( 0 ), PDGCharge(0.0), ParticleName(""), initialKineticEnergy( 0. ), initialMomentum( G4ThreeVector() ) {;} G4SmoothTrajectory::G4SmoothTrajectory(const G4Track* aTrack) { G4ParticleDefinition * fpParticleDefinition = aTrack->GetDefinition(); ParticleName = fpParticleDefinition->GetParticleName(); PDGCharge = fpParticleDefinition->GetPDGCharge(); PDGEncoding = fpParticleDefinition->GetPDGEncoding(); fTrackID = aTrack->GetTrackID(); fParentID = aTrack->GetParentID(); initialKineticEnergy = aTrack->GetKineticEnergy(); initialMomentum = aTrack->GetMomentum(); positionRecord = new TrajectoryPointContainer(); // Following is for the first trajectory point positionRecord->push_back(new G4SmoothTrajectoryPoint(aTrack->GetPosition())); // The first point has no auxiliary points, so set the auxiliary // points vector to NULL positionRecord->push_back(new G4SmoothTrajectoryPoint(aTrack->GetPosition(), 0)); } G4SmoothTrajectory::G4SmoothTrajectory(G4SmoothTrajectory & right):G4VTrajectory() { ParticleName = right.ParticleName; PDGCharge = right.PDGCharge; PDGEncoding = right.PDGEncoding; fTrackID = right.fTrackID; fParentID = right.fParentID; initialKineticEnergy = right.initialKineticEnergy; initialMomentum = right.initialMomentum; positionRecord = new TrajectoryPointContainer(); for(size_t i=0;isize();i++) { G4SmoothTrajectoryPoint* rightPoint = (G4SmoothTrajectoryPoint*)((*(right.positionRecord))[i]); positionRecord->push_back(new G4SmoothTrajectoryPoint(*rightPoint)); } } G4SmoothTrajectory::~G4SmoothTrajectory() { if (positionRecord) { size_t i; for(i=0;isize();i++) { delete (*positionRecord)[i]; } positionRecord->clear(); delete positionRecord; } } void G4SmoothTrajectory::ShowTrajectory(std::ostream& os) const { // Invoke the default implementation in G4VTrajectory... G4VTrajectory::ShowTrajectory(os); // ... or override with your own code here. } void G4SmoothTrajectory::DrawTrajectory() const { // Invoke the default implementation in G4VTrajectory... G4VTrajectory::DrawTrajectory(); // ... or override with your own code here. } const std::map* G4SmoothTrajectory::GetAttDefs() const { G4bool isNew; std::map* store = G4AttDefStore::GetInstance("G4SmoothTrajectory",isNew); if (isNew) { G4String ID("ID"); (*store)[ID] = G4AttDef(ID,"Track ID","Physics","","G4int"); G4String PID("PID"); (*store)[PID] = G4AttDef(PID,"Parent ID","Physics","","G4int"); G4String PN("PN"); (*store)[PN] = G4AttDef(PN,"Particle Name","Physics","","G4String"); G4String Ch("Ch"); (*store)[Ch] = G4AttDef(Ch,"Charge","Physics","e+","G4double"); G4String PDG("PDG"); (*store)[PDG] = G4AttDef(PDG,"PDG Encoding","Physics","","G4int"); G4String IKE("IKE"); (*store)[IKE] = G4AttDef(IKE, "Initial kinetic energy", "Physics","G4BestUnit","G4double"); G4String IMom("IMom"); (*store)[IMom] = G4AttDef(IMom, "Initial momentum", "Physics","G4BestUnit","G4ThreeVector"); G4String IMag("IMag"); (*store)[IMag] = G4AttDef (IMag, "Initial momentum magnitude", "Physics","G4BestUnit","G4double"); G4String NTP("NTP"); (*store)[NTP] = G4AttDef(NTP,"No. of points","Physics","","G4int"); } return store; } std::vector* G4SmoothTrajectory::CreateAttValues() const { std::vector* values = new std::vector; values->push_back (G4AttValue("ID",G4UIcommand::ConvertToString(fTrackID),"")); values->push_back (G4AttValue("PID",G4UIcommand::ConvertToString(fParentID),"")); values->push_back(G4AttValue("PN",ParticleName,"")); values->push_back (G4AttValue("Ch",G4UIcommand::ConvertToString(PDGCharge),"")); values->push_back (G4AttValue("PDG",G4UIcommand::ConvertToString(PDGEncoding),"")); values->push_back (G4AttValue("IKE",G4BestUnit(initialKineticEnergy,"Energy"),"")); values->push_back (G4AttValue("IMom",G4BestUnit(initialMomentum,"Energy"),"")); values->push_back (G4AttValue("IMag",G4BestUnit(initialMomentum.mag(),"Energy"),"")); values->push_back (G4AttValue("NTP",G4UIcommand::ConvertToString(GetPointEntries()),"")); #ifdef G4ATTDEBUG G4cout << G4AttCheck(values,GetAttDefs()); #endif return values; } void G4SmoothTrajectory::AppendStep(const G4Step* aStep) { positionRecord->push_back( new G4SmoothTrajectoryPoint(aStep->GetPostStepPoint()->GetPosition(), aStep->GetPointerToVectorOfAuxiliaryPoints())); } G4ParticleDefinition* G4SmoothTrajectory::GetParticleDefinition() { return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName)); } void G4SmoothTrajectory::MergeTrajectory(G4VTrajectory* secondTrajectory) { if(!secondTrajectory) return; G4SmoothTrajectory* seco = (G4SmoothTrajectory*)secondTrajectory; G4int ent = seco->GetPointEntries(); for(G4int i=1;ipush_back((*(seco->positionRecord))[i]); } delete (*seco->positionRecord)[0]; seco->positionRecord->clear(); }