00001
00002
00003
00004
00005 #include "MollerDetectorHit.hh"
00006 #include "G4AttCheck.hh"
00007
00008 #include "G4ios.hh"
00009 #include "G4VVisManager.hh"
00010 #include "G4Circle.hh"
00011 #include "G4Colour.hh"
00012 #include "G4AttDefStore.hh"
00013 #include "G4AttDef.hh"
00014 #include "G4AttValue.hh"
00015 #include "G4UIcommand.hh"
00016 #include "G4UnitsTable.hh"
00017 #include "G4VisAttributes.hh"
00018 #include "G4LogicalVolume.hh"
00019
00020 G4Allocator<MollerDetectorHit> MollerDetectorHitAllocator;
00021
00022 MollerDetectorHit::MollerDetectorHit()
00023 {
00024 layerID = -1;
00025 trackID = -1;
00026
00027 worldPos = G4ThreeVector(0,0,0);
00028 vertexPos = G4ThreeVector(0,0,0);
00029
00030 kineE=0;
00031 kineE0=0;
00032
00033 creatorProcess = "";
00034 particleName = "";
00035 ion = 0;
00036 }
00037
00038 MollerDetectorHit::MollerDetectorHit(G4int z)
00039 {
00040 layerID = z;
00041 trackID = -1;
00042
00043 worldPos = G4ThreeVector(0,0,0);
00044 vertexPos = G4ThreeVector(0,0,0);
00045
00046 kineE=0;
00047 kineE0=0;
00048
00049 creatorProcess = "";
00050 particleName = "";
00051 ion = 0;
00052 }
00053
00054 MollerDetectorHit::~MollerDetectorHit()
00055 {;}
00056
00057 MollerDetectorHit::MollerDetectorHit(const MollerDetectorHit &right)
00058 : G4VHit()
00059 {
00060
00061 layerID = right.layerID;
00062 trackID = right.trackID;
00063 worldPos = right.worldPos;
00064 vertexPos = right.vertexPos;
00065
00066 kineE=right.kineE;
00067 kineE0=right.kineE0;
00068
00069 creatorProcess = right.creatorProcess;
00070 particleName = right.particleName;
00071 ion = right.ion;
00072 }
00073
00074 const MollerDetectorHit& MollerDetectorHit::operator=(const MollerDetectorHit &right)
00075 {
00076
00077 layerID = right.layerID;
00078 trackID = right.trackID;
00079 worldPos = right.worldPos;
00080 vertexPos = right.vertexPos;
00081
00082 kineE=right.kineE;
00083 kineE0=right.kineE0;
00084
00085 creatorProcess = right.creatorProcess;
00086 particleName = right.particleName;
00087 ion = right.ion;
00088
00089 return *this;
00090 }
00091
00092 int MollerDetectorHit::operator==(const MollerDetectorHit &) const
00093 {
00094 return 0;
00095 }
00096
00097
00098
00099 void MollerDetectorHit::Draw()
00100 {
00101 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
00102 if(pVVisManager)
00103 {
00104 G4Circle circle(worldPos);
00105 circle.SetScreenSize(0.04);
00106 circle.SetFillStyle(G4Circle::filled);
00107 G4Colour colour(1.,0.,0.);
00108 G4VisAttributes attribs(colour);
00109 circle.SetVisAttributes(attribs);
00110 pVVisManager->Draw(circle);
00111 }
00112 }
00113
00114 const std::map<G4String,G4AttDef>* MollerDetectorHit::GetAttDefs() const
00115 {
00116 G4bool isNew;
00117 std::map<G4String,G4AttDef>* store
00118 = G4AttDefStore::GetInstance("MollerDetectorHit",isNew);
00119
00120 if (isNew) {
00121 G4String HitType("HitType");
00122 (*store)[HitType] = G4AttDef(HitType,"Hit Type","Physics","","G4String");
00123
00124 G4String ID("ID");
00125 (*store)[ID] = G4AttDef(ID,"ID","Physics","","G4int");
00126
00127 G4String trackID("trackID");
00128 (*store)[trackID] = G4AttDef(trackID,"trackID","Physics","","G4int");
00129
00130 G4String Pos("Pos");
00131 (*store)[Pos] = G4AttDef(Pos, "Position",
00132 "Physics","G4BestUnit","G4ThreeVector");
00133
00134 G4String Creator("Creator");
00135 (*store)[Creator] = G4AttDef(Creator, "CreatorProcess",
00136 "Physics","","G4String");
00137
00138 G4String Name("Name");
00139 (*store)[Name] = G4AttDef(Name, "ParticleName",
00140 "Physics","","G4String");
00141
00142 G4String Ion("Ion");
00143 (*store)[Ion] = G4AttDef(Ion, "IsMoller",
00144 "Physics","","G4bool");
00145 }
00146 return store;
00147 }
00148
00149 std::vector<G4AttValue>* MollerDetectorHit::CreateAttValues() const
00150 {
00151 std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
00152
00153 values->push_back(G4AttValue("HitType","MolleDetectorHit",""));
00154
00155 values->push_back
00156 (G4AttValue("ID",G4UIcommand::ConvertToString(layerID),""));
00157
00158 values->push_back
00159 (G4AttValue("trackID",G4UIcommand::ConvertToString(trackID),""));
00160
00161 values->push_back
00162 (G4AttValue("Pos",G4BestUnit(worldPos,"Length"),""));
00163
00164 values->push_back
00165 (G4AttValue("Creator",creatorProcess,""));
00166
00167 values->push_back
00168 (G4AttValue("Name",particleName,""));
00169
00170 values->push_back
00171 (G4AttValue("Ion",ion,""));
00172
00173 return values;
00174
00175 }
00176
00177 void MollerDetectorHit::Print()
00178 {
00179 G4cout << " Layer[" << layerID << "] : Particle type: " << particleName
00180 << " created by " << creatorProcess << " \nat (x,y,z) "
00181 << vertexPos.x() << ", " << vertexPos.y() << ", " << vertexPos.z()
00182 << " with KineE0= " << kineE0/GeV << "[GeV]"
00183 << " hit detector \nat (x,y,z) " << worldPos.x()
00184 << ", " << worldPos.y() << ", " << worldPos.z()
00185 << " with KineE= " << kineE/GeV << "[GeV]\n"
00186 << G4endl;
00187 }
00188
00189