G4PSSphereSurfaceFlux.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id$
00028 //
00029 // G4PSSphereSurfaceFlux
00030 #include "G4PSSphereSurfaceFlux.hh"
00031 
00032 #include "G4SystemOfUnits.hh"
00033 #include "G4StepStatus.hh"
00034 #include "G4Track.hh"
00035 #include "G4VSolid.hh"
00036 #include "G4VPhysicalVolume.hh"
00037 #include "G4VPVParameterisation.hh"
00038 #include "G4UnitsTable.hh"
00039 #include "G4GeometryTolerance.hh"
00041 // (Description)
00042 //   This is a primitive scorer class for scoring only Surface Flux.
00043 //  Flux version assumes only for G4Sphere shape. 
00044 //
00045 // Surface is defined  at the inside of sphere.
00046 // Direction                  -Rmin   +Rmax
00047 //   0  IN || OUT            ->|<-     |
00048 //   1  IN                   ->|       |
00049 //   2  OUT                    |<-     |
00050 //
00051 // Created: 2005-11-14  Tsukasa ASO, Akinori Kimura.
00052 // 29-Mar-2007  T.Aso,  Bug fix for momentum direction at outgoing flux.
00053 // 2010-07-22   Introduce Unit specification.
00054 // 2010-07-22   Add weighted and divideByAre options
00055 // 2011-02-21   Get correct momentum direction in Flux_Out. 
00056 // 2011-09-09   Modify comment in PrintAll().
00058 
00059 G4PSSphereSurfaceFlux::G4PSSphereSurfaceFlux(G4String name, 
00060                                          G4int direction, G4int depth)
00061   : G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),
00062     weighted(true),divideByArea(true)
00063 {
00064     DefineUnitAndCategory();
00065     SetUnit("percm2");
00066 }
00067 
00068 G4PSSphereSurfaceFlux::G4PSSphereSurfaceFlux(G4String name, 
00069                                              G4int direction,
00070                                              const G4String& unit,
00071                                              G4int depth)
00072   : G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),
00073     weighted(true),divideByArea(true)
00074 {
00075     DefineUnitAndCategory();
00076     SetUnit(unit);
00077 }
00078 
00079 G4PSSphereSurfaceFlux::~G4PSSphereSurfaceFlux()
00080 {;}
00081 
00082 G4bool G4PSSphereSurfaceFlux::ProcessHits(G4Step* aStep,G4TouchableHistory*)
00083 {
00084   G4StepPoint* preStep = aStep->GetPreStepPoint();
00085 
00086   G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
00087   G4VPVParameterisation* physParam = physVol->GetParameterisation();
00088   G4VSolid * solid = 0;
00089   if(physParam)
00090   { // for parameterized volume
00091     G4int idx = ((G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable()))
00092                 ->GetReplicaNumber(indexDepth);
00093     solid = physParam->ComputeSolid(idx, physVol);
00094     solid->ComputeDimensions(physParam,idx,physVol);
00095   }
00096   else
00097   { // for ordinary volume
00098     solid = physVol->GetLogicalVolume()->GetSolid();
00099   }
00100 
00101   G4Sphere* sphereSolid = (G4Sphere*)(solid);
00102 
00103   G4int dirFlag =IsSelectedSurface(aStep,sphereSolid);
00104   if ( dirFlag > 0 ) {
00105     if ( fDirection == fFlux_InOut || fDirection == dirFlag ){
00106 
00107       G4StepPoint* thisStep=0;
00108       if ( dirFlag == fFlux_In ){
00109         thisStep = preStep;
00110       }else if ( dirFlag == fFlux_Out ){
00111         thisStep = aStep->GetPostStepPoint();
00112       }else{
00113         return FALSE;
00114       }
00115 
00116       G4TouchableHandle theTouchable = thisStep->GetTouchableHandle();
00117       G4ThreeVector pdirection = thisStep->GetMomentumDirection();
00118       G4ThreeVector localdir  = 
00119         theTouchable->GetHistory()->GetTopTransform().TransformAxis(pdirection);
00120       G4double localdirL2 = localdir.x()*localdir.x()
00121         +localdir.y()*localdir.y()
00122         +localdir.z()*localdir.z();
00123       G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
00124       G4ThreeVector localpos1 = 
00125         theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
00126       G4double localR2 = localpos1.x()*localpos1.x()
00127         +localpos1.y()*localpos1.y()
00128         +localpos1.z()*localpos1.z();
00129       G4double anglefactor = (localdir.x()*localpos1.x()
00130                               +localdir.y()*localpos1.y()
00131                               +localdir.z()*localpos1.z())
00132         /std::sqrt(localdirL2)/std::sqrt(localR2);
00133 
00134       G4double radi   = sphereSolid->GetInsideRadius();
00135       G4double dph    = sphereSolid->GetDeltaPhiAngle()/radian;
00136       G4double stth   = sphereSolid->GetStartThetaAngle()/radian;
00137       G4double enth   = stth+sphereSolid->GetDeltaThetaAngle()/radian;
00138       G4double square = radi*radi*dph*( -std::cos(enth) + std::cos(stth) );
00139 
00140       G4double current = 1.0;
00141       if ( weighted ) thisStep->GetWeight(); // Flux (Particle Weight)
00142       if ( divideByArea ) current = current/square;  // Flux with angle.
00143 
00144       current /= anglefactor;
00145 
00146       G4int index = GetIndex(aStep);
00147       EvtMap->add(index,current);
00148     }
00149   }
00150 
00151   return TRUE;
00152 }
00153 
00154 G4int G4PSSphereSurfaceFlux::IsSelectedSurface(G4Step* aStep, G4Sphere* sphereSolid){
00155 
00156   G4TouchableHandle theTouchable = 
00157     aStep->GetPreStepPoint()->GetTouchableHandle();
00158   G4double kCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
00159   
00160   if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
00161     // Entering Geometry
00162     G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
00163     G4ThreeVector localpos1 = 
00164       theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
00165     G4double localR2 = localpos1.x()*localpos1.x()
00166                       +localpos1.y()*localpos1.y()
00167                       +localpos1.z()*localpos1.z();
00168     //G4double InsideRadius2 = 
00169     //  sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
00170     //if(std::fabs( localR2 - InsideRadius2 ) < kCarTolerance ){
00171     G4double InsideRadius = sphereSolid->GetInsideRadius();
00172     if ( localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
00173          &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
00174       return fFlux_In;
00175     }
00176   }
00177 
00178   if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
00179     // Exiting Geometry
00180     G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
00181     G4ThreeVector localpos2 = 
00182       theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
00183     G4double localR2 = localpos2.x()*localpos2.x()
00184                       +localpos2.y()*localpos2.y()
00185                       +localpos2.z()*localpos2.z();
00186     //G4double InsideRadius2 = 
00187     //  sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
00188     //if(std::facb(localR2 - InsideRadius2) ) < kCarTolerance ){
00189     G4double InsideRadius = sphereSolid->GetInsideRadius();
00190     if ( localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
00191          &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
00192       return fFlux_Out;
00193     }
00194   }
00195 
00196   return -1;
00197 }
00198 
00199 void G4PSSphereSurfaceFlux::Initialize(G4HCofThisEvent* HCE)
00200 {
00201   EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName());
00202   if ( HCID < 0 ) HCID = GetCollectionID(0);
00203   HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
00204 }
00205 
00206 void G4PSSphereSurfaceFlux::EndOfEvent(G4HCofThisEvent*)
00207 {;}
00208 
00209 void G4PSSphereSurfaceFlux::clear(){
00210   EvtMap->clear();
00211 }
00212 
00213 void G4PSSphereSurfaceFlux::DrawAll()
00214 {;}
00215 
00216 void G4PSSphereSurfaceFlux::PrintAll()
00217 {
00218   G4cout << " MultiFunctionalDet  " << detector->GetName() << G4endl;
00219   G4cout << " PrimitiveScorer " << GetName() <<G4endl; 
00220   G4cout << " Number of entries " << EvtMap->entries() << G4endl;
00221   std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
00222   for(; itr != EvtMap->GetMap()->end(); itr++) {
00223     G4cout << "  copy no.: " << itr->first
00224            << "  Flux  : " << *(itr->second)/GetUnitValue()
00225            << " ["<<GetUnit()<<"]"
00226            << G4endl;
00227   }
00228 }
00229 
00230 void G4PSSphereSurfaceFlux::SetUnit(const G4String& unit)
00231 {
00232     if ( divideByArea ) {
00233         CheckAndSetUnit(unit,"Per Unit Surface");
00234     } else {
00235         if (unit == "" ){
00236             unitName = unit;
00237             unitValue = 1.0;
00238         }else{
00239             G4String msg = "Invalid unit ["+unit+"] (Current  unit is [" +GetUnit()+"] ) for " + GetName();
00240             G4Exception("G4PSSphereSurfaceFlux::SetUnit","DetPS0016",JustWarning,msg);
00241         }
00242     }
00243 }
00244 
00245 void G4PSSphereSurfaceFlux::DefineUnitAndCategory(){
00246    // Per Unit Surface
00247    new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
00248    new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
00249    new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
00250 }
00251 

Generated on Mon May 27 17:49:29 2013 for Geant4 by  doxygen 1.4.7