G4VScoringMesh Class Reference

#include <G4VScoringMesh.hh>

Inheritance diagram for G4VScoringMesh:

G4ScoringBox G4ScoringCylinder

Public Member Functions

 G4VScoringMesh (const G4String &wName)
virtual ~G4VScoringMesh ()
virtual void Construct (G4VPhysicalVolume *fWorldPhys)=0
virtual void List () const
const G4StringGetWorldName () const
G4bool IsActive () const
void Activate (G4bool vl=true)
MeshShape GetShape () const
void Accumulate (G4THitsMap< G4double > *map)
void Dump ()
void DrawMesh (const G4String &psName, G4VScoreColorMap *colorMap, G4int axflg=111)
void DrawMesh (const G4String &psName, G4int idxPlane, G4int iColumn, G4VScoreColorMap *colorMap)
virtual void Draw (std::map< G4int, G4double * > *map, G4VScoreColorMap *colorMap, G4int axflg=111)=0
virtual void DrawColumn (std::map< G4int, G4double * > *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn)=0
void ResetScore ()
void SetSize (G4double size[3])
G4ThreeVector GetSize () const
void SetCenterPosition (G4double centerPosition[3])
G4ThreeVector GetTranslation () const
void RotateX (G4double delta)
void RotateY (G4double delta)
void RotateZ (G4double delta)
G4RotationMatrix GetRotationMatrix () const
void SetNumberOfSegments (G4int nSegment[3])
void GetNumberOfSegments (G4int nSegment[3])
void SetPrimitiveScorer (G4VPrimitiveScorer *ps)
void SetFilter (G4VSDFilter *filter)
void SetCurrentPrimitiveScorer (const G4String &name)
G4bool FindPrimitiveScorer (const G4String &psname)
G4bool IsCurrentPrimitiveScorerNull ()
G4String GetPSUnit (const G4String &psname)
G4String GetCurrentPSUnit ()
void SetCurrentPSUnit (const G4String &unit)
G4double GetPSUnitValue (const G4String &psname)
void SetDrawPSName (const G4String &psname)
void GetDivisionAxisNames (G4String divisionAxisNames[3])
void SetNullToCurrentPrimitiveScorer ()
void SetVerboseLevel (G4int vl)
MeshScoreMap GetScoreMap ()
G4bool ReadyForQuantity () const

Protected Member Functions

G4VPrimitiveScorerGetPrimitiveScorer (const G4String &name)

Protected Attributes

G4String fWorldName
G4VPrimitiveScorerfCurrentPS
G4bool fConstructed
G4bool fActive
MeshShape fShape
G4double fSize [3]
G4ThreeVector fCenterPosition
G4RotationMatrixfRotationMatrix
G4int fNSegment [3]
std::map< G4String, G4THitsMap<
G4double > * > 
fMap
G4MultiFunctionalDetectorfMFD
G4int verboseLevel
G4bool sizeIsSet
G4bool nMeshIsSet
G4String fDrawUnit
G4double fDrawUnitValue
G4String fDrawPSName
G4String fDivisionAxisNames [3]

Detailed Description

Definition at line 52 of file G4VScoringMesh.hh.


Constructor & Destructor Documentation

G4VScoringMesh::G4VScoringMesh ( const G4String wName  ) 

Definition at line 46 of file G4VScoringMesh.cc.

References G4SDManager::AddNewDetector(), fDivisionAxisNames, fMFD, fNSegment, fSize, and G4SDManager::GetSDMpointer().

00047   : fWorldName(wName),fCurrentPS(0),fConstructed(false),fActive(true),
00048     fRotationMatrix(0), fMFD(new G4MultiFunctionalDetector(wName)),
00049     verboseLevel(0),sizeIsSet(false),nMeshIsSet(false),
00050     fDrawUnit(""), fDrawUnitValue(1.)
00051 {
00052   G4SDManager::GetSDMpointer()->AddNewDetector(fMFD);
00053 
00054   fSize[0] = fSize[1] = fSize[2] = 0.;
00055   fNSegment[0] = fNSegment[1] = fNSegment[2] = 1;
00056   fDivisionAxisNames[0] = fDivisionAxisNames[1] = fDivisionAxisNames[2] = "";
00057 }

G4VScoringMesh::~G4VScoringMesh (  )  [virtual]

Definition at line 59 of file G4VScoringMesh.cc.

00060 {
00061   ;
00062 }


Member Function Documentation

void G4VScoringMesh::Accumulate ( G4THitsMap< G4double > *  map  )  [inline]

Definition at line 186 of file G4VScoringMesh.hh.

References fMap, G4cout, G4endl, G4VHitsCollection::GetName(), G4THitsMap< T >::GetSize(), G4THitsMap< T >::PrintAllHits(), and verboseLevel.

00187 {
00188   G4String psName = map->GetName();
00189   std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
00190   *(fMapItr->second) += *map;
00191 
00192   if(verboseLevel > 9) {
00193     G4cout << G4endl;
00194     G4cout << "G4VScoringMesh::Accumulate()" << G4endl;
00195     G4cout << "  PS name : " << psName << G4endl;
00196     if(fMapItr == fMap.end()) {
00197       G4cout << "  "
00198              << psName << " was not found." << G4endl;
00199     } else {
00200       G4cout << "  map size : " << map->GetSize() << G4endl;
00201       map->PrintAllHits();
00202     }
00203     G4cout << G4endl;
00204   }
00205 }

void G4VScoringMesh::Activate ( G4bool  vl = true  )  [inline]

Definition at line 72 of file G4VScoringMesh.hh.

References fActive.

00073   { fActive = vl; }

virtual void G4VScoringMesh::Construct ( G4VPhysicalVolume fWorldPhys  )  [pure virtual]

Implemented in G4ScoringBox, and G4ScoringCylinder.

virtual void G4VScoringMesh::Draw ( std::map< G4int, G4double * > *  map,
G4VScoreColorMap colorMap,
G4int  axflg = 111 
) [pure virtual]

Implemented in G4ScoringBox, and G4ScoringCylinder.

Referenced by DrawMesh().

virtual void G4VScoringMesh::DrawColumn ( std::map< G4int, G4double * > *  map,
G4VScoreColorMap colorMap,
G4int  idxProj,
G4int  idxColumn 
) [pure virtual]

Implemented in G4ScoringBox, and G4ScoringCylinder.

Referenced by DrawMesh().

void G4VScoringMesh::DrawMesh ( const G4String psName,
G4int  idxPlane,
G4int  iColumn,
G4VScoreColorMap colorMap 
)

Definition at line 301 of file G4VScoringMesh.cc.

References DrawColumn(), fDrawPSName, fDrawUnit, fDrawUnitValue, fMap, G4cerr, G4endl, GetPSUnit(), and GetPSUnitValue().

00302 {
00303   fDrawPSName = psName;
00304   std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
00305   if(fMapItr!=fMap.end()) {
00306     fDrawUnit = GetPSUnit(psName);
00307     fDrawUnitValue = GetPSUnitValue(psName);
00308     DrawColumn(fMapItr->second->GetMap(),colorMap,idxPlane,iColumn);
00309   } else {
00310     G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl;
00311   }
00312 }

void G4VScoringMesh::DrawMesh ( const G4String psName,
G4VScoreColorMap colorMap,
G4int  axflg = 111 
)

Definition at line 288 of file G4VScoringMesh.cc.

References Draw(), fDrawPSName, fDrawUnit, fDrawUnitValue, fMap, G4cerr, G4endl, GetPSUnit(), and GetPSUnitValue().

Referenced by G4ScoringManager::DrawMesh().

00289 {
00290   fDrawPSName = psName;
00291   std::map<G4String, G4THitsMap<G4double>* >::const_iterator fMapItr = fMap.find(psName);
00292   if(fMapItr!=fMap.end()) {
00293     fDrawUnit = GetPSUnit(psName);
00294     fDrawUnitValue = GetPSUnitValue(psName);
00295     Draw(fMapItr->second->GetMap(), colorMap,axflg);
00296   } else {
00297     G4cerr << "Scorer <" << psName << "> is not defined. Method ignored." << G4endl;
00298   }
00299 }

void G4VScoringMesh::Dump (  ) 

Definition at line 274 of file G4VScoringMesh.cc.

References fMap, fWorldName, G4cout, and G4endl.

00274                           {
00275   G4cout << "scoring mesh name: " << fWorldName << G4endl;
00276 
00277   G4cout << "# of G4THitsMap : " << fMap.size() << G4endl;
00278   std::map<G4String, G4THitsMap<G4double>* >::iterator itr = fMap.begin();
00279   for(; itr != fMap.end(); itr++) {
00280     G4cout << "[" << itr->first << "]" << G4endl;
00281     itr->second->PrintAllHits();
00282   }
00283   G4cout << G4endl;
00284 
00285 }

G4bool G4VScoringMesh::FindPrimitiveScorer ( const G4String psname  ) 

Definition at line 172 of file G4VScoringMesh.cc.

References fMap.

Referenced by G4ScoreQuantityMessenger::CheckMeshPS().

00172                                                                   {
00173   std::map<G4String, G4THitsMap<G4double>* >::iterator itr = fMap.find(psname);;
00174   if(itr == fMap.end()) return false;
00175   return true;
00176 }

G4String G4VScoringMesh::GetCurrentPSUnit (  ) 

Definition at line 187 of file G4VScoringMesh.cc.

References fCurrentPS, G4cerr, G4endl, and G4VPrimitiveScorer::GetUnit().

Referenced by G4ScoreQuantityMessenger::SetNewValue().

00187                                          {
00188     G4String unit = "";
00189   if(fCurrentPS == 0) {
00190       G4String msg = "ERROR : G4VScoringMesh::GetCurrentPSUnit() : ";
00191       msg += " Current primitive scorer is null.";
00192       G4cerr << msg << G4endl;
00193   }else{
00194      unit =  fCurrentPS->GetUnit();
00195   }
00196   return unit;
00197 }

void G4VScoringMesh::GetDivisionAxisNames ( G4String  divisionAxisNames[3]  ) 

Definition at line 218 of file G4VScoringMesh.cc.

References fDivisionAxisNames.

Referenced by G4VScoreWriter::DumpAllQuantitiesToFile(), and G4VScoreWriter::DumpQuantityToFile().

00218                                                                        {
00219   for(int i = 0; i < 3; i++) divisionAxisNames[i] = fDivisionAxisNames[i];
00220 }

void G4VScoringMesh::GetNumberOfSegments ( G4int  nSegment[3]  ) 

Definition at line 103 of file G4VScoringMesh.cc.

References fNSegment.

Referenced by G4GMocrenFileSceneHandler::AddSolid(), and G4VScoreWriter::SetScoringMesh().

00103                                                           {
00104   for(int i = 0; i < 3; i++) nSegment[i] = fNSegment[i];
00105 }

G4VPrimitiveScorer * G4VScoringMesh::GetPrimitiveScorer ( const G4String name  )  [protected]

Definition at line 222 of file G4VScoringMesh.cc.

References fMFD, G4MultiFunctionalDetector::GetNumberOfPrimitives(), and G4MultiFunctionalDetector::GetPrimitive().

Referenced by GetPSUnit(), GetPSUnitValue(), and SetCurrentPrimitiveScorer().

00222                                                                              {
00223   if(fMFD == 0) return 0;
00224 
00225   G4int nps = fMFD->GetNumberOfPrimitives();
00226   for(G4int i = 0; i < nps; i++) {
00227     G4VPrimitiveScorer * ps = fMFD->GetPrimitive(i);
00228     if(name == ps->GetName()) return ps;
00229   }
00230 
00231   return 0;
00232 }

G4String G4VScoringMesh::GetPSUnit ( const G4String psname  ) 

Definition at line 178 of file G4VScoringMesh.cc.

References fMap, GetPrimitiveScorer(), and G4VPrimitiveScorer::GetUnit().

Referenced by DrawMesh(), G4VScoreWriter::DumpAllQuantitiesToFile(), and G4VScoreWriter::DumpQuantityToFile().

00178                                                           {
00179   std::map<G4String, G4THitsMap<G4double>* >::iterator itr = fMap.find(psname);;
00180   if(itr == fMap.end()) {
00181     return G4String("");
00182   } else {
00183     return GetPrimitiveScorer(psname)->GetUnit();
00184   }
00185 }

G4double G4VScoringMesh::GetPSUnitValue ( const G4String psname  ) 

Definition at line 209 of file G4VScoringMesh.cc.

References fMap, GetPrimitiveScorer(), and G4VPrimitiveScorer::GetUnitValue().

Referenced by DrawMesh(), G4VScoreWriter::DumpAllQuantitiesToFile(), and G4VScoreWriter::DumpQuantityToFile().

00209                                                                {
00210   std::map<G4String, G4THitsMap<G4double>* >::iterator itr = fMap.find(psname);;
00211   if(itr == fMap.end()) {
00212     return 1.;
00213   } else {
00214     return GetPrimitiveScorer(psname)->GetUnitValue();
00215   }
00216 }

G4RotationMatrix G4VScoringMesh::GetRotationMatrix (  )  const [inline]

Definition at line 108 of file G4VScoringMesh.hh.

References fRotationMatrix.

Referenced by G4GMocrenFileSceneHandler::AddSolid().

00108                                              {
00109     if(fRotationMatrix) return *fRotationMatrix;
00110     else return G4RotationMatrix::IDENTITY;
00111   }

MeshScoreMap G4VScoringMesh::GetScoreMap (  )  [inline]

Definition at line 150 of file G4VScoringMesh.hh.

References fMap.

Referenced by G4VScoreWriter::DumpAllQuantitiesToFile(), and G4VScoreWriter::DumpQuantityToFile().

00150 {return fMap;}

MeshShape G4VScoringMesh::GetShape (  )  const [inline]

Definition at line 75 of file G4VScoringMesh.hh.

References fShape.

Referenced by G4ScoringMessenger::SetNewValue().

00076   { return fShape; }

G4ThreeVector G4VScoringMesh::GetSize (  )  const

Definition at line 83 of file G4VScoringMesh.cc.

References fSize, and sizeIsSet.

Referenced by G4GMocrenFileSceneHandler::AddSolid().

00083                                             {
00084   if(sizeIsSet)
00085     return G4ThreeVector(fSize[0], fSize[1], fSize[2]);
00086   else
00087     return G4ThreeVector(0., 0., 0.);
00088 }

G4ThreeVector G4VScoringMesh::GetTranslation (  )  const [inline]

Definition at line 100 of file G4VScoringMesh.hh.

References fCenterPosition.

Referenced by G4GMocrenFileSceneHandler::AddSolid().

00100 {return fCenterPosition;}

const G4String& G4VScoringMesh::GetWorldName (  )  const [inline]

Definition at line 66 of file G4VScoringMesh.hh.

References fWorldName.

Referenced by G4VScoreWriter::DumpAllQuantitiesToFile(), G4VScoreWriter::DumpQuantityToFile(), and G4ScoringMessenger::SetNewValue().

00067   { return fWorldName; }

G4bool G4VScoringMesh::IsActive (  )  const [inline]

Definition at line 69 of file G4VScoringMesh.hh.

References fActive.

00070   { return fActive; }

G4bool G4VScoringMesh::IsCurrentPrimitiveScorerNull (  )  [inline]

Definition at line 126 of file G4VScoringMesh.hh.

References fCurrentPS.

Referenced by G4ScoreQuantityMessenger::SetNewValue().

00126                                         {
00127     if(fCurrentPS == NULL) return true;
00128     else return false;
00129   }

void G4VScoringMesh::List (  )  const [virtual]

Reimplemented in G4ScoringBox, and G4ScoringCylinder.

Definition at line 233 of file G4VScoringMesh.cc.

References fCenterPosition, fMFD, fNSegment, fRotationMatrix, G4cout, G4endl, G4VPrimitiveScorer::GetFilter(), G4VSDFilter::GetName(), G4VPrimitiveScorer::GetName(), G4MultiFunctionalDetector::GetNumberOfPrimitives(), and G4MultiFunctionalDetector::GetPrimitive().

Referenced by G4ScoringCylinder::List(), and G4ScoringBox::List().

00233                                 {
00234 
00235   G4cout << " # of segments: ("
00236          << fNSegment[0] << ", "
00237          << fNSegment[1] << ", "
00238          << fNSegment[2] << ")"
00239          << G4endl;
00240   G4cout << " displacement: ("
00241          << fCenterPosition.x()/cm << ", "
00242          << fCenterPosition.y()/cm << ", "
00243          << fCenterPosition.z()/cm << ") [cm]"
00244          << G4endl;
00245   if(fRotationMatrix != 0) {
00246     G4cout << " rotation matrix: "
00247            << fRotationMatrix->xx() << "  "
00248            << fRotationMatrix->xy() << "  "
00249            << fRotationMatrix->xz() << G4endl
00250            << "                  "
00251            << fRotationMatrix->yx() << "  "
00252            << fRotationMatrix->yy() << "  "
00253            << fRotationMatrix->yz() << G4endl
00254            << "                  "
00255            << fRotationMatrix->zx() << "  "
00256            << fRotationMatrix->zy() << "  "
00257            << fRotationMatrix->zz() << G4endl;
00258   }
00259 
00260 
00261   G4cout << " registered primitve scorers : " << G4endl;
00262   G4int nps = fMFD->GetNumberOfPrimitives();
00263   G4VPrimitiveScorer * ps;
00264   for(int i = 0; i < nps; i++) {
00265     ps = fMFD->GetPrimitive(i);
00266     G4cout << "   " << i << "  " << ps->GetName();
00267     if(ps->GetFilter() != 0) G4cout << "     with  " << ps->GetFilter()->GetName();
00268     G4cout << G4endl;
00269   }
00270 
00271 
00272 }

G4bool G4VScoringMesh::ReadyForQuantity (  )  const [inline]

Definition at line 152 of file G4VScoringMesh.hh.

References nMeshIsSet, and sizeIsSet.

Referenced by SetPrimitiveScorer().

00153   { return (sizeIsSet && nMeshIsSet); }

void G4VScoringMesh::ResetScore (  ) 

Definition at line 64 of file G4VScoringMesh.cc.

References fMap, G4cout, G4endl, and verboseLevel.

Referenced by G4ScoringCylinder::Construct(), and G4ScoringBox::Construct().

00064                                 {
00065   if(verboseLevel > 9) G4cout << "G4VScoringMesh::ResetScore() is called." << G4endl;
00066   std::map<G4String, G4THitsMap<G4double>* >::iterator itr = fMap.begin();
00067   for(; itr != fMap.end(); itr++) {
00068     if(verboseLevel > 9) G4cout << "G4VScoringMesh::ResetScore()" << itr->first << G4endl;
00069     itr->second->clear();
00070   }
00071 }

void G4VScoringMesh::RotateX ( G4double  delta  ) 

Definition at line 106 of file G4VScoringMesh.cc.

References fRotationMatrix.

Referenced by G4ScoringMessenger::SetNewValue().

00106                                            {
00107   if(fRotationMatrix == 0) fRotationMatrix = new G4RotationMatrix();
00108   fRotationMatrix->rotateX(delta);
00109 }

void G4VScoringMesh::RotateY ( G4double  delta  ) 

Definition at line 111 of file G4VScoringMesh.cc.

References fRotationMatrix.

Referenced by G4ScoringMessenger::SetNewValue().

00111                                            {
00112   if(fRotationMatrix == 0) fRotationMatrix = new G4RotationMatrix();
00113   fRotationMatrix->rotateY(delta);
00114 }

void G4VScoringMesh::RotateZ ( G4double  delta  ) 

Definition at line 116 of file G4VScoringMesh.cc.

References fRotationMatrix.

Referenced by G4ScoringMessenger::SetNewValue().

00116                                            {
00117   if(fRotationMatrix == 0) fRotationMatrix = new G4RotationMatrix();
00118   fRotationMatrix->rotateZ(delta);
00119 }

void G4VScoringMesh::SetCenterPosition ( G4double  centerPosition[3]  ) 

Definition at line 89 of file G4VScoringMesh.cc.

References fCenterPosition.

Referenced by G4ScoringMessenger::SetNewValue().

00089                                                                  {
00090   fCenterPosition = G4ThreeVector(centerPosition[0], centerPosition[1], centerPosition[2]);
00091 }

void G4VScoringMesh::SetCurrentPrimitiveScorer ( const G4String name  ) 

Definition at line 164 of file G4VScoringMesh.cc.

References fCurrentPS, G4cerr, G4endl, and GetPrimitiveScorer().

Referenced by G4ScoreQuantityMessenger::SetNewValue().

00164                                                                     {
00165   fCurrentPS = GetPrimitiveScorer(name);
00166   if(fCurrentPS == 0) {
00167     G4cerr << "ERROR : G4VScoringMesh::SetCurrentPrimitiveScorer() : The primitive scorer <"
00168            << name << "> does not found." << G4endl;
00169   }
00170 }

void G4VScoringMesh::SetCurrentPSUnit ( const G4String unit  ) 

Definition at line 199 of file G4VScoringMesh.cc.

References fCurrentPS, G4cerr, G4endl, and G4VPrimitiveScorer::SetUnit().

Referenced by G4ScoreQuantityMessenger::SetNewValue().

00199                                                           {
00200   if(fCurrentPS == 0) {
00201       G4String msg = "ERROR : G4VScoringMesh::GetCurrentPSUnit() : ";
00202       msg += " Current primitive scorer is null.";
00203       G4cerr << msg << G4endl;
00204   }else{
00205       fCurrentPS->SetUnit(unit);
00206   }
00207 }

void G4VScoringMesh::SetDrawPSName ( const G4String psname  )  [inline]

Definition at line 139 of file G4VScoringMesh.hh.

References fDrawPSName.

00139 {fDrawPSName = psname;}

void G4VScoringMesh::SetFilter ( G4VSDFilter filter  ) 

Definition at line 144 of file G4VScoringMesh.cc.

References fCurrentPS, G4cerr, G4cout, G4endl, G4VPrimitiveScorer::GetFilter(), G4VPrimitiveScorer::GetName(), G4VSDFilter::GetName(), G4VPrimitiveScorer::SetFilter(), and verboseLevel.

Referenced by G4ScoreQuantityMessenger::FParticleCommand(), G4ScoreQuantityMessenger::FParticleWithEnergyCommand(), and G4ScoreQuantityMessenger::SetNewValue().

00144                                                    {
00145 
00146   if(fCurrentPS == 0) {
00147     G4cerr << "ERROR : G4VScoringMesh::SetSDFilter() : a quantity must be defined first. This method is ignored." << G4endl;
00148     return;
00149   }
00150   if(verboseLevel > 0) G4cout << "G4VScoringMesh::SetFilter() : "
00151                               << filter->GetName()
00152                               << " is set to "
00153                               << fCurrentPS->GetName() << G4endl;
00154 
00155   G4VSDFilter* oldFilter = fCurrentPS->GetFilter();
00156   if(oldFilter)
00157   {
00158     G4cout << "WARNING : G4VScoringMesh::SetFilter() : " << oldFilter->GetName() 
00159            << " is overwritten by " << filter->GetName() << G4endl;
00160   }
00161   fCurrentPS->SetFilter(filter);
00162 }

void G4VScoringMesh::SetNullToCurrentPrimitiveScorer (  )  [inline]

Definition at line 145 of file G4VScoringMesh.hh.

References fCurrentPS.

Referenced by G4ScoreQuantityMessenger::CheckMeshPS().

00145 {fCurrentPS = NULL;}

void G4VScoringMesh::SetNumberOfSegments ( G4int  nSegment[3]  ) 

Definition at line 92 of file G4VScoringMesh.cc.

References fNSegment, G4Exception(), JustWarning, and nMeshIsSet.

00092                                                           {
00093   if ( !nMeshIsSet ){
00094     for(int i = 0; i < 3; i++) fNSegment[i] = nSegment[i];
00095     nMeshIsSet = true;
00096   } else {
00097     G4String message = "   The size of scoring segments can not be changed.";
00098     G4Exception("G4VScoringMesh::SetNumberOfSegments()",
00099                 "DigiHitsUtilsScoreVScoringMesh000", JustWarning,
00100                 message);
00101   }
00102 }

void G4VScoringMesh::SetPrimitiveScorer ( G4VPrimitiveScorer ps  ) 

Definition at line 121 of file G4VScoringMesh.cc.

References fCurrentPS, fMap, fMFD, fNSegment, fWorldName, G4cerr, G4cout, G4endl, G4VPrimitiveScorer::GetName(), ReadyForQuantity(), G4MultiFunctionalDetector::RegisterPrimitive(), G4VPrimitiveScorer::SetNijk(), and verboseLevel.

00121                                                                {
00122 
00123   if(!ReadyForQuantity())
00124   {
00125     G4cerr << "ERROR : G4VScoringMesh::SetPrimitiveScorer() : " << ps->GetName() 
00126            << " does not yet have mesh size or number of bins. Set them first." << G4endl
00127            << "This Method is ignored." << G4endl;
00128     return;
00129   }
00130   if(verboseLevel > 0) G4cout << "G4VScoringMesh::SetPrimitiveScorer() : "
00131                               << ps->GetName() << " is registered."
00132                               << " 3D size: ("
00133                               << fNSegment[0] << ", "
00134                               << fNSegment[1] << ", "
00135                               << fNSegment[2] << ")" << G4endl;
00136 
00137   ps->SetNijk(fNSegment[0], fNSegment[1], fNSegment[2]);
00138   fCurrentPS = ps;
00139   fMFD->RegisterPrimitive(ps);
00140   G4THitsMap<G4double> * map = new G4THitsMap<G4double>(fWorldName, ps->GetName());
00141   fMap[ps->GetName()] = map;
00142 }

void G4VScoringMesh::SetSize ( G4double  size[3]  ) 

Definition at line 72 of file G4VScoringMesh.cc.

References fSize, G4Exception(), JustWarning, and sizeIsSet.

Referenced by G4ScoringMessenger::SetNewValue().

00072                                              {
00073   if ( !sizeIsSet ){
00074     for(int i = 0; i < 3; i++) fSize[i] = size[i];
00075     sizeIsSet = true;
00076   }else{
00077     G4String message = "   The size of scoring mesh can not be changed.";
00078     G4Exception("G4VScoringMesh::SetSize()",
00079                 "DigiHitsUtilsScoreVScoringMesh000", JustWarning,
00080                 message);
00081   }
00082 }

void G4VScoringMesh::SetVerboseLevel ( G4int  vl  )  [inline]

Definition at line 147 of file G4VScoringMesh.hh.

References verboseLevel.

Referenced by G4ScoringManager::RegisterScoringMesh().

00148   { verboseLevel = vl; }


Field Documentation

G4bool G4VScoringMesh::fActive [protected]

Definition at line 163 of file G4VScoringMesh.hh.

Referenced by Activate(), and IsActive().

G4ThreeVector G4VScoringMesh::fCenterPosition [protected]

Definition at line 167 of file G4VScoringMesh.hh.

Referenced by G4ScoringCylinder::Draw(), G4ScoringBox::Draw(), G4ScoringCylinder::DrawColumn(), G4ScoringBox::DrawColumn(), GetTranslation(), List(), and SetCenterPosition().

G4bool G4VScoringMesh::fConstructed [protected]

Definition at line 162 of file G4VScoringMesh.hh.

Referenced by G4ScoringCylinder::Construct(), and G4ScoringBox::Construct().

G4VPrimitiveScorer* G4VScoringMesh::fCurrentPS [protected]

Definition at line 161 of file G4VScoringMesh.hh.

Referenced by GetCurrentPSUnit(), IsCurrentPrimitiveScorerNull(), SetCurrentPrimitiveScorer(), SetCurrentPSUnit(), SetFilter(), SetNullToCurrentPrimitiveScorer(), and SetPrimitiveScorer().

G4String G4VScoringMesh::fDivisionAxisNames[3] [protected]

Definition at line 183 of file G4VScoringMesh.hh.

Referenced by G4ScoringBox::G4ScoringBox(), G4ScoringCylinder::G4ScoringCylinder(), G4VScoringMesh(), and GetDivisionAxisNames().

G4String G4VScoringMesh::fDrawPSName [protected]

Definition at line 181 of file G4VScoringMesh.hh.

Referenced by G4ScoringCylinder::Draw(), G4ScoringBox::Draw(), G4ScoringCylinder::DrawColumn(), G4ScoringBox::DrawColumn(), DrawMesh(), and SetDrawPSName().

G4String G4VScoringMesh::fDrawUnit [protected]

Definition at line 179 of file G4VScoringMesh.hh.

Referenced by G4ScoringCylinder::Draw(), G4ScoringBox::Draw(), G4ScoringCylinder::DrawColumn(), G4ScoringBox::DrawColumn(), and DrawMesh().

G4double G4VScoringMesh::fDrawUnitValue [protected]

Definition at line 180 of file G4VScoringMesh.hh.

Referenced by G4ScoringCylinder::Draw(), G4ScoringBox::Draw(), G4ScoringCylinder::DrawColumn(), G4ScoringBox::DrawColumn(), and DrawMesh().

std::map<G4String, G4THitsMap<G4double>* > G4VScoringMesh::fMap [protected]

Definition at line 171 of file G4VScoringMesh.hh.

Referenced by Accumulate(), DrawMesh(), Dump(), FindPrimitiveScorer(), GetPSUnit(), GetPSUnitValue(), GetScoreMap(), ResetScore(), and SetPrimitiveScorer().

G4MultiFunctionalDetector* G4VScoringMesh::fMFD [protected]

Definition at line 172 of file G4VScoringMesh.hh.

Referenced by G4VScoringMesh(), GetPrimitiveScorer(), List(), and SetPrimitiveScorer().

G4int G4VScoringMesh::fNSegment[3] [protected]

Definition at line 169 of file G4VScoringMesh.hh.

Referenced by G4ScoringCylinder::Draw(), G4ScoringBox::Draw(), G4ScoringCylinder::DrawColumn(), G4ScoringBox::DrawColumn(), G4VScoringMesh(), GetNumberOfSegments(), List(), SetNumberOfSegments(), and SetPrimitiveScorer().

G4RotationMatrix* G4VScoringMesh::fRotationMatrix [protected]

Definition at line 168 of file G4VScoringMesh.hh.

Referenced by G4ScoringCylinder::Draw(), G4ScoringBox::Draw(), G4ScoringCylinder::DrawColumn(), G4ScoringBox::DrawColumn(), GetRotationMatrix(), List(), RotateX(), RotateY(), and RotateZ().

MeshShape G4VScoringMesh::fShape [protected]

Definition at line 164 of file G4VScoringMesh.hh.

Referenced by G4ScoringBox::G4ScoringBox(), G4ScoringCylinder::G4ScoringCylinder(), and GetShape().

G4double G4VScoringMesh::fSize[3] [protected]

Definition at line 166 of file G4VScoringMesh.hh.

Referenced by G4ScoringCylinder::Draw(), G4ScoringBox::Draw(), G4ScoringCylinder::DrawColumn(), G4ScoringBox::DrawColumn(), G4VScoringMesh(), GetSize(), G4ScoringCylinder::List(), G4ScoringBox::List(), G4ScoringCylinder::SetRMax(), SetSize(), and G4ScoringCylinder::SetZSize().

G4String G4VScoringMesh::fWorldName [protected]

Definition at line 160 of file G4VScoringMesh.hh.

Referenced by Dump(), GetWorldName(), G4ScoringCylinder::List(), G4ScoringBox::List(), and SetPrimitiveScorer().

G4bool G4VScoringMesh::nMeshIsSet [protected]

Definition at line 177 of file G4VScoringMesh.hh.

Referenced by ReadyForQuantity(), and SetNumberOfSegments().

G4bool G4VScoringMesh::sizeIsSet [protected]

Definition at line 176 of file G4VScoringMesh.hh.

Referenced by GetSize(), ReadyForQuantity(), and SetSize().

G4int G4VScoringMesh::verboseLevel [protected]

Definition at line 174 of file G4VScoringMesh.hh.

Referenced by Accumulate(), G4ScoringCylinder::Construct(), G4ScoringBox::Construct(), ResetScore(), SetFilter(), SetPrimitiveScorer(), and SetVerboseLevel().


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:55 2013 for Geant4 by  doxygen 1.4.7