Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Types | Public Member Functions
G4ScoringCylinder Class Reference

#include <G4ScoringCylinder.hh>

Inheritance diagram for G4ScoringCylinder:
G4VScoringMesh

Public Types

enum  IDX { IZ, IPHI, IR }
 

Public Member Functions

 G4ScoringCylinder (G4String wName)
 
 ~G4ScoringCylinder ()
 
virtual void Construct (G4VPhysicalVolume *fWorldPhys)
 
virtual void List () const
 
virtual void Draw (std::map< G4int, G4double * > *map, G4VScoreColorMap *colorMap, G4int axflg=111)
 
virtual void DrawColumn (std::map< G4int, G4double * > *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn)
 
void SetRMax (G4double rMax)
 
void SetZSize (G4double zSize)
 
void RegisterPrimitives (std::vector< G4VPrimitiveScorer * > &vps)
 
void GetRZPhi (G4int index, G4int q[3]) const
 
- Public Member Functions inherited from G4VScoringMesh
 G4VScoringMesh (const G4String &wName)
 
virtual ~G4VScoringMesh ()
 
virtual void WorkerConstruct (G4VPhysicalVolume *fWorldPhys)
 
const G4StringGetWorldName () const
 
G4bool IsActive () const
 
void Activate (G4bool vl=true)
 
MeshShape GetShape () const
 
void Accumulate (G4THitsMap< G4double > *map)
 
void Merge (const G4VScoringMesh *scMesh)
 
void Dump ()
 
void DrawMesh (const G4String &psName, G4VScoreColorMap *colorMap, G4int axflg=111)
 
void DrawMesh (const G4String &psName, G4int idxPlane, G4int iColumn, G4VScoreColorMap *colorMap)
 
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 () const
 
G4bool ReadyForQuantity () const
 
void SetMeshElementLogical (G4LogicalVolume *val)
 
G4LogicalVolumeGetMeshElementLogical () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VScoringMesh
G4VPrimitiveScorerGetPrimitiveScorer (const G4String &name)
 
- Protected Attributes inherited from G4VScoringMesh
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]
 
G4LogicalVolumefMeshElementLogical
 

Detailed Description

Definition at line 41 of file G4ScoringCylinder.hh.

Member Enumeration Documentation

Constructor & Destructor Documentation

G4ScoringCylinder::G4ScoringCylinder ( G4String  wName)

Definition at line 55 of file G4ScoringCylinder.cc.

References cylinderMesh, G4VScoringMesh::fDivisionAxisNames, and G4VScoringMesh::fShape.

56  :G4VScoringMesh(wName)
57 {
59 
60  fDivisionAxisNames[0] = "Z";
61  fDivisionAxisNames[1] = "PHI";
62  fDivisionAxisNames[2] = "R";
63 }
G4VScoringMesh(const G4String &wName)
G4String fDivisionAxisNames[3]
G4ScoringCylinder::~G4ScoringCylinder ( )

Definition at line 65 of file G4ScoringCylinder.cc.

66 {;}

Member Function Documentation

void G4ScoringCylinder::Construct ( G4VPhysicalVolume fWorldPhys)
virtual

Implements G4VScoringMesh.

Definition at line 68 of file G4ScoringCylinder.cc.

References G4VScoringMesh::fConstructed, G4cout, G4endl, G4VPhysicalVolume::GetName(), G4VScoringMesh::ResetScore(), and G4VScoringMesh::verboseLevel.

69 {
70  if(fConstructed) {
71 
72  if(verboseLevel > 0)
73  G4cout << fWorldPhys->GetName() << " --- All quantities are reset." << G4endl;
74  ResetScore();
75 
76  } else {
77  fConstructed = true;
78  SetupGeometry(fWorldPhys);
79  }
80 }
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
#define G4endl
Definition: G4ios.hh:61
void G4ScoringCylinder::Draw ( std::map< G4int, G4double * > *  map,
G4VScoreColorMap colorMap,
G4int  axflg = 111 
)
virtual

Implements G4VScoringMesh.

Definition at line 223 of file G4ScoringCylinder.cc.

References test::c, DBL_MAX, G4VVisManager::Draw(), G4VScoreColorMap::DrawColorChart(), G4VScoringMesh::fCenterPosition, G4VScoringMesh::fDrawPSName, G4VScoringMesh::fDrawUnit, G4VScoringMesh::fDrawUnitValue, G4VScoringMesh::fNSegment, G4VScoringMesh::fRotationMatrix, G4VScoringMesh::fSize, G4cout, G4endl, G4VVisManager::GetConcreteInstance(), G4VScoreColorMap::GetMapColor(), GetRZPhi(), G4VScoreColorMap::IfFloatMinMax(), HepGeom::Transform3D::inverse(), IPHI, IR, iz, IZ, G4VisAttributes::SetColour(), G4VisAttributes::SetForceAuxEdgeVisible(), G4VisAttributes::SetForceSolid(), G4VScoreColorMap::SetMinMax(), G4VScoreColorMap::SetPSName(), G4VScoreColorMap::SetPSUnit(), python.hepunit::twopi, and z.

223  {
224 
226  if(pVisManager) {
227 
228  // cell vectors
229  std::vector<double> ephi;
230  for(int phi = 0; phi < fNSegment[IPHI]; phi++) ephi.push_back(0.);
231  //-
232  std::vector<std::vector<double> > zphicell; // zphicell[Z][PHI]
233  for(int z = 0; z < fNSegment[IZ]; z++) zphicell.push_back(ephi);
234  //-
235  std::vector<std::vector<double> > rphicell; // rphicell[R][PHI]
236  for(int r = 0; r < fNSegment[IR]; r++) rphicell.push_back(ephi);
237 
238  // projections
239  G4int q[3];
240  std::map<G4int, G4double*>::iterator itr = map->begin();
241  for(; itr != map->end(); itr++) {
242  if(itr->first < 0) {
243  G4cout << itr->first << G4endl;
244  continue;
245  }
246  GetRZPhi(itr->first, q);
247 
248  zphicell[q[IZ]][q[IPHI]] += *(itr->second)/fDrawUnitValue;
249  rphicell[q[IR]][q[IPHI]] += *(itr->second)/fDrawUnitValue;
250  }
251 
252  // search min./max. values
253  G4double zphimin = DBL_MAX, rphimin = DBL_MAX;
254  G4double zphimax = 0., rphimax = 0.;
255  for(int iphi = 0; iphi < fNSegment[IPHI]; iphi++) {
256  for(int iz = 0; iz < fNSegment[IZ]; iz++) {
257  if(zphimin > zphicell[iz][iphi]) zphimin = zphicell[iz][iphi];
258  if(zphimax < zphicell[iz][iphi]) zphimax = zphicell[iz][iphi];
259  }
260  for(int ir = 0; ir < fNSegment[IR]; ir++) {
261  if(rphimin > rphicell[ir][iphi]) rphimin = rphicell[ir][iphi];
262  if(rphimax < rphicell[ir][iphi]) rphimax = rphicell[ir][iphi];
263  }
264  }
265 
266  G4VisAttributes att;
267  att.SetForceSolid(true);
268  att.SetForceAuxEdgeVisible(true);
269 
270 
271  G4Scale3D scale;
272  if(axflg/100==1) {
273  // rz plane
274  }
275  axflg = axflg%100;
276  if(axflg/10==1) {
277  // z-phi plane
278  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(zphimin, zphimax); }
279 
280  G4double zhalf = fSize[1]/fNSegment[IZ];
281  for(int phi = 0; phi < fNSegment[IPHI]; phi++) {
282  for(int z = 0; z < fNSegment[IZ]; z++) {
283  //-
284  G4double angle = twopi/fNSegment[IPHI]*phi;
285  G4double dphi = twopi/fNSegment[IPHI];
286  G4Tubs cylinder("z-phi", // name
287  fSize[0]*0.99, fSize[0], // inner radius, outer radius
288  zhalf, // half length in z
289  angle, dphi*0.99999); // starting phi angle, delta angle
290  //-
291  G4ThreeVector zpos(0., 0., -fSize[1] + fSize[1]/fNSegment[IZ]*(1 + 2.*z));
292  G4Transform3D trans;
293  if(fRotationMatrix) {
294  trans = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(zpos);
295  trans = G4Translate3D(fCenterPosition)*trans;
296  } else {
298  }
299  G4double c[4];
300  colorMap->GetMapColor(zphicell[z][phi], c);
301  att.SetColour(c[0], c[1], c[2]);//, c[3]);
302  //-
303  pVisManager->Draw(cylinder, att, trans);
304  }
305  }
306  }
307  axflg = axflg%10;
308  if(axflg==1) {
309  // r-phi plane
310  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(rphimin, rphimax); }
311 
312  G4double rsize = fSize[0]/fNSegment[IR];
313  for(int phi = 0; phi < fNSegment[IPHI]; phi++) {
314  for(int r = 0; r < fNSegment[IR]; r++) {
315 
316  G4double rs[2] = {rsize*r, rsize*(r+1)};
317  G4double angle = twopi/fNSegment[IPHI]*phi;
318  G4double dphi = twopi/fNSegment[IPHI];
319  G4Tubs cylinder("z-phi", rs[0], rs[1], 0.001,
320  angle, dphi*0.99999);
321  /*
322  G4cout << ">>>> "
323  << rs[0] << " - " << rs[1] << " : "
324  << angle << " - " << angle + dphi
325  << G4endl;
326  */
327 
328  G4ThreeVector zposn(0., 0., -fSize[1]);
329  G4ThreeVector zposp(0., 0., fSize[1]);
330  G4Transform3D transn, transp;
331  if(fRotationMatrix) {
332  transn = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(zposn);
333  transn = G4Translate3D(fCenterPosition)*transn;
334  transp = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(zposp);
335  transp = G4Translate3D(fCenterPosition)*transp;
336  } else {
339  }
340  G4double c[4];
341  colorMap->GetMapColor(rphicell[r][phi], c);
342  att.SetColour(c[0], c[1], c[2]);//, c[3]);
343  /*
344  G4cout << " " << c[0] << ", "
345  << c[1] << ", " << c[2] << G4endl;
346  */
347  pVisManager->Draw(cylinder, att, transn);
348  pVisManager->Draw(cylinder, att, transp);
349  }
350  }
351  }
352 
353  colorMap->SetPSUnit(fDrawUnit);
354  colorMap->SetPSName(fDrawPSName);
355  colorMap->DrawColorChart();
356 
357  }
358 }
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
void SetColour(const G4Colour &)
void SetPSName(G4String &psName)
G4ThreeVector fCenterPosition
static G4VVisManager * GetConcreteInstance()
G4double z
Definition: TRTMaterials.hh:39
Definition: G4Tubs.hh:84
void SetPSUnit(G4String &unit)
virtual void GetMapColor(G4double val, G4double color[4])=0
Transform3D inverse() const
Definition: Transform3D.cc:142
G4double fSize[3]
void SetForceSolid(G4bool)
int G4int
Definition: G4Types.hh:78
void SetMinMax(G4double minVal, G4double maxVal)
G4double fDrawUnitValue
G4GLOB_DLL std::ostream G4cout
G4RotationMatrix * fRotationMatrix
G4double iz
Definition: TRTMaterials.hh:39
void SetForceAuxEdgeVisible(G4bool)
HepGeom::Rotate3D G4Rotate3D
HepGeom::Translate3D G4Translate3D
G4bool IfFloatMinMax() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
G4String fDrawPSName
virtual void DrawColorChart(G4int nPoint=5)
void GetRZPhi(G4int index, G4int q[3]) const
void G4ScoringCylinder::DrawColumn ( std::map< G4int, G4double * > *  map,
G4VScoreColorMap colorMap,
G4int  idxProj,
G4int  idxColumn 
)
virtual

Implements G4VScoringMesh.

Definition at line 360 of file G4ScoringCylinder.cc.

References test::c, DBL_MAX, G4VVisManager::Draw(), G4VScoreColorMap::DrawColorChart(), G4VScoringMesh::fCenterPosition, G4VScoringMesh::fDrawPSName, G4VScoringMesh::fDrawUnit, G4VScoringMesh::fDrawUnitValue, G4VScoringMesh::fNSegment, G4VScoringMesh::fRotationMatrix, G4VScoringMesh::fSize, G4cerr, G4cout, G4endl, G4VVisManager::GetConcreteInstance(), G4VScoreColorMap::GetMapColor(), GetRZPhi(), G4VScoreColorMap::IfFloatMinMax(), HepGeom::Transform3D::inverse(), IPHI, IR, IZ, python.hepunit::radian, G4VisAttributes::SetColour(), G4VisAttributes::SetForceAuxEdgeVisible(), G4VisAttributes::SetForceSolid(), G4VScoreColorMap::SetMinMax(), G4VScoreColorMap::SetPSName(), G4VScoreColorMap::SetPSUnit(), python.hepunit::twopi, and z.

362 {
363  G4int projAxis = 0;
364  switch(idxProj) {
365  case 0:
366  projAxis = IR;
367  break;
368  case 1:
369  projAxis = IZ;
370  break;
371  case 2:
372  projAxis = IPHI;
373  break;
374  }
375 
376  if(idxColumn<0 || idxColumn>=fNSegment[projAxis])
377  {
378  G4cerr << "Warning : Column number " << idxColumn << " is out of scoring mesh [0," << fNSegment[projAxis]-1 <<
379  "]. Method ignored." << G4endl;
380  return;
381  }
383  if(pVisManager) {
384 
385  // cell vectors
386  std::vector<std::vector<std::vector<double> > > cell; // cell[R][Z][PHI]
387  std::vector<double> ephi;
388  for(int phi = 0; phi < fNSegment[IPHI]; phi++) ephi.push_back(0.);
389  std::vector<std::vector<double> > ezphi;
390  for(int z = 0; z < fNSegment[IZ]; z++) ezphi.push_back(ephi);
391  for(int r = 0; r < fNSegment[IR]; r++) cell.push_back(ezphi);
392 
393  std::vector<std::vector<double> > rzcell; // rzcell[R][Z]
394  std::vector<double> ez;
395  for(int z = 0; z < fNSegment[IZ]; z++) ez.push_back(0.);
396  for(int r = 0; r < fNSegment[IR]; r++) rzcell.push_back(ez);
397 
398  std::vector<std::vector<double> > zphicell; // zphicell[Z][PHI]
399  for(int z = 0; z < fNSegment[IZ]; z++) zphicell.push_back(ephi);
400 
401  std::vector<std::vector<double> > rphicell; // rphicell[R][PHI]
402  for(int r = 0; r < fNSegment[IR]; r++) rphicell.push_back(ephi);
403 
404  // projections
405  G4int q[3];
406  std::map<G4int, G4double*>::iterator itr = map->begin();
407  for(; itr != map->end(); itr++) {
408  if(itr->first < 0) {
409  G4cout << itr->first << G4endl;
410  continue;
411  }
412  GetRZPhi(itr->first, q);
413 
414  if(projAxis == IR && q[IR] == idxColumn) { // zphi plane
415  zphicell[q[IZ]][q[IPHI]] += *(itr->second)/fDrawUnitValue;
416  }
417  if(projAxis == IZ && q[IZ] == idxColumn) { // rphi plane
418  rphicell[q[IR]][q[IPHI]] += *(itr->second)/fDrawUnitValue;
419  }
420  if(projAxis == IPHI && q[IPHI] == idxColumn) { // rz plane
421  rzcell[q[IR]][q[IZ]] += *(itr->second)/fDrawUnitValue;
422  }
423  }
424 
425  // search min./max. values
426  G4double rzmin = DBL_MAX, zphimin = DBL_MAX, rphimin = DBL_MAX;
427  G4double rzmax = 0., zphimax = 0., rphimax = 0.;
428  for(int r = 0; r < fNSegment[IR]; r++) {
429  for(int phi = 0; phi < fNSegment[IPHI]; phi++) {
430  if(rphimin > rphicell[r][phi]) rphimin = rphicell[r][phi];
431  if(rphimax < rphicell[r][phi]) rphimax = rphicell[r][phi];
432  }
433  for(int z = 0; z < fNSegment[IZ]; z++) {
434  if(rzmin > rzcell[r][z]) rzmin = rzcell[r][z];
435  if(rzmax < rzcell[r][z]) rzmax = rzcell[r][z];
436  }
437  }
438  for(int z = 0; z < fNSegment[IZ]; z++) {
439  for(int phi = 0; phi < fNSegment[IPHI]; phi++) {
440  if(zphimin > zphicell[z][phi]) zphimin = zphicell[z][phi];
441  if(zphimax < zphicell[z][phi]) zphimax = zphicell[z][phi];
442  }
443  }
444 
445 
446  G4VisAttributes att;
447  att.SetForceSolid(true);
448  att.SetForceAuxEdgeVisible(true);
449 
450 
451  G4Scale3D scale;
452  // z-phi plane
453  if(projAxis == IR) {
454  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(zphimin,zphimax); }
455 
456  G4double zhalf = fSize[1]/fNSegment[IZ];
457  G4double rsize[2] = {fSize[0]/fNSegment[IR]*idxColumn,
458  fSize[0]/fNSegment[IR]*(idxColumn+1)};
459  for(int phi = 0; phi < fNSegment[IPHI]; phi++) {
460  for(int z = 0; z < fNSegment[IZ]; z++) {
461 
462  G4double angle = twopi/fNSegment[IPHI]*phi*radian;
463  G4double dphi = twopi/fNSegment[IPHI]*radian;
464  G4Tubs cylinder("z-phi", rsize[0], rsize[1], zhalf,
465  angle, dphi*0.99999);
466 
467  G4ThreeVector zpos(0., 0., -fSize[1] + fSize[1]/fNSegment[IZ]*(1 + 2.*z));
468  G4Transform3D trans;
469  if(fRotationMatrix) {
470  trans = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(zpos);
471  trans = G4Translate3D(fCenterPosition)*trans;
472  } else {
474  }
475  G4double c[4];
476  colorMap->GetMapColor(zphicell[z][phi], c);
477  att.SetColour(c[0], c[1], c[2]);//, c[3]);
478  pVisManager->Draw(cylinder, att, trans);
479  }
480  }
481 
482  // r-phi plane
483  } else if(projAxis == IZ) {
484  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(rphimin,rphimax); }
485 
486  G4double rsize = fSize[0]/fNSegment[IR];
487  for(int phi = 0; phi < fNSegment[IPHI]; phi++) {
488  for(int r = 0; r < fNSegment[IR]; r++) {
489 
490  G4double rs[2] = {rsize*r, rsize*(r+1)};
491  G4double angle = twopi/fNSegment[IPHI]*phi*radian;
492  G4double dz = fSize[1]/fNSegment[IZ];
493  G4double dphi = twopi/fNSegment[IPHI]*radian;
494  G4Tubs cylinder("r-phi", rs[0], rs[1], dz,
495  angle, dphi*0.99999);
496  G4ThreeVector zpos(0., 0.,
497  -fSize[1]+fSize[1]/fNSegment[IZ]*(idxColumn*2+1));
498  G4Transform3D trans;
499  if(fRotationMatrix) {
500  trans = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(zpos);
501  trans = G4Translate3D(fCenterPosition)*trans;
502  } else {
504  }
505  G4double c[4];
506  colorMap->GetMapColor(rphicell[r][phi], c);
507  att.SetColour(c[0], c[1], c[2]);//, c[3]);
508  pVisManager->Draw(cylinder, att, trans);
509  }
510  }
511 
512  // r-z plane
513  } else if(projAxis == IPHI) {
514  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(rzmin,rzmax); }
515 
516  G4double rsize = fSize[0]/fNSegment[IR];
517  G4double zhalf = fSize[1]/fNSegment[IZ];
518  G4double angle = twopi/fNSegment[IPHI]*idxColumn*radian;
519  G4double dphi = twopi/fNSegment[IPHI]*radian;
520  for(int z = 0; z < fNSegment[IZ]; z++) {
521  for(int r = 0; r < fNSegment[IR]; r++) {
522 
523  G4double rs[2] = {rsize*r, rsize*(r+1)};
524  G4Tubs cylinder("z-phi", rs[0], rs[1], zhalf,
525  angle, dphi);
526 
527  G4ThreeVector zpos(0., 0.,
528  -fSize[1]+fSize[1]/fNSegment[IZ]*(2.*z+1));
529  G4Transform3D trans;
530  if(fRotationMatrix) {
531  trans = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(zpos);
532  trans = G4Translate3D(fCenterPosition)*trans;
533  } else {
535  }
536  G4double c[4];
537  colorMap->GetMapColor(rzcell[r][z], c);
538  att.SetColour(c[0], c[1], c[2]);//, c[3]);
539  pVisManager->Draw(cylinder, att, trans);
540  }
541  }
542  }
543  }
544 
545  colorMap->SetPSUnit(fDrawUnit);
546  colorMap->SetPSName(fDrawPSName);
547  colorMap->DrawColorChart();
548 
549 }
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
void SetColour(const G4Colour &)
void SetPSName(G4String &psName)
G4ThreeVector fCenterPosition
static G4VVisManager * GetConcreteInstance()
G4double z
Definition: TRTMaterials.hh:39
Definition: G4Tubs.hh:84
void SetPSUnit(G4String &unit)
virtual void GetMapColor(G4double val, G4double color[4])=0
Transform3D inverse() const
Definition: Transform3D.cc:142
G4double fSize[3]
void SetForceSolid(G4bool)
int G4int
Definition: G4Types.hh:78
void SetMinMax(G4double minVal, G4double maxVal)
G4double fDrawUnitValue
G4GLOB_DLL std::ostream G4cout
G4RotationMatrix * fRotationMatrix
void SetForceAuxEdgeVisible(G4bool)
HepGeom::Rotate3D G4Rotate3D
HepGeom::Translate3D G4Translate3D
G4bool IfFloatMinMax() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
G4String fDrawPSName
G4GLOB_DLL std::ostream G4cerr
virtual void DrawColorChart(G4int nPoint=5)
void GetRZPhi(G4int index, G4int q[3]) const
void G4ScoringCylinder::GetRZPhi ( G4int  index,
G4int  q[3] 
) const

Definition at line 551 of file G4ScoringCylinder.cc.

References IPHI, IR, and IZ.

Referenced by Draw(), and DrawColumn().

551  {
552  // index = k + j * k-size + i * jk-plane-size
553 
554  // nested : z -> phi -> r
555  G4int i = IZ;
556  G4int j = IPHI;
557  G4int k = IR;
558  G4int jk = fNSegment[j]*fNSegment[k];
559  q[i] = index/jk;
560  q[j] = (index - q[i]*jk)/fNSegment[k];
561  q[k] = index - q[j]*fNSegment[k] - q[i]*jk;
562 }
int G4int
Definition: G4Types.hh:78
void G4ScoringCylinder::List ( ) const
virtual

Reimplemented from G4VScoringMesh.

Definition at line 211 of file G4ScoringCylinder.cc.

References python.hepunit::cm, G4VScoringMesh::fSize, G4VScoringMesh::fWorldName, G4cout, G4endl, and G4VScoringMesh::List().

211  {
212  G4cout << "G4ScoringCylinder : " << fWorldName << " --- Shape: Cylindrical mesh" << G4endl;
213 
214  G4cout << " Size (R, Dz): ("
215  << fSize[0]/cm << ", "
216  << fSize[1]/cm << ") [cm]"
217  << G4endl;
218 
220 }
G4double fSize[3]
G4GLOB_DLL std::ostream G4cout
virtual void List() const
#define G4endl
Definition: G4ios.hh:61
void G4ScoringCylinder::RegisterPrimitives ( std::vector< G4VPrimitiveScorer * > &  vps)
void G4ScoringCylinder::SetRMax ( G4double  rMax)
inline

Definition at line 55 of file G4ScoringCylinder.hh.

References G4VScoringMesh::fSize.

55 {fSize[0] = rMax;}
G4double fSize[3]
void G4ScoringCylinder::SetZSize ( G4double  zSize)
inline

Definition at line 56 of file G4ScoringCylinder.hh.

References G4VScoringMesh::fSize.

56 {fSize[1] = zSize;} // half height
G4double fSize[3]

The documentation for this class was generated from the following files: