00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030 #include "G4PSFlatSurfaceFlux.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
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00059
00060 G4PSFlatSurfaceFlux::G4PSFlatSurfaceFlux(G4String name,
00061 G4int direction, G4int depth)
00062 : G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),
00063 weighted(true),divideByArea(true)
00064 {
00065 DefineUnitAndCategory();
00066 SetUnit("percm2");
00067 }
00068
00069 G4PSFlatSurfaceFlux::G4PSFlatSurfaceFlux(G4String name,
00070 G4int direction,
00071 const G4String& unit,
00072 G4int depth)
00073 : G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),
00074 weighted(true),divideByArea(true)
00075 {
00076 DefineUnitAndCategory();
00077 SetUnit(unit);
00078 }
00079
00080 G4PSFlatSurfaceFlux::~G4PSFlatSurfaceFlux()
00081 {;}
00082
00083 G4bool G4PSFlatSurfaceFlux::ProcessHits(G4Step* aStep,G4TouchableHistory*)
00084 {
00085 G4StepPoint* preStep = aStep->GetPreStepPoint();
00086 G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
00087 G4VPVParameterisation* physParam = physVol->GetParameterisation();
00088 G4VSolid * solid = 0;
00089 if(physParam)
00090 {
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 {
00098 solid = physVol->GetLogicalVolume()->GetSolid();
00099 }
00100
00101 G4Box* boxSolid = (G4Box*)(solid);
00102
00103 G4int dirFlag =IsSelectedSurface(aStep,boxSolid);
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
00121 G4double angleFactor = localdir.z();
00122 if ( angleFactor < 0 ) angleFactor *= -1.;
00123 G4double flux = 1.0;
00124 if ( weighted ) flux *=preStep->GetWeight();
00125
00126 G4double square = 4.*boxSolid->GetXHalfLength()*boxSolid->GetYHalfLength();
00127
00128 flux = flux/angleFactor;
00129 if ( divideByArea ) flux /= square;
00130
00131 G4int index = GetIndex(aStep);
00132 EvtMap->add(index,flux);
00133 }
00134 }
00135 #ifdef debug
00136 G4cout << " PASSED vol "
00137 << index << " trk "<<trkid<<" len " << fFlatSurfaceFlux<<G4endl;
00138 #endif
00139
00140 return TRUE;
00141 }
00142
00143 G4int G4PSFlatSurfaceFlux::IsSelectedSurface(G4Step* aStep, G4Box* boxSolid){
00144
00145 G4TouchableHandle theTouchable =
00146 aStep->GetPreStepPoint()->GetTouchableHandle();
00147 G4double kCarTolerance=G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
00148
00149 if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
00150
00151 G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
00152 G4ThreeVector localpos1 =
00153 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
00154 if(std::fabs( localpos1.z() + boxSolid->GetZHalfLength())<kCarTolerance ){
00155 return fFlux_In;
00156 }
00157 }
00158
00159 if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
00160
00161 G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
00162 G4ThreeVector localpos2 =
00163 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
00164 if(std::fabs( localpos2.z() + boxSolid->GetZHalfLength())<kCarTolerance ){
00165 return fFlux_Out;
00166 }
00167 }
00168
00169 return -1;
00170 }
00171
00172 void G4PSFlatSurfaceFlux::Initialize(G4HCofThisEvent* HCE)
00173 {
00174 EvtMap = new G4THitsMap<G4double>(GetMultiFunctionalDetector()->GetName(),
00175 GetName());
00176 if ( HCID < 0 ) HCID = GetCollectionID(0);
00177 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
00178 }
00179
00180 void G4PSFlatSurfaceFlux::EndOfEvent(G4HCofThisEvent*)
00181 {;}
00182
00183 void G4PSFlatSurfaceFlux::clear(){
00184 EvtMap->clear();
00185 }
00186
00187 void G4PSFlatSurfaceFlux::DrawAll()
00188 {;}
00189
00190 void G4PSFlatSurfaceFlux::PrintAll()
00191 {
00192 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
00193 G4cout << " PrimitiveScorer" << GetName() <<G4endl;
00194 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
00195 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
00196 for(; itr != EvtMap->GetMap()->end(); itr++) {
00197 G4cout << " copy no.: " << itr->first
00198 << " flux : " << *(itr->second)/GetUnitValue()
00199 << " [" << GetUnit() <<"]"
00200 << G4endl;
00201 }
00202 }
00203
00204 void G4PSFlatSurfaceFlux::SetUnit(const G4String& unit)
00205 {
00206 if ( divideByArea ) {
00207 CheckAndSetUnit(unit,"Per Unit Surface");
00208 } else {
00209 if (unit == "" ){
00210 unitName = unit;
00211 unitValue = 1.0;
00212 }else{
00213 G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] ) for " + GetName();
00214 G4Exception("G4PSFlatSurfaceFlux::SetUnit","DetPS0008",JustWarning,msg);
00215 }
00216 }
00217 }
00218
00219 void G4PSFlatSurfaceFlux::DefineUnitAndCategory(){
00220
00221 new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
00222 new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
00223 new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
00224 }
00225
00226