Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ScoringBox Class Reference

#include <G4ScoringBox.hh>

+ Inheritance diagram for G4ScoringBox:

Public Member Functions

 G4ScoringBox (G4String wName)
 
 ~G4ScoringBox () override=default
 
void List () const override
 
void Draw (RunScore *map, G4VScoreColorMap *colorMap, G4int axflg=111) override
 
void DrawColumn (RunScore *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn) override
 
void SetSegmentDirection (G4int dir)
 
- Public Member Functions inherited from G4VScoringMesh
 G4VScoringMesh (const G4String &wName)
 
virtual ~G4VScoringMesh ()=default
 
virtual void Construct (G4VPhysicalVolume *fWorldPhys)
 
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 Accumulate (G4THitsMap< G4StatDouble > *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 SetAngles (G4double, G4double)
 
G4double GetStartAngle () const
 
G4double GetAngleSpan () 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
 
G4VPrimitiveScorerGetPrimitiveScorer (const G4String &name)
 
void SetMeshElementLogical (G4LogicalVolume *val)
 
G4LogicalVolumeGetMeshElementLogical () const
 
void SetParallelWorldProcess (G4ParallelWorldProcess *proc)
 
G4ParallelWorldProcessGetParallelWorldProcess () const
 
void GeometryHasBeenDestroyed ()
 
void SetCopyNumberLevel (G4int val)
 
G4int GetCopyNumberLevel () const
 
G4bool LayeredMassFlg ()
 

Protected Member Functions

void SetupGeometry (G4VPhysicalVolume *fWorldPhys) override
 
- Protected Member Functions inherited from G4VScoringMesh

Additional Inherited Members

- Public Types inherited from G4VScoringMesh
enum class  MeshShape {
  box , cylinder , sphere , realWorldLogVol ,
  probe , undefined = -1
}
 
using EventScore = G4THitsMap<G4double>
 
using RunScore = G4THitsMap<G4StatDouble>
 
using MeshScoreMap = std::map<G4String, RunScore*>
 
- Protected Attributes inherited from G4VScoringMesh
G4String fWorldName
 
G4VPrimitiveScorerfCurrentPS
 
G4bool fConstructed
 
G4bool fActive
 
MeshShape fShape
 
G4double fSize [3]
 
G4double fAngle [2]
 
G4ThreeVector fCenterPosition
 
G4RotationMatrixfRotationMatrix
 
G4int fNSegment [3]
 
MeshScoreMap fMap
 
G4MultiFunctionalDetectorfMFD
 
G4int verboseLevel
 
G4bool sizeIsSet
 
G4bool nMeshIsSet
 
G4String fDrawUnit
 
G4double fDrawUnitValue
 
G4String fDrawPSName
 
G4String fDivisionAxisNames [3]
 
G4LogicalVolumefMeshElementLogical
 
G4ParallelWorldProcessfParallelWorldProcess
 
G4bool fGeometryHasBeenDestroyed
 
G4int copyNumberLevel
 
G4bool layeredMassFlg
 

Detailed Description

Definition at line 42 of file G4ScoringBox.hh.

Constructor & Destructor Documentation

◆ G4ScoringBox()

G4ScoringBox::G4ScoringBox ( G4String wName)

Definition at line 53 of file G4ScoringBox.cc.

54 : G4VScoringMesh(wName)
55 , fSegmentDirection(-1)
56{
58 fDivisionAxisNames[0] = "X";
59 fDivisionAxisNames[1] = "Y";
60 fDivisionAxisNames[2] = "Z";
61}
G4String fDivisionAxisNames[3]
G4VScoringMesh(const G4String &wName)

◆ ~G4ScoringBox()

G4ScoringBox::~G4ScoringBox ( )
overridedefault

Member Function Documentation

◆ Draw()

void G4ScoringBox::Draw ( RunScore * map,
G4VScoreColorMap * colorMap,
G4int axflg = 111 )
overridevirtual

Implements G4VScoringMesh.

Definition at line 239 of file G4ScoringBox.cc.

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

◆ DrawColumn()

void G4ScoringBox::DrawColumn ( RunScore * map,
G4VScoreColorMap * colorMap,
G4int idxProj,
G4int idxColumn )
overridevirtual

Implements G4VScoringMesh.

Definition at line 504 of file G4ScoringBox.cc.

506{
507 G4int iColumn[3] = { 2, 0, 1 };
508 if(idxColumn < 0 || idxColumn >= fNSegment[iColumn[idxProj]])
509 {
510 G4cerr << "ERROR : Column number " << idxColumn
511 << " is out of scoring mesh [0," << fNSegment[iColumn[idxProj]] - 1
512 << "]. Method ignored." << G4endl;
513 return;
514 }
516 if(pVisManager != nullptr)
517 {
518 pVisManager->BeginDraw();
519
520 // cell vectors
521 std::vector<std::vector<std::vector<double>>> cell; // cell[X][Y][Z]
522 std::vector<double> ez;
523 for(int z = 0; z < fNSegment[2]; z++)
524 ez.push_back(0.);
525 std::vector<std::vector<double>> eyz;
526 for(int y = 0; y < fNSegment[1]; y++)
527 eyz.push_back(ez);
528 for(int x = 0; x < fNSegment[0]; x++)
529 cell.push_back(eyz);
530
531 std::vector<std::vector<double>> xycell; // xycell[X][Y]
532 std::vector<double> ey;
533 for(int y = 0; y < fNSegment[1]; y++)
534 ey.push_back(0.);
535 for(int x = 0; x < fNSegment[0]; x++)
536 xycell.push_back(ey);
537
538 std::vector<std::vector<double>> yzcell; // yzcell[Y][Z]
539 for(int y = 0; y < fNSegment[1]; y++)
540 yzcell.push_back(ez);
541
542 std::vector<std::vector<double>> xzcell; // xzcell[X][Z]
543 for(int x = 0; x < fNSegment[0]; x++)
544 xzcell.push_back(ez);
545
546 // projections
547 G4int q[3];
548 auto itr = map->GetMap()->begin();
549 for(; itr != map->GetMap()->end(); itr++)
550 {
551 GetXYZ(itr->first, q);
552
553 if(idxProj == 0 && q[2] == idxColumn)
554 { // xy plane
555 xycell[q[0]][q[1]] += (itr->second->sum_wx()) / fDrawUnitValue;
556 }
557 if(idxProj == 1 && q[0] == idxColumn)
558 { // yz plane
559 yzcell[q[1]][q[2]] += (itr->second->sum_wx()) / fDrawUnitValue;
560 }
561 if(idxProj == 2 && q[1] == idxColumn)
562 { // zx plane
563 xzcell[q[0]][q[2]] += (itr->second->sum_wx()) / fDrawUnitValue;
564 }
565 }
566
567 // search max. & min. values in each slice
568 G4double xymin = DBL_MAX, yzmin = DBL_MAX, xzmin = DBL_MAX;
569 G4double xymax = 0., yzmax = 0., xzmax = 0.;
570 for(int x = 0; x < fNSegment[0]; x++)
571 {
572 for(int y = 0; y < fNSegment[1]; y++)
573 {
574 if(xymin > xycell[x][y])
575 xymin = xycell[x][y];
576 if(xymax < xycell[x][y])
577 xymax = xycell[x][y];
578 }
579 for(int z = 0; z < fNSegment[2]; z++)
580 {
581 if(xzmin > xzcell[x][z])
582 xzmin = xzcell[x][z];
583 if(xzmax < xzcell[x][z])
584 xzmax = xzcell[x][z];
585 }
586 }
587 for(int y = 0; y < fNSegment[1]; y++)
588 {
589 for(int z = 0; z < fNSegment[2]; z++)
590 {
591 if(yzmin > yzcell[y][z])
592 yzmin = yzcell[y][z];
593 if(yzmax < yzcell[y][z])
594 yzmax = yzcell[y][z];
595 }
596 }
597
598 G4VisAttributes att;
599 att.SetForceSolid(true);
600 att.SetForceAuxEdgeVisible(true);
601
602 G4Scale3D scale;
603 // xy plane
604 if(idxProj == 0)
605 {
606 if(colorMap->IfFloatMinMax())
607 {
608 colorMap->SetMinMax(xymin, xymax);
609 }
610 for(int x = 0; x < fNSegment[0]; x++)
611 {
612 for(int y = 0; y < fNSegment[1]; y++)
613 {
614 G4Box xyplate("xy", fSize[0] / fNSegment[0], fSize[1] / fNSegment[1],
615 fSize[2] / fNSegment[2]);
616
617 G4ThreeVector pos(GetReplicaPosition(x, y, idxColumn));
618 G4Transform3D trans;
619 if(fRotationMatrix != nullptr)
620 {
621 trans = G4Rotate3D(*fRotationMatrix).inverse() * G4Translate3D(pos);
622 trans = G4Translate3D(fCenterPosition) * trans;
623 }
624 else
625 {
627 }
628 G4double c[4];
629 colorMap->GetMapColor(xycell[x][y], c);
630 att.SetColour(c[0], c[1], c[2]);
631
632 G4Polyhedron* poly = xyplate.GetPolyhedron();
633 poly->Transform(trans);
634 poly->SetVisAttributes(att);
635 pVisManager->Draw(*poly);
636 }
637 }
638 }
639 else
640 // yz plane
641 if(idxProj == 1)
642 {
643 if(colorMap->IfFloatMinMax())
644 {
645 colorMap->SetMinMax(yzmin, yzmax);
646 }
647 for(int y = 0; y < fNSegment[1]; y++)
648 {
649 for(int z = 0; z < fNSegment[2]; z++)
650 {
651 G4Box yzplate("yz", fSize[0] / fNSegment[0], fSize[1] / fNSegment[1],
652 fSize[2] / fNSegment[2]);
653
654 G4ThreeVector pos(GetReplicaPosition(idxColumn, y, z));
655 G4Transform3D trans;
656 if(fRotationMatrix != nullptr)
657 {
658 trans = G4Rotate3D(*fRotationMatrix).inverse() * G4Translate3D(pos);
659 trans = G4Translate3D(fCenterPosition) * trans;
660 }
661 else
662 {
664 }
665 G4double c[4];
666 colorMap->GetMapColor(yzcell[y][z], c);
667 att.SetColour(c[0], c[1], c[2]); //, c[3]);
668
669 G4Polyhedron* poly = yzplate.GetPolyhedron();
670 poly->Transform(trans);
671 poly->SetVisAttributes(att);
672 pVisManager->Draw(*poly);
673 }
674 }
675 }
676 else
677 // xz plane
678 if(idxProj == 2)
679 {
680 if(colorMap->IfFloatMinMax())
681 {
682 colorMap->SetMinMax(xzmin, xzmax);
683 }
684 for(int x = 0; x < fNSegment[0]; x++)
685 {
686 for(int z = 0; z < fNSegment[2]; z++)
687 {
688 G4Box xzplate("xz", fSize[0] / fNSegment[0], fSize[1] / fNSegment[1],
689 fSize[2] / fNSegment[2]);
690
691 G4ThreeVector pos(GetReplicaPosition(x, idxColumn, z));
692 G4Transform3D trans;
693 if(fRotationMatrix != nullptr)
694 {
695 trans = G4Rotate3D(*fRotationMatrix).inverse() * G4Translate3D(pos);
696 trans = G4Translate3D(fCenterPosition) * trans;
697 }
698 else
699 {
701 }
702 G4double c[4];
703 colorMap->GetMapColor(xzcell[x][z], c);
704 att.SetColour(c[0], c[1], c[2]); //, c[3]);
705
706 G4Polyhedron* poly = xzplate.GetPolyhedron();
707 poly->Transform(trans);
708 poly->SetVisAttributes(att);
709 pVisManager->Draw(*poly);
710 }
711 }
712 }
713 pVisManager->EndDraw();
714 }
715
716 colorMap->SetPSUnit(fDrawUnit);
717 colorMap->SetPSName(fDrawPSName);
718 colorMap->DrawColorChart();
719}
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition G4ios.hh:67

◆ List()

void G4ScoringBox::List ( ) const
overridevirtual

Reimplemented from G4VScoringMesh.

Definition at line 230 of file G4ScoringBox.cc.

231{
232 G4cout << "G4ScoringBox : " << fWorldName << " --- Shape: Box mesh" << G4endl;
233 G4cout << " Size (x, y, z): (" << fSize[0] / cm << ", " << fSize[1] / cm
234 << ", " << fSize[2] / cm << ") [cm]" << G4endl;
235
237}
G4GLOB_DLL std::ostream G4cout
virtual void List() const

◆ SetSegmentDirection()

void G4ScoringBox::SetSegmentDirection ( G4int dir)
inline

Definition at line 55 of file G4ScoringBox.hh.

55{ fSegmentDirection = dir; }

◆ SetupGeometry()

void G4ScoringBox::SetupGeometry ( G4VPhysicalVolume * fWorldPhys)
overrideprotectedvirtual

Implements G4VScoringMesh.

Definition at line 63 of file G4ScoringBox.cc.

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

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