source: trunk/source/digits_hits/scorer/src/G4PSSphereSurfaceFlux.cc@ 1290

Last change on this file since 1290 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 8.2 KB
Line 
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: G4PSSphereSurfaceFlux.cc,v 1.3 2009/11/14 00:01:13 asaim Exp $
28// GEANT4 tag $Name: geant4-09-03 $
29//
30// G4PSSphereSurfaceFlux
31#include "G4PSSphereSurfaceFlux.hh"
32#include "G4StepStatus.hh"
33#include "G4Track.hh"
34#include "G4VSolid.hh"
35#include "G4VPhysicalVolume.hh"
36#include "G4VPVParameterisation.hh"
37#include "G4UnitsTable.hh"
38#include "G4GeometryTolerance.hh"
39////////////////////////////////////////////////////////////////////////////////
40// (Description)
41// This is a primitive scorer class for scoring only Surface Flux.
42// Flux version assumes only for G4Sphere shape.
43//
44// Surface is defined at the inside of sphere.
45// Direction -Rmin +Rmax
46// 0 IN || OUT ->|<- |
47// 1 IN ->| |
48// 2 OUT |<- |
49//
50// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
51// 29-Mar-2007 T.Aso, Bug fix for momentum direction at outgoing flux.
52//
53///////////////////////////////////////////////////////////////////////////////
54
55G4PSSphereSurfaceFlux::G4PSSphereSurfaceFlux(G4String name,
56 G4int direction, G4int depth)
57 :G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction)
58{;}
59
60G4PSSphereSurfaceFlux::~G4PSSphereSurfaceFlux()
61{;}
62
63G4bool G4PSSphereSurfaceFlux::ProcessHits(G4Step* aStep,G4TouchableHistory*)
64{
65 G4StepPoint* preStep = aStep->GetPreStepPoint();
66 G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
67 G4VPVParameterisation* physParam = physVol->GetParameterisation();
68 G4VSolid * solid = 0;
69 if(physParam)
70 { // for parameterized volume
71 G4int idx = ((G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable()))
72 ->GetReplicaNumber(indexDepth);
73 solid = physParam->ComputeSolid(idx, physVol);
74 solid->ComputeDimensions(physParam,idx,physVol);
75 }
76 else
77 { // for ordinary volume
78 solid = physVol->GetLogicalVolume()->GetSolid();
79 }
80
81// if( solid->GetEntityType() != "G4Sphere" ){
82// G4Exception("G4PSSphereSurfaceFluxScorer. - Solid type is not supported.");
83// return FALSE;
84// }
85 G4Sphere* sphereSolid = (G4Sphere*)(solid);
86
87 G4int dirFlag =IsSelectedSurface(aStep,sphereSolid);
88 if ( dirFlag > 0 ) {
89 if ( fDirection == fFlux_InOut || fDirection == dirFlag ){
90
91 G4StepPoint* thisStep=0;
92 if ( dirFlag == fFlux_In ){
93 thisStep = preStep;
94 }else if ( dirFlag == fFlux_Out ){
95 thisStep = aStep->GetPreStepPoint();
96 }else{
97 return FALSE;
98 }
99
100 G4TouchableHandle theTouchable = thisStep->GetTouchableHandle();
101 G4ThreeVector pdirection = thisStep->GetMomentumDirection();
102 G4ThreeVector localdir =
103 theTouchable->GetHistory()->GetTopTransform().TransformAxis(pdirection);
104 G4double localdirL2 = localdir.x()*localdir.x()
105 +localdir.y()*localdir.y()
106 +localdir.z()*localdir.z();
107 G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
108 G4ThreeVector localpos1 =
109 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
110 G4double localR2 = localpos1.x()*localpos1.x()
111 +localpos1.y()*localpos1.y()
112 +localpos1.z()*localpos1.z();
113 G4double anglefactor = (localdir.x()*localpos1.x()
114 +localdir.y()*localpos1.y()
115 +localdir.z()*localpos1.z())
116 /std::sqrt(localdirL2)/std::sqrt(localR2);
117
118 G4double radi = sphereSolid->GetInsideRadius();
119 G4double dph = sphereSolid->GetDeltaPhiAngle()/radian;
120 G4double stth = sphereSolid->GetStartThetaAngle()/radian;
121 G4double enth = stth+sphereSolid->GetDeltaThetaAngle()/radian;
122 G4double square = radi*radi*dph*( -std::cos(enth) + std::cos(stth) );
123
124 G4double current = thisStep->GetWeight(); // Flux (Particle Weight)
125 current = current/square; // Flux with angle.
126
127 current /= anglefactor;
128
129 G4int index = GetIndex(aStep);
130 EvtMap->add(index,current);
131 }
132 }
133
134 return TRUE;
135}
136
137G4int G4PSSphereSurfaceFlux::IsSelectedSurface(G4Step* aStep, G4Sphere* sphereSolid){
138
139 G4TouchableHandle theTouchable =
140 aStep->GetPreStepPoint()->GetTouchableHandle();
141 G4double kCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
142
143 if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
144 // Entering Geometry
145 G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
146 G4ThreeVector localpos1 =
147 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
148 G4double localR2 = localpos1.x()*localpos1.x()
149 +localpos1.y()*localpos1.y()
150 +localpos1.z()*localpos1.z();
151 //G4double InsideRadius2 =
152 // sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
153 //if(std::fabs( localR2 - InsideRadius2 ) < kCarTolerance ){
154 G4double InsideRadius = sphereSolid->GetInsideRadius();
155 if ( localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
156 &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
157 return fFlux_In;
158 }
159 }
160
161 if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
162 // Exiting Geometry
163 G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
164 G4ThreeVector localpos2 =
165 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
166 G4double localR2 = localpos2.x()*localpos2.x()
167 +localpos2.y()*localpos2.y()
168 +localpos2.z()*localpos2.z();
169 //G4double InsideRadius2 =
170 // sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
171 //if(std::facb(localR2 - InsideRadius2) ) < kCarTolerance ){
172 G4double InsideRadius = sphereSolid->GetInsideRadius();
173 if ( localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
174 &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
175 return fFlux_Out;
176 }
177 }
178
179 return -1;
180}
181
182void G4PSSphereSurfaceFlux::Initialize(G4HCofThisEvent* HCE)
183{
184 EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName());
185 if ( HCID < 0 ) HCID = GetCollectionID(0);
186 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
187}
188
189void G4PSSphereSurfaceFlux::EndOfEvent(G4HCofThisEvent*)
190{;}
191
192void G4PSSphereSurfaceFlux::clear(){
193 EvtMap->clear();
194}
195
196void G4PSSphereSurfaceFlux::DrawAll()
197{;}
198
199void G4PSSphereSurfaceFlux::PrintAll()
200{
201 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
202 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
203 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
204 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
205 for(; itr != EvtMap->GetMap()->end(); itr++) {
206 G4cout << " copy no.: " << itr->first
207 << " current : " << *(itr->second)
208 << G4endl;
209 }
210}
211
Note: See TracBrowser for help on using the repository browser.