Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4PSDoseDeposit.cc
이 파일의 문서화 페이지로 가기
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4PSDoseDeposit.cc 81087 2014-05-20 15:44:27Z gcosmo $
28 // GEANT4 tag $Name: geant4-09-04 $
29 //
30 // G4PSDoseDeposit
31 #include "G4PSDoseDeposit.hh"
32 #include "G4VSolid.hh"
33 #include "G4VPhysicalVolume.hh"
34 #include "G4VPVParameterisation.hh"
35 #include "G4UnitsTable.hh"
36 
38 // (Description)
39 // This is a primitive scorer class for scoring only energy deposit.
40 //
41 //
42 // Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
43 // 2010-07-22 Introduce Unit specification.
44 //
46 
48  :G4VPrimitiveScorer(name,depth),HCID(-1),EvtMap(0)
49 {
50  SetUnit("Gy");
51 }
52 
54  G4int depth)
55  :G4VPrimitiveScorer(name,depth),HCID(-1),EvtMap(0)
56 {
57  SetUnit(unit);
58 }
59 
61 {;}
62 
64 {
66  if ( edep == 0. ) return FALSE;
67 
68  G4int idx = ((G4TouchableHistory*)
69  (aStep->GetPreStepPoint()->GetTouchable()))
70  ->GetReplicaNumber(indexDepth);
71  G4double cubicVolume = ComputeVolume(aStep, idx);
72 
73 
74  G4double density = aStep->GetTrack()->GetStep()->GetPreStepPoint()->GetMaterial()->GetDensity();
75  G4double dose = edep / ( density * cubicVolume );
76  dose *= aStep->GetPreStepPoint()->GetWeight();
77  G4int index = GetIndex(aStep);
78  EvtMap->add(index,dose);
79  return TRUE;
80 }
81 
83 {
85  GetName());
86  if(HCID < 0) {HCID = GetCollectionID(0);}
88 }
89 
91 {;}
92 
94 {
95  EvtMap->clear();
96 }
97 
99 {;}
100 
102 {
103  G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
104  G4cout << " PrimitiveScorer " << GetName() << G4endl;
105  G4cout << " Number of entries " << EvtMap->entries() << G4endl;
106  std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
107  for(; itr != EvtMap->GetMap()->end(); itr++) {
108  G4cout << " copy no.: " << itr->first
109  << " dose deposit: "
110  << *(itr->second)/GetUnitValue()
111  << " ["<<GetUnit() <<"]"
112  << G4endl;
113  }
114 }
115 
117 {
118  CheckAndSetUnit(unit,"Dose");
119 }
120 
122 
123  G4VPhysicalVolume* physVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
124  G4VPVParameterisation* physParam = physVol->GetParameterisation();
125  G4VSolid* solid = 0;
126  if(physParam)
127  { // for parameterized volume
128  if(idx<0)
129  {
131  ED << "Incorrect replica number --- GetReplicaNumber : " << idx << G4endl;
132  G4Exception("G4PSDoseDeposit::ComputeVolume","DetPS0004",JustWarning,ED);
133  }
134  solid = physParam->ComputeSolid(idx, physVol);
135  solid->ComputeDimensions(physParam,idx,physVol);
136  }
137  else
138  { // for ordinary volume
139  solid = physVol->GetLogicalVolume()->GetSolid();
140  }
141 
142  return solid->GetCubicVolume();
143 }
144 
145 
146 
const G4String & GetUnit() const
const XML_Char * name
Definition: expat.h:151
const G4VTouchable * GetTouchable() const
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
G4int GetCollectionID(G4int)
G4LogicalVolume * GetLogicalVolume() const
G4String GetName() const
virtual void EndOfEvent(G4HCofThisEvent *)
G4StepPoint * GetPreStepPoint() const
Float_t dose
G4THitsMap< G4double > * EvtMap
G4int add(const G4int &key, U *&aHit) const
Definition: G4THitsMap.hh:150
virtual G4VPVParameterisation * GetParameterisation() const =0
#define G4endl
Definition: G4ios.hh:61
virtual ~G4PSDoseDeposit()
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4VPhysicalVolume * GetPhysicalVolume() const
virtual G4double ComputeVolume(G4Step *, G4int idx)
virtual G4int GetIndex(G4Step *)
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
G4PSDoseDeposit(G4String name, G4int depth=0)
virtual void PrintAll()
G4MultiFunctionalDetector * GetMultiFunctionalDetector() const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4Material * GetMaterial() const
virtual void SetUnit(const G4String &unit)
#define FALSE
Definition: globals.hh:52
virtual void DrawAll()
void clear()
Definition: G4THitsMap.hh:510
G4Track * GetTrack() const
#define TRUE
Definition: globals.hh:55
G4double GetTotalEnergyDeposit() const
Definition: G4Step.hh:76
const G4Step * GetStep() const
Double_t edep
void CheckAndSetUnit(const G4String &unit, const G4String &category)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
int G4int
Definition: G4Types.hh:78
virtual void Initialize(G4HCofThisEvent *)
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
G4GLOB_DLL std::ostream G4cout
G4VSolid * GetSolid() const
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137
G4double GetWeight() const
G4MultiFunctionalDetector * detector
G4double GetUnitValue() const
G4int entries() const
Definition: G4THitsMap.hh:121
Map_t * GetMap() const
Definition: G4THitsMap.hh:117
G4double GetDensity() const
Definition: G4Material.hh:181
virtual void clear()