source: trunk/source/digits_hits/scorer/src/G4PSSphereSurfaceCurrent.cc@ 952

Last change on this file since 952 was 850, checked in by garnier, 17 years ago

geant4.8.2 beta

File size: 6.3 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: G4PSSphereSurfaceCurrent.cc,v 1.2 2007/08/14 21:23:52 taso Exp $
28// GEANT4 tag $Name: HEAD $
29//
30// G4PSSphereSurfaceCurrent
31#include "G4PSSphereSurfaceCurrent.hh"
32#include "G4StepStatus.hh"
33#include "G4Track.hh"
34#include "G4UnitsTable.hh"
35#include "G4GeometryTolerance.hh"
36////////////////////////////////////////////////////////////////////////////////
37// (Description)
38// This is a primitive scorer class for scoring only Surface Current.
39// Current version assumes only for G4Sphere shape.
40//
41// Surface is defined at the inside of sphere.
42// Direction -Rmin +Rmax
43// 0 IN || OUT ->|<- |
44// 1 IN ->| |
45// 2 OUT |<- |
46//
47// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
48//
49///////////////////////////////////////////////////////////////////////////////
50
51G4PSSphereSurfaceCurrent::G4PSSphereSurfaceCurrent(G4String name,
52 G4int direction, G4int depth)
53 :G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),
54 weighted(true),divideByArea(true)
55{;}
56
57G4PSSphereSurfaceCurrent::~G4PSSphereSurfaceCurrent()
58{;}
59
60G4bool G4PSSphereSurfaceCurrent::ProcessHits(G4Step* aStep,G4TouchableHistory*)
61{
62 G4StepPoint* preStep = aStep->GetPreStepPoint();
63 G4VSolid * solid =
64 preStep->GetPhysicalVolume()->GetLogicalVolume()->GetSolid();
65 if( solid->GetEntityType() != "G4Sphere" ){
66 G4Exception("G4PSSphereSurfaceCurrentScorer. - Solid type is not supported.");
67 return FALSE;
68 }
69 G4Sphere* sphereSolid = (G4Sphere*)(solid);
70
71 G4int dirFlag =IsSelectedSurface(aStep,sphereSolid);
72 if ( dirFlag > 0 ) {
73 if ( fDirection == fCurrent_InOut || fDirection == dirFlag ){
74 G4double radi = sphereSolid->GetInsideRadius();
75 G4double dph = sphereSolid->GetDeltaPhiAngle()/radian;
76 G4double stth = sphereSolid->GetStartThetaAngle()/radian;
77 G4double enth = stth+sphereSolid->GetDeltaThetaAngle()/radian;
78 G4double current = 1.0;
79 if ( weighted) current = preStep->GetWeight(); // Current (Particle Weight)
80 if ( divideByArea ){
81 G4double square = radi*radi*dph*( -std::cos(enth) + std::cos(stth) );
82 current = current/square; // Current with angle.
83 }
84
85 G4int index = GetIndex(aStep);
86 EvtMap->add(index,current);
87 }
88 }
89
90 return TRUE;
91}
92
93G4int G4PSSphereSurfaceCurrent::IsSelectedSurface(G4Step* aStep, G4Sphere* sphereSolid){
94
95 G4TouchableHandle theTouchable =
96 aStep->GetPreStepPoint()->GetTouchableHandle();
97 G4double kCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
98
99 if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
100 // Entering Geometry
101 G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
102 G4ThreeVector localpos1 =
103 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
104 G4double localR2 = localpos1.x()*localpos1.x()
105 +localpos1.y()*localpos1.y()
106 +localpos1.z()*localpos1.z();
107 G4double InsideRadius2 =
108 sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
109 if(std::fabs( localR2 - InsideRadius2 ) < kCarTolerance ){
110 return fCurrent_In;
111 }
112 }
113
114 if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
115 // Exiting Geometry
116 G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
117 G4ThreeVector localpos2 =
118 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
119 G4double localR2 = localpos2.x()*localpos2.x()
120 +localpos2.y()*localpos2.y()
121 +localpos2.z()*localpos2.z();
122 G4double InsideRadius2 =
123 sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
124 if(std::fabs( localR2 - InsideRadius2 ) < kCarTolerance ){
125 return fCurrent_Out;
126 }
127 }
128
129 return -1;
130}
131
132void G4PSSphereSurfaceCurrent::Initialize(G4HCofThisEvent* HCE)
133{
134 EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName());
135 if ( HCID < 0 ) HCID = GetCollectionID(0);
136 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
137}
138
139void G4PSSphereSurfaceCurrent::EndOfEvent(G4HCofThisEvent*)
140{;}
141
142void G4PSSphereSurfaceCurrent::clear(){
143 EvtMap->clear();
144}
145
146void G4PSSphereSurfaceCurrent::DrawAll()
147{;}
148
149void G4PSSphereSurfaceCurrent::PrintAll()
150{
151 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
152 G4cout << " PrimitiveScorer " << GetName() <<G4endl;
153 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
154 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
155 for(; itr != EvtMap->GetMap()->end(); itr++) {
156 G4cout << " copy no.: " << itr->first << " current : " ;
157 if ( divideByArea ) G4cout << *(itr->second)*cm*cm << " [/cm2]" ;
158 else G4cout << *(itr->second)*cm*cm << " [track]" ;
159 G4cout << G4endl;
160 }
161}
162
Note: See TracBrowser for help on using the repository browser.