Geant4  10.02.p03
G4ScoringBox Class Reference

#include <G4ScoringBox.hh>

Inheritance diagram for G4ScoringBox:
Collaboration diagram for G4ScoringBox:

Public Member Functions

 G4ScoringBox (G4String wName)
 
 ~G4ScoringBox ()
 
void List () const
 
void Draw (std::map< G4int, G4double *> *map, G4VScoreColorMap *colorMap, G4int axflg=111)
 
void DrawColumn (std::map< G4int, G4double *> *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn)
 
void SetSegmentDirection (G4int dir)
 
- Public Member Functions inherited from G4VScoringMesh
 G4VScoringMesh (const G4String &wName)
 
virtual ~G4VScoringMesh ()
 
void Construct (G4VPhysicalVolume *fWorldPhys)
 
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
 
void SetParallelWorldProcess (G4ParallelWorldProcess *proc)
 
G4ParallelWorldProcessGetParallelWorldProcess () const
 
void GeometryHasBeenDestroyed ()
 

Protected Member Functions

virtual void SetupGeometry (G4VPhysicalVolume *fWorldPhys)
 
- Protected Member Functions inherited from G4VScoringMesh
G4VPrimitiveScorerGetPrimitiveScorer (const G4String &name)
 

Private Member Functions

G4ThreeVector GetReplicaPosition (G4int x, G4int y, G4int z)
 
void GetXYZ (G4int index, G4int q[3]) const
 
G4int GetIndex (G4int x, G4int y, G4int z) const
 

Private Attributes

G4int fSegmentDirection
 

Additional Inherited Members

- 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
 
G4ParallelWorldProcessfParallelWorldProcess
 
G4bool fGeometryHasBeenDestroyed
 

Detailed Description

Definition at line 43 of file G4ScoringBox.hh.

Constructor & Destructor Documentation

◆ G4ScoringBox()

G4ScoringBox::G4ScoringBox ( G4String  wName)

Definition at line 53 of file G4ScoringBox.cc.

55 {
56  fShape = boxMesh;
57  fDivisionAxisNames[0] = "X";
58  fDivisionAxisNames[1] = "Y";
59  fDivisionAxisNames[2] = "Z";
60 }
G4VScoringMesh(const G4String &wName)
G4int fSegmentDirection
Definition: G4ScoringBox.hh:71
G4String fDivisionAxisNames[3]

◆ ~G4ScoringBox()

G4ScoringBox::~G4ScoringBox ( )

Definition at line 62 of file G4ScoringBox.cc.

63 {
64 }

Member Function Documentation

◆ Draw()

void G4ScoringBox::Draw ( std::map< G4int, G4double *> *  map,
G4VScoreColorMap colorMap,
G4int  axflg = 111 
)
virtual

Implements G4VScoringMesh.

Definition at line 239 of file G4ScoringBox.cc.

240  {
241 
243  if(pVisManager) {
244 
245  // cell vectors
246  std::vector<std::vector<std::vector<double> > > cell; // cell[X][Y][Z]
247  std::vector<double> ez;
248  for(int z = 0; z < fNSegment[2]; z++) ez.push_back(0.);
249  std::vector<std::vector<double> > eyz;
250  for(int y = 0; y < fNSegment[1]; y++) eyz.push_back(ez);
251  for(int x = 0; x < fNSegment[0]; x++) cell.push_back(eyz);
252 
253  std::vector<std::vector<double> > xycell; // xycell[X][Y]
254  std::vector<double> ey;
255  for(int y = 0; y < fNSegment[1]; y++) ey.push_back(0.);
256  for(int x = 0; x < fNSegment[0]; x++) xycell.push_back(ey);
257 
258  std::vector<std::vector<double> > yzcell; // yzcell[Y][Z]
259  for(int y = 0; y < fNSegment[1]; y++) yzcell.push_back(ez);
260 
261  std::vector<std::vector<double> > xzcell; // xzcell[X][Z]
262  for(int x = 0; x < fNSegment[0]; x++) xzcell.push_back(ez);
263 
264  // projections
265  G4int q[3];
266  std::map<G4int, G4double*>::iterator itr = map->begin();
267  for(; itr != map->end(); itr++) {
268  GetXYZ(itr->first, q);
269 
270  xycell[q[0]][q[1]] += *(itr->second)/fDrawUnitValue;
271  yzcell[q[1]][q[2]] += *(itr->second)/fDrawUnitValue;
272  xzcell[q[0]][q[2]] += *(itr->second)/fDrawUnitValue;
273  }
274 
275  // search max. & min. values in each slice
276  G4double xymin = DBL_MAX, yzmin = DBL_MAX, xzmin = DBL_MAX;
277  G4double xymax = 0., yzmax = 0., xzmax = 0.;
278  for(int x = 0; x < fNSegment[0]; x++) {
279  for(int y = 0; y < fNSegment[1]; y++) {
280  if(xymin > xycell[x][y]) xymin = xycell[x][y];
281  if(xymax < xycell[x][y]) xymax = xycell[x][y];
282  }
283  for(int z = 0; z < fNSegment[2]; z++) {
284  if(xzmin > xzcell[x][z]) xzmin = xzcell[x][z];
285  if(xzmax < xzcell[x][z]) xzmax = xzcell[x][z];
286  }
287  }
288  for(int y = 0; y < fNSegment[1]; y++) {
289  for(int z = 0; z < fNSegment[2]; z++) {
290  if(yzmin > yzcell[y][z]) yzmin = yzcell[y][z];
291  if(yzmax < yzcell[y][z]) yzmax = yzcell[y][z];
292  }
293  }
294 
295 
296  G4VisAttributes att;
297  att.SetForceSolid(true);
298  att.SetForceAuxEdgeVisible(true);
299  G4double thick = 0.01;
300 
302  if(axflg/100==1) {
303  pVisManager->BeginDraw();
304 
305  // xy plane
306  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(xymin ,xymax); }
307  G4ThreeVector zhalf(0., 0., fSize[2]/fNSegment[2]-thick);
308  for(int x = 0; x < fNSegment[0]; x++) {
309  for(int y = 0; y < fNSegment[1]; y++) {
310 
311  G4ThreeVector pos(GetReplicaPosition(x, y, 0) - zhalf);
312  G4ThreeVector pos2(GetReplicaPosition(x, y, fNSegment[2]-1) + zhalf);
313  G4Transform3D trans, trans2;
314  if(fRotationMatrix) {
315  trans = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(pos);
316  trans = G4Translate3D(fCenterPosition)*trans;
317  trans2 = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(pos2);
318  trans2 = G4Translate3D(fCenterPosition)*trans2;
319  } else {
322  }
323  G4double c[4];
324  colorMap->GetMapColor(xycell[x][y], c);
325  att.SetColour(c[0], c[1], c[2]);//, c[3]);
326 
327  G4Box xyplate("xy", fSize[0]/fNSegment[0], fSize[1]/fNSegment[1],
328  thick);
329  G4Polyhedron * poly = xyplate.GetPolyhedron();
330  poly->Transform(trans);
331  poly->SetVisAttributes(&att);
332  pVisManager->Draw(*poly);
333 
334  G4Box xyplate2 = xyplate;
335  G4Polyhedron * poly2 = xyplate2.GetPolyhedron();
336  poly2->Transform(trans2);
337  poly2->SetVisAttributes(&att);
338  pVisManager->Draw(*poly2);
339 
340  /*
341  G4double nodes[][3] =
342  {{-fSize[0]/fNSegment[0], -fSize[1]/fNSegment[1], 0.},
343  { fSize[0]/fNSegment[0], -fSize[1]/fNSegment[1], 0.},
344  { fSize[0]/fNSegment[0], fSize[1]/fNSegment[1], 0.},
345  {-fSize[0]/fNSegment[0], fSize[1]/fNSegment[1], 0.}};
346  G4int facets[][4] = {{4, 3, 2, 1}};
347  G4int facets2[][4] = {{1, 2, 3, 4}};
348 
349  G4Polyhedron poly, poly2;
350  poly.createPolyhedron(4, 1, nodes, facets);
351  poly.Transform(trans);
352  poly.SetVisAttributes(att);
353  pVisManager->Draw(poly);
354 
355  poly2.createPolyhedron(4, 1, nodes, facets2);
356  poly2.Transform(trans2);
357  poly2.SetVisAttributes(att);
358  pVisManager->Draw(poly2);
359  */
360  }
361  }
362  pVisManager->EndDraw();
363  }
364  axflg = axflg%100;
365  if(axflg/10==1) {
366  pVisManager->BeginDraw();
367 
368  // yz plane
369  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(yzmin, yzmax); }
370  G4ThreeVector xhalf(fSize[0]/fNSegment[0]-thick, 0., 0.);
371  for(int y = 0; y < fNSegment[1]; y++) {
372  for(int z = 0; z < fNSegment[2]; z++) {
373 
374  G4ThreeVector pos(GetReplicaPosition(0, y, z) - xhalf);
375  G4ThreeVector pos2(GetReplicaPosition(fNSegment[0]-1, y, z) + xhalf);
376  G4Transform3D trans, trans2;
377  if(fRotationMatrix) {
378  trans = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(pos);
379  trans = G4Translate3D(fCenterPosition)*trans;
380  trans2 = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(pos2);
381  trans2 = G4Translate3D(fCenterPosition)*trans2;
382  } else {
385  }
386  G4double c[4];
387  colorMap->GetMapColor(yzcell[y][z], c);
388  att.SetColour(c[0], c[1], c[2]);//, c[3]);
389 
390  G4Box yzplate("yz", thick,//fSize[0]/fNSegment[0]*0.001,
391  fSize[1]/fNSegment[1],
392  fSize[2]/fNSegment[2]);
393  G4Polyhedron * poly = yzplate.GetPolyhedron();
394  poly->Transform(trans);
395  poly->SetVisAttributes(&att);
396  pVisManager->Draw(*poly);
397 
398  G4Box yzplate2 = yzplate;
399  G4Polyhedron * poly2 = yzplate2.GetPolyhedron();
400  poly2->Transform(trans2);
401  poly2->SetVisAttributes(&att);
402  pVisManager->Draw(*poly2);
403 
404  /*
405  G4double nodes[][3] =
406  {{0., -fSize[1]/fNSegment[1], -fSize[2]/fNSegment[2]},
407  {0., fSize[1]/fNSegment[1], -fSize[2]/fNSegment[2]},
408  {0., fSize[1]/fNSegment[1], fSize[2]/fNSegment[2]},
409  {0., -fSize[1]/fNSegment[1], fSize[2]/fNSegment[2]}};
410  G4int facets[][4] = {{4, 3, 2, 1}};
411  G4int facets2[][4] = {{1, 2, 3, 4}};
412 
413  G4Polyhedron poly, poly2;
414  poly.createPolyhedron(4, 1, nodes, facets);
415  poly.Transform(trans);
416  poly.SetVisAttributes(att);
417  pVisManager->Draw(poly);
418 
419  poly2.createPolyhedron(4, 1, nodes, facets2);
420  poly2.Transform(trans2);
421  poly2.SetVisAttributes(att);
422  pVisManager->Draw(poly2);
423  */
424  }
425  }
426  pVisManager->EndDraw();
427  }
428  axflg = axflg%10;
429  if(axflg==1) {
430  pVisManager->BeginDraw();
431 
432  // xz plane
433  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(xzmin,xzmax); }
434  G4ThreeVector yhalf(0., fSize[1]/fNSegment[1]-thick, 0.);
435  for(int x = 0; x < fNSegment[0]; x++) {
436  for(int z = 0; z < fNSegment[2]; z++) {
437 
438  G4ThreeVector pos(GetReplicaPosition(x, 0, z) - yhalf);
439  G4ThreeVector pos2(GetReplicaPosition(x, fNSegment[1]-1, z) + yhalf);
440  G4Transform3D trans, trans2;
441  if(fRotationMatrix) {
442  trans = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(pos);
443  trans = G4Translate3D(fCenterPosition)*trans;
444  trans2 = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(pos2);
445  trans2 = G4Translate3D(fCenterPosition)*trans2;
446  } else {
449  }
450  G4double c[4];
451  colorMap->GetMapColor(xzcell[x][z], c);
452  att.SetColour(c[0], c[1], c[2]);//, c[3]);
453 
454  G4Box xzplate("xz", fSize[0]/fNSegment[0], thick,//fSize[1]/fNSegment[1]*0.001,
455  fSize[2]/fNSegment[2]);
456  G4Polyhedron * poly = xzplate.GetPolyhedron();
457  poly->Transform(trans);
458  poly->SetVisAttributes(&att);
459  pVisManager->Draw(*poly);
460 
461  G4Box xzplate2 = xzplate;
462  G4Polyhedron * poly2 = xzplate2.GetPolyhedron();
463  poly2->Transform(trans2);
464  poly2->SetVisAttributes(&att);
465  pVisManager->Draw(*poly2);
466 
467 
468  /*
469  G4double nodes[][3] =
470  {{-fSize[1]/fNSegment[1], 0., -fSize[2]/fNSegment[2]},
471  { fSize[1]/fNSegment[1], 0., -fSize[2]/fNSegment[2]},
472  { fSize[1]/fNSegment[1], 0., fSize[2]/fNSegment[2]},
473  {-fSize[1]/fNSegment[1], 0., fSize[2]/fNSegment[2]}};
474  G4int facets[][4] = {{1, 2, 3, 4}};
475  G4int facets2[][4] = {{4, 3, 2, 1}};
476 
477  G4Polyhedron poly, poly2;
478  poly.createPolyhedron(4, 1, nodes, facets);
479  poly.Transform(trans);
480  poly.SetVisAttributes(att);
481  pVisManager->Draw(poly);
482 
483  poly2.createPolyhedron(4, 1, nodes, facets2);
484  poly2.Transform(trans2);
485  poly2.SetVisAttributes(att);
486  pVisManager->Draw(poly2);
487  */
488  }
489  }
490  pVisManager->EndDraw();
491  }
492  }
493  colorMap->SetPSUnit(fDrawUnit);
494  colorMap->SetPSName(fDrawPSName);
495  colorMap->DrawColorChart();
496 }
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
void SetColour(const G4Colour &)
void SetPSName(G4String &psName)
G4ThreeVector fCenterPosition
G4ThreeVector GetReplicaPosition(G4int x, G4int y, G4int z)
static G4VVisManager * GetConcreteInstance()
Definition: G4Box.hh:64
void SetPSUnit(G4String &unit)
virtual void GetMapColor(G4double val, G4double color[4])=0
G4double fSize[3]
void SetForceSolid(G4bool)
HepPolyhedron & Transform(const G4Transform3D &t)
int G4int
Definition: G4Types.hh:78
Double_t y
void SetMinMax(G4double minVal, G4double maxVal)
G4bool IfFloatMinMax() const
G4double fDrawUnitValue
G4RotationMatrix * fRotationMatrix
virtual void EndDraw()=0
Double_t scale
virtual void BeginDraw(const G4Transform3D &objectTransformation=G4Transform3D())=0
void SetForceAuxEdgeVisible(G4bool)
void SetVisAttributes(const G4VisAttributes *)
Definition: G4Visible.cc:80
HepGeom::Rotate3D G4Rotate3D
HepGeom::Translate3D G4Translate3D
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4CSGSolid.cc:134
void GetXYZ(G4int index, G4int q[3]) const
static const G4double pos
G4String fDrawPSName
Transform3D inverse() const
Definition: Transform3D.cc:142
virtual void DrawColorChart(G4int nPoint=5)
Here is the call graph for this function:

◆ DrawColumn()

void G4ScoringBox::DrawColumn ( std::map< G4int, G4double *> *  map,
G4VScoreColorMap colorMap,
G4int  idxProj,
G4int  idxColumn 
)
virtual

Implements G4VScoringMesh.

Definition at line 521 of file G4ScoringBox.cc.

524 {
525  G4int iColumn[3] = {2, 0, 1};
526  if(idxColumn<0 || idxColumn>=fNSegment[iColumn[idxProj]])
527  {
528  G4cerr << "ERROR : Column number " << idxColumn
529  << " is out of scoring mesh [0," << fNSegment[iColumn[idxProj]]-1
530  << "]. Method ignored." << G4endl;
531  return;
532  }
534  if(pVisManager) {
535  pVisManager->BeginDraw();
536 
537  // cell vectors
538  std::vector<std::vector<std::vector<double> > > cell; // cell[X][Y][Z]
539  std::vector<double> ez;
540  for(int z = 0; z < fNSegment[2]; z++) ez.push_back(0.);
541  std::vector<std::vector<double> > eyz;
542  for(int y = 0; y < fNSegment[1]; y++) eyz.push_back(ez);
543  for(int x = 0; x < fNSegment[0]; x++) cell.push_back(eyz);
544 
545  std::vector<std::vector<double> > xycell; // xycell[X][Y]
546  std::vector<double> ey;
547  for(int y = 0; y < fNSegment[1]; y++) ey.push_back(0.);
548  for(int x = 0; x < fNSegment[0]; x++) xycell.push_back(ey);
549 
550  std::vector<std::vector<double> > yzcell; // yzcell[Y][Z]
551  for(int y = 0; y < fNSegment[1]; y++) yzcell.push_back(ez);
552 
553  std::vector<std::vector<double> > xzcell; // xzcell[X][Z]
554  for(int x = 0; x < fNSegment[0]; x++) xzcell.push_back(ez);
555 
556  // projections
557  G4int q[3];
558  std::map<G4int, G4double*>::iterator itr = map->begin();
559  for(; itr != map->end(); itr++) {
560  GetXYZ(itr->first, q);
561 
562  if(idxProj == 0 && q[2] == idxColumn) { // xy plane
563  xycell[q[0]][q[1]] += *(itr->second)/fDrawUnitValue;
564  }
565  if(idxProj == 1 && q[0] == idxColumn) { // yz plane
566  yzcell[q[1]][q[2]] += *(itr->second)/fDrawUnitValue;
567  }
568  if(idxProj == 2 && q[1] == idxColumn) { // zx plane
569  xzcell[q[0]][q[2]] += *(itr->second)/fDrawUnitValue;
570  }
571  }
572 
573  // search max. & min. values in each slice
574  G4double xymin = DBL_MAX, yzmin = DBL_MAX, xzmin = DBL_MAX;
575  G4double xymax = 0., yzmax = 0., xzmax = 0.;
576  for(int x = 0; x < fNSegment[0]; x++) {
577  for(int y = 0; y < fNSegment[1]; y++) {
578  if(xymin > xycell[x][y]) xymin = xycell[x][y];
579  if(xymax < xycell[x][y]) xymax = xycell[x][y];
580  }
581  for(int z = 0; z < fNSegment[2]; z++) {
582  if(xzmin > xzcell[x][z]) xzmin = xzcell[x][z];
583  if(xzmax < xzcell[x][z]) xzmax = xzcell[x][z];
584  }
585  }
586  for(int y = 0; y < fNSegment[1]; y++) {
587  for(int z = 0; z < fNSegment[2]; z++) {
588  if(yzmin > yzcell[y][z]) yzmin = yzcell[y][z];
589  if(yzmax < yzcell[y][z]) yzmax = yzcell[y][z];
590  }
591  }
592 
593 
594  G4VisAttributes att;
595  att.SetForceSolid(true);
596  att.SetForceAuxEdgeVisible(true);
597 
599  // xy plane
600  if(idxProj == 0) {
601  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(xymin,xymax); }
602  for(int x = 0; x < fNSegment[0]; x++) {
603  for(int y = 0; y < fNSegment[1]; y++) {
604  G4Box xyplate("xy", fSize[0]/fNSegment[0], fSize[1]/fNSegment[1],
605  fSize[2]/fNSegment[2]);
606 
607  G4ThreeVector pos(GetReplicaPosition(x, y, idxColumn));
608  G4Transform3D trans;
609  if(fRotationMatrix) {
610  trans = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(pos);
611  trans = G4Translate3D(fCenterPosition)*trans;
612  } else {
614  }
615  G4double c[4];
616  colorMap->GetMapColor(xycell[x][y], c);
617  att.SetColour(c[0], c[1], c[2]);
618 
619  G4Polyhedron * poly = xyplate.GetPolyhedron();
620  poly->Transform(trans);
621  poly->SetVisAttributes(att);
622  pVisManager->Draw(*poly);
623  }
624  }
625 
626  } else
627  // yz plane
628  if(idxProj == 1) {
629  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(yzmin,yzmax); }
630  for(int y = 0; y < fNSegment[1]; y++) {
631  for(int z = 0; z < fNSegment[2]; z++) {
632  G4Box yzplate("yz", fSize[0]/fNSegment[0], fSize[1]/fNSegment[1],
633  fSize[2]/fNSegment[2]);
634 
635  G4ThreeVector pos(GetReplicaPosition(idxColumn, y, z));
636  G4Transform3D trans;
637  if(fRotationMatrix) {
638  trans = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(pos);
639  trans = G4Translate3D(fCenterPosition)*trans;
640  } else {
642  }
643  G4double c[4];
644  colorMap->GetMapColor(yzcell[y][z], c);
645  att.SetColour(c[0], c[1], c[2]);//, c[3]);
646 
647  G4Polyhedron * poly = yzplate.GetPolyhedron();
648  poly->Transform(trans);
649  poly->SetVisAttributes(att);
650  pVisManager->Draw(*poly);
651  }
652  }
653  } else
654  // xz plane
655  if(idxProj == 2) {
656  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(xzmin,xzmax);}
657  for(int x = 0; x < fNSegment[0]; x++) {
658  for(int z = 0; z < fNSegment[2]; z++) {
659  G4Box xzplate("xz", fSize[0]/fNSegment[0], fSize[1]/fNSegment[1],
660  fSize[2]/fNSegment[2]);
661 
662  G4ThreeVector pos(GetReplicaPosition(x, idxColumn, z));
663  G4Transform3D trans;
664  if(fRotationMatrix) {
665  trans = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(pos);
666  trans = G4Translate3D(fCenterPosition)*trans;
667  } else {
669  }
670  G4double c[4];
671  colorMap->GetMapColor(xzcell[x][z], c);
672  att.SetColour(c[0], c[1], c[2]);//, c[3]);
673 
674  G4Polyhedron * poly = xzplate.GetPolyhedron();
675  poly->Transform(trans);
676  poly->SetVisAttributes(att);
677  pVisManager->Draw(*poly);
678  }
679  }
680  }
681  pVisManager->EndDraw();
682 
683  }
684 
685  colorMap->SetPSUnit(fDrawUnit);
686  colorMap->SetPSName(fDrawPSName);
687  colorMap->DrawColorChart();
688 }
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
void SetColour(const G4Colour &)
void SetPSName(G4String &psName)
G4ThreeVector fCenterPosition
G4ThreeVector GetReplicaPosition(G4int x, G4int y, G4int z)
static G4VVisManager * GetConcreteInstance()
Definition: G4Box.hh:64
void SetPSUnit(G4String &unit)
virtual void GetMapColor(G4double val, G4double color[4])=0
G4double fSize[3]
void SetForceSolid(G4bool)
HepPolyhedron & Transform(const G4Transform3D &t)
int G4int
Definition: G4Types.hh:78
Double_t y
void SetMinMax(G4double minVal, G4double maxVal)
G4bool IfFloatMinMax() const
G4double fDrawUnitValue
G4RotationMatrix * fRotationMatrix
virtual void EndDraw()=0
Double_t scale
virtual void BeginDraw(const G4Transform3D &objectTransformation=G4Transform3D())=0
void SetForceAuxEdgeVisible(G4bool)
void SetVisAttributes(const G4VisAttributes *)
Definition: G4Visible.cc:80
HepGeom::Rotate3D G4Rotate3D
HepGeom::Translate3D G4Translate3D
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
void GetXYZ(G4int index, G4int q[3]) const
static const G4double pos
G4String fDrawPSName
G4GLOB_DLL std::ostream G4cerr
Transform3D inverse() const
Definition: Transform3D.cc:142
virtual void DrawColorChart(G4int nPoint=5)
Here is the call graph for this function:

◆ GetIndex()

G4int G4ScoringBox::GetIndex ( G4int  x,
G4int  y,
G4int  z 
) const
private

Definition at line 517 of file G4ScoringBox.cc.

517  {
518  return x + y*fNSegment[0] + z*fNSegment[0]*fNSegment[1];
519 }
Double_t y
Here is the caller graph for this function:

◆ GetReplicaPosition()

G4ThreeVector G4ScoringBox::GetReplicaPosition ( G4int  x,
G4int  y,
G4int  z 
)
private

Definition at line 498 of file G4ScoringBox.cc.

498  {
499 
501  fSize[2]/fNSegment[2]);
502  G4ThreeVector pos(-fSize[0] + 2*(x+0.5)*width.x(),
503  -fSize[1] + 2*(y+0.5)*width.y(),
504  -fSize[2] + 2*(z+0.5)*width.z());
505 
506  return pos;
507 }
#define width
G4double fSize[3]
Double_t y
static const G4double pos
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetXYZ()

void G4ScoringBox::GetXYZ ( G4int  index,
G4int  q[3] 
) const
private

Definition at line 509 of file G4ScoringBox.cc.

509  {
510 
511  q[0] = index/(fNSegment[2]*fNSegment[1]);
512  q[1] = (index - q[0]*fNSegment[2]*fNSegment[1])/fNSegment[2];
513  q[2] = index - q[1]*fNSegment[2] - q[0]*fNSegment[2]*fNSegment[1];
514 
515 }
Int_t index
Here is the caller graph for this function:

◆ List()

void G4ScoringBox::List ( ) const
virtual

Reimplemented from G4VScoringMesh.

Definition at line 228 of file G4ScoringBox.cc.

228  {
229  G4cout << "G4ScoringBox : " << fWorldName << " --- Shape: Box mesh" << G4endl;
230  G4cout << " Size (x, y, z): ("
231  << fSize[0]/cm << ", "
232  << fSize[1]/cm << ", "
233  << fSize[2]/cm << ") [cm]"
234  << G4endl;
235 
237 }
static const double cm
Definition: G4SIunits.hh:118
G4double fSize[3]
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
virtual void List() const
Here is the call graph for this function:

◆ SetSegmentDirection()

void G4ScoringBox::SetSegmentDirection ( G4int  dir)
inline

Definition at line 60 of file G4ScoringBox.hh.

G4int fSegmentDirection
Definition: G4ScoringBox.hh:71
TDirectory * dir
Here is the call graph for this function:

◆ SetupGeometry()

void G4ScoringBox::SetupGeometry ( G4VPhysicalVolume fWorldPhys)
protectedvirtual

Implements G4VScoringMesh.

Definition at line 66 of file G4ScoringBox.cc.

66  {
67 
68  if(verboseLevel > 9) G4cout << "G4ScoringBox::SetupGeometry() ..." << G4endl;
69 
70  // World
71  G4VPhysicalVolume * scoringWorld = fWorldPhys;
72  G4LogicalVolume * worldLogical = scoringWorld->GetLogicalVolume();
73 
74  // Scoring Mesh
75  if(verboseLevel > 9) G4cout << fWorldName << G4endl;
76  G4String boxName = fWorldName;
77 
78  if(verboseLevel > 9)
79  G4cout << fSize[0] << ", " << fSize[1] << ", " << fSize[2] << G4endl;
80  G4VSolid * boxSolid = new G4Box(boxName+"0", fSize[0], fSize[1], fSize[2]);
81  G4LogicalVolume * boxLogical = new G4LogicalVolume(boxSolid, 0, boxName+"_0");
83  boxLogical, boxName+"0", worldLogical, false, 0);
84 
85  //G4double fsegment[3][3];
86  //G4int segOrder[3];
87  //GetSegmentOrder(fSegmentDirection, fNSegment, segOrder, fsegment);
88  //EAxis axis[3] = {kXAxis, kYAxis, kZAxis};
89 
90  G4String layerName[2] = {boxName + "_1", boxName + "_2"};
91  G4VSolid * layerSolid[2];
92  G4LogicalVolume * layerLogical[2];
93 
94  //-- fisrt nested layer (replicated to x direction)
95  if(verboseLevel > 9) G4cout << "layer 1 :" << G4endl;
96  layerSolid[0] = new G4Box(layerName[0],
97  fSize[0]/fNSegment[0],
98  fSize[1],
99  fSize[2]);
100  layerLogical[0] = new G4LogicalVolume(layerSolid[0], 0, layerName[0]);
101  if(fNSegment[0] > 1) {
102  if(verboseLevel > 9)
103  G4cout << "G4ScoringBox::Construct() : Replicate to x direction" << G4endl;
105  {
106  new G4PVReplica(layerName[0], layerLogical[0], boxLogical, kXAxis,
107  fNSegment[0], fSize[0]/fNSegment[0]*2.);
108  }
109  else
110  {
111  new G4PVDivision(layerName[0], layerLogical[0], boxLogical, kXAxis,
112  fNSegment[0], 0.);
113  }
114  } else if(fNSegment[0] == 1) {
115  if(verboseLevel > 9)
116  G4cout << "G4ScoringBox::Construct() : Placement" << G4endl;
117  new G4PVPlacement(0, G4ThreeVector(0.,0.,0.), layerLogical[0], layerName[0],
118  boxLogical, false, 0);
119  } else
120  G4cerr << "ERROR : G4ScoringBox::SetupGeometry() : invalid parameter ("
121  << fNSegment[0] << ") "
122  << "in placement of the first nested layer." << G4endl;
123 
124  if(verboseLevel > 9) {
125  G4cout << fSize[0]/fNSegment[0] << ", "
126  << fSize[1] << ", "
127  << fSize[2] << G4endl;
128  G4cout << layerName[0] << ": kXAxis, "
129  << fNSegment[0] << ", "
130  << 2.*fSize[0]/fNSegment[0] << G4endl;
131  }
132 
133  // second nested layer (replicated to y direction)
134  if(verboseLevel > 9) G4cout << "layer 2 :" << G4endl;
135  layerSolid[1] = new G4Box(layerName[1],
136  fSize[0]/fNSegment[0],
137  fSize[1]/fNSegment[1],
138  fSize[2]);
139  layerLogical[1] = new G4LogicalVolume(layerSolid[1], 0, layerName[1]);
140  if(fNSegment[1] > 1) {
141  if(verboseLevel > 9)
142  G4cout << "G4ScoringBox::Construct() : Replicate to y direction" << G4endl;
144  {
145  new G4PVReplica(layerName[1], layerLogical[1], layerLogical[0], kYAxis,
146  fNSegment[1], fSize[1]/fNSegment[1]*2.);
147  }
148  else
149  {
150  new G4PVDivision(layerName[1], layerLogical[1], layerLogical[0], kYAxis,
151  fNSegment[1], 0.);
152  }
153  } else if(fNSegment[1] == 1) {
154  if(verboseLevel > 9)
155  G4cout << "G4ScoringBox::Construct() : Placement" << G4endl;
156  new G4PVPlacement(0, G4ThreeVector(0.,0.,0.), layerLogical[1], layerName[1],
157  layerLogical[0], false, 0);
158  } else
159  G4cerr << "ERROR : G4ScoringBox::SetupGeometry() : invalid parameter ("
160  << fNSegment[1] << ") "
161  << "in placement of the second nested layer." << G4endl;
162 
163  if(verboseLevel > 9) {
164  G4cout << fSize[0]/fNSegment[0] << ", "
165  << fSize[1]/fNSegment[1] << ", "
166  << fSize[2] << G4endl;
167  G4cout << layerName[1] << ": kYAxis, "
168  << fNSegment[1] << ", "
169  << 2.*fSize[1]/fNSegment[1] << G4endl;
170  }
171 
172  // mesh elements (replicated to z direction)
173  if(verboseLevel > 9) G4cout << "mesh elements :" << G4endl;
174  G4String elementName = boxName +"_3";
175  G4VSolid * elementSolid = new G4Box(elementName,
176  fSize[0]/fNSegment[0],
177  fSize[1]/fNSegment[1],
178  fSize[2]/fNSegment[2]);
179  fMeshElementLogical = new G4LogicalVolume(elementSolid, 0, elementName);
180  if(fNSegment[2] > 1) {
181  if(verboseLevel > 9)
182  G4cout << "G4ScoringBox::Construct() : Replicate to z direction" << G4endl;
183 
185  {
186  new G4PVReplica(elementName, fMeshElementLogical, layerLogical[1], kZAxis,
187  fNSegment[2], 2.*fSize[2]/fNSegment[2]);
188  }
189  else
190  {
191  new G4PVDivision(elementName, fMeshElementLogical, layerLogical[1], kZAxis,
192  fNSegment[2], 0.);
193  }
194  } else if(fNSegment[2] == 1) {
195  if(verboseLevel > 9)
196  G4cout << "G4ScoringBox::Construct() : Placement" << G4endl;
198  elementName, layerLogical[1], false, 0);
199  } else
200  G4cerr << "ERROR : G4ScoringBox::SetupGeometry() : "
201  << "invalid parameter (" << fNSegment[2] << ") "
202  << "in mesh element placement." << G4endl;
203 
204  if(verboseLevel > 9) {
205  G4cout << fSize[0]/fNSegment[0] << ", "
206  << fSize[1]/fNSegment[1] << ", "
207  << fSize[2]/fNSegment[2] << G4endl;
208  G4cout << elementName << ": kZAxis, "
209  << fNSegment[2] << ", "
210  << 2.*fSize[2]/fNSegment[2] << G4endl;
211  }
212 
213 
214  // set the sensitive detector
216 
217 
218  // vis. attributes
219  G4VisAttributes * visatt = new G4VisAttributes(G4Colour(.5,.5,.5));
220  visatt->SetVisibility(false);
221  layerLogical[0]->SetVisAttributes(visatt);
222  layerLogical[1]->SetVisAttributes(visatt);
223  visatt->SetVisibility(true);
225 }
G4ThreeVector fCenterPosition
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Box.hh:64
G4MultiFunctionalDetector * fMFD
void SetVisibility(G4bool)
G4double fSize[3]
G4GLOB_DLL std::ostream G4cout
static G4int GetReplicaLevel()
G4RotationMatrix * fRotationMatrix
G4LogicalVolume * fMeshElementLogical
#define G4endl
Definition: G4ios.hh:61
G4LogicalVolume * GetLogicalVolume() const
void SetVisAttributes(const G4VisAttributes *pVA)
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
G4GLOB_DLL std::ostream G4cerr
Here is the call graph for this function:

Member Data Documentation

◆ fSegmentDirection

G4int G4ScoringBox::fSegmentDirection
private

Definition at line 71 of file G4ScoringBox.hh.


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