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

#include <G4ScoringCylinder.hh>

+ Inheritance diagram for G4ScoringCylinder:

Public Member Functions

 G4ScoringCylinder (const G4String &wName)
 
 ~G4ScoringCylinder () 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 SetRMin (G4double rMin)
 
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 ()=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
 

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 G4ScoringCylinder.hh.

Constructor & Destructor Documentation

◆ G4ScoringCylinder()

G4ScoringCylinder::G4ScoringCylinder ( const G4String & wName)

Definition at line 57 of file G4ScoringCylinder.cc.

58 : G4VScoringMesh(wName)
59{
61
62 fDivisionAxisNames[0] = "Z";
63 fDivisionAxisNames[1] = "PHI";
64 fDivisionAxisNames[2] = "R";
65}
G4String fDivisionAxisNames[3]
G4VScoringMesh(const G4String &wName)

◆ ~G4ScoringCylinder()

G4ScoringCylinder::~G4ScoringCylinder ( )
overridedefault

Member Function Documentation

◆ Draw()

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

Implements G4VScoringMesh.

Definition at line 253 of file G4ScoringCylinder.cc.

255{
256 G4VVisManager* pVisManager = G4VVisManager::GetConcreteInstance();
257 if(pVisManager != nullptr)
258 {
259 // cell vectors
260 std::vector<G4double> ephi;
261 for(G4int phi = 0; phi < fNSegment[IPHI]; phi++)
262 ephi.push_back(0.);
263 //-
264 std::vector<std::vector<G4double>> zphicell; // zphicell[Z][PHI]
265 for(G4int z = 0; z < fNSegment[IZ]; z++)
266 zphicell.push_back(ephi);
267 //-
268 std::vector<std::vector<G4double>> rphicell; // rphicell[R][PHI]
269 for(G4int r = 0; r < fNSegment[IR]; r++)
270 rphicell.push_back(ephi);
271
272 // projections
273 G4int q[3];
274 auto itr = map->GetMap()->begin();
275 for(; itr != map->GetMap()->end(); itr++)
276 {
277 if(itr->first < 0)
278 {
279 G4cout << itr->first << G4endl;
280 continue;
281 }
282 GetRZPhi(itr->first, q);
283
284 zphicell[q[IZ]][q[IPHI]] += (itr->second->sum_wx()) / fDrawUnitValue;
285 rphicell[q[IR]][q[IPHI]] += (itr->second->sum_wx()) / fDrawUnitValue;
286 }
287
288 // search min./max. values
289 G4double zphimin = DBL_MAX, rphimin = DBL_MAX;
290 G4double zphimax = 0., rphimax = 0.;
291 for(G4int iphi = 0; iphi < fNSegment[IPHI]; iphi++)
292 {
293 for(G4int iz = 0; iz < fNSegment[IZ]; iz++)
294 {
295 if(zphimin > zphicell[iz][iphi])
296 zphimin = zphicell[iz][iphi];
297 if(zphimax < zphicell[iz][iphi])
298 zphimax = zphicell[iz][iphi];
299 }
300 for(G4int ir = 0; ir < fNSegment[IR]; ir++)
301 {
302 if(rphimin > rphicell[ir][iphi])
303 rphimin = rphicell[ir][iphi];
304 if(rphimax < rphicell[ir][iphi])
305 rphimax = rphicell[ir][iphi];
306 }
307 }
308
309 G4VisAttributes att;
310 att.SetForceSolid(true);
311 att.SetForceAuxEdgeVisible(true);
312
313 G4Scale3D scale;
314 if(axflg / 100 == 1)
315 {
316 // rz plane
317 }
318 axflg = axflg % 100;
319 if(axflg / 10 == 1)
320 {
321 pVisManager->BeginDraw();
322
323 // z-phi plane
324 if(colorMap->IfFloatMinMax())
325 {
326 colorMap->SetMinMax(zphimin, zphimax);
327 }
328
329 G4double zhalf = fSize[2] / fNSegment[IZ];
330 for(G4int phi = 0; phi < fNSegment[IPHI]; phi++)
331 {
332 for(G4int z = 0; z < fNSegment[IZ]; z++)
333 {
334 //-
335 G4double angle = fAngle[0] + fAngle[1] / fNSegment[IPHI] * phi;
336 G4double dphi = fAngle[1] / fNSegment[IPHI];
337 G4Tubs cylinder(
338 "z-phi", // name
339 fSize[1] * 0.99, fSize[1], // inner radius, outer radius
340 zhalf, // half length in z
341 angle, dphi * 0.99999); // starting phi angle, delta angle
342 //-
343 G4ThreeVector zpos(
344 0., 0., -fSize[2] + fSize[2] / fNSegment[IZ] * (1 + 2. * z));
345 G4Transform3D trans;
346 if(fRotationMatrix != nullptr)
347 {
348 trans =
349 G4Rotate3D(*fRotationMatrix).inverse() * G4Translate3D(zpos);
350 trans = G4Translate3D(fCenterPosition) * trans;
351 }
352 else
353 {
355 }
356 G4double c[4];
357 colorMap->GetMapColor(zphicell[z][phi], c);
358 att.SetColour(c[0], c[1], c[2]); //, c[3]);
359 //-
360 G4Polyhedron* poly = cylinder.GetPolyhedron();
361 poly->Transform(trans);
362 poly->SetVisAttributes(att);
363 pVisManager->Draw(*poly);
364 }
365 }
366 pVisManager->EndDraw();
367 }
368 axflg = axflg % 10;
369 if(axflg == 1)
370 {
371 pVisManager->BeginDraw();
372
373 // r-phi plane
374 if(colorMap->IfFloatMinMax())
375 {
376 colorMap->SetMinMax(rphimin, rphimax);
377 }
378
379 G4double rsize = (fSize[1] - fSize[0]) / fNSegment[IR];
380 for(G4int phi = 0; phi < fNSegment[IPHI]; phi++)
381 {
382 for(G4int r = 0; r < fNSegment[IR]; r++)
383 {
384 G4double rs[2] = { fSize[0] + rsize * r, fSize[0] + rsize * (r + 1) };
385 G4double angle = fAngle[0] + fAngle[1] / fNSegment[IPHI] * phi;
386 G4double dphi = fAngle[1] / fNSegment[IPHI];
387 G4Tubs cylindern("z-phi", rs[0], rs[1], 0.001, angle, dphi * 0.99999);
388 G4Tubs cylinderp = cylindern;
389
390 G4ThreeVector zposn(0., 0., -fSize[2]);
391 G4ThreeVector zposp(0., 0., fSize[2]);
392 G4Transform3D transn, transp;
393 if(fRotationMatrix != nullptr)
394 {
395 transn =
396 G4Rotate3D(*fRotationMatrix).inverse() * G4Translate3D(zposn);
397 transn = G4Translate3D(fCenterPosition) * transn;
398 transp =
399 G4Rotate3D(*fRotationMatrix).inverse() * G4Translate3D(zposp);
400 transp = G4Translate3D(fCenterPosition) * transp;
401 }
402 else
403 {
406 }
407 G4double c[4];
408 colorMap->GetMapColor(rphicell[r][phi], c);
409 att.SetColour(c[0], c[1], c[2]); //, c[3]);
410
411 G4Polyhedron* polyn = cylindern.GetPolyhedron();
412 polyn->Transform(transn);
413 polyn->SetVisAttributes(att);
414 pVisManager->Draw(*polyn);
415
416 G4Polyhedron* polyp = cylinderp.GetPolyhedron();
417 polyp->Transform(transp);
418 polyp->SetVisAttributes(att);
419 pVisManager->Draw(*polyp);
420 }
421 }
422
423 pVisManager->EndDraw();
424 }
425
426 colorMap->SetPSUnit(fDrawUnit);
427 colorMap->SetPSName(fDrawPSName);
428 colorMap->DrawColorChart();
429 }
430}
CLHEP::Hep3Vector G4ThreeVector
HepGeom::Transform3D G4Transform3D
HepGeom::Translate3D G4Translate3D
HepGeom::Scale3D G4Scale3D
HepGeom::Rotate3D G4Rotate3D
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4Polyhedron * GetPolyhedron() const override
void GetRZPhi(G4int index, G4int q[3]) const
virtual void DrawColorChart(G4int nPoint=5)
G4bool IfFloatMinMax() const
void SetPSName(const G4String &psName)
void SetMinMax(G4double minVal, G4double maxVal)
virtual void GetMapColor(G4double val, G4double color[4])=0
void SetPSUnit(const G4String &unit)
G4RotationMatrix * fRotationMatrix
G4double fDrawUnitValue
G4double fAngle[2]
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 G4ScoringCylinder::DrawColumn ( RunScore * map,
G4VScoreColorMap * colorMap,
G4int idxProj,
G4int idxColumn )
overridevirtual

Implements G4VScoringMesh.

Definition at line 432 of file G4ScoringCylinder.cc.

434{
435 G4int projAxis = 0;
436 switch(idxProj)
437 {
438 case 0:
439 projAxis = IR;
440 break;
441 case 1:
442 projAxis = IZ;
443 break;
444 case 2:
445 projAxis = IPHI;
446 break;
447 }
448
449 if(idxColumn < 0 || idxColumn >= fNSegment[projAxis])
450 {
451 G4cerr << "Warning : Column number " << idxColumn
452 << " is out of scoring mesh [0," << fNSegment[projAxis] - 1
453 << "]. Method ignored." << G4endl;
454 return;
455 }
456 G4VVisManager* pVisManager = G4VVisManager::GetConcreteInstance();
457 if(pVisManager != nullptr)
458 {
459 // cell vectors
460 std::vector<std::vector<std::vector<G4double>>> cell; // cell[R][Z][PHI]
461 std::vector<G4double> ephi;
462 for(G4int phi = 0; phi < fNSegment[IPHI]; phi++)
463 ephi.push_back(0.);
464 std::vector<std::vector<G4double>> ezphi;
465 for(G4int z = 0; z < fNSegment[IZ]; z++)
466 ezphi.push_back(ephi);
467 for(G4int r = 0; r < fNSegment[IR]; r++)
468 cell.push_back(ezphi);
469
470 std::vector<std::vector<G4double>> rzcell; // rzcell[R][Z]
471 std::vector<G4double> ez;
472 for(G4int z = 0; z < fNSegment[IZ]; z++)
473 ez.push_back(0.);
474 for(G4int r = 0; r < fNSegment[IR]; r++)
475 rzcell.push_back(ez);
476
477 std::vector<std::vector<G4double>> zphicell; // zphicell[Z][PHI]
478 for(G4int z = 0; z < fNSegment[IZ]; z++)
479 zphicell.push_back(ephi);
480
481 std::vector<std::vector<G4double>> rphicell; // rphicell[R][PHI]
482 for(G4int r = 0; r < fNSegment[IR]; r++)
483 rphicell.push_back(ephi);
484
485 // projections
486 G4int q[3];
487 auto itr = map->GetMap()->begin();
488 for(; itr != map->GetMap()->end(); itr++)
489 {
490 if(itr->first < 0)
491 {
492 G4cout << itr->first << G4endl;
493 continue;
494 }
495 GetRZPhi(itr->first, q);
496
497 if(projAxis == IR && q[IR] == idxColumn)
498 { // zphi plane
499 zphicell[q[IZ]][q[IPHI]] += (itr->second->sum_wx()) / fDrawUnitValue;
500 }
501 if(projAxis == IZ && q[IZ] == idxColumn)
502 { // rphi plane
503 rphicell[q[IR]][q[IPHI]] += (itr->second->sum_wx()) / fDrawUnitValue;
504 }
505 if(projAxis == IPHI && q[IPHI] == idxColumn)
506 { // rz plane
507 rzcell[q[IR]][q[IZ]] += (itr->second->sum_wx()) / fDrawUnitValue;
508 }
509 }
510
511 // search min./max. values
512 G4double rzmin = DBL_MAX, zphimin = DBL_MAX, rphimin = DBL_MAX;
513 G4double rzmax = 0., zphimax = 0., rphimax = 0.;
514 for(G4int r = 0; r < fNSegment[IR]; r++)
515 {
516 for(G4int phi = 0; phi < fNSegment[IPHI]; phi++)
517 {
518 if(rphimin > rphicell[r][phi])
519 rphimin = rphicell[r][phi];
520 if(rphimax < rphicell[r][phi])
521 rphimax = rphicell[r][phi];
522 }
523 for(G4int z = 0; z < fNSegment[IZ]; z++)
524 {
525 if(rzmin > rzcell[r][z])
526 rzmin = rzcell[r][z];
527 if(rzmax < rzcell[r][z])
528 rzmax = rzcell[r][z];
529 }
530 }
531 for(G4int z = 0; z < fNSegment[IZ]; z++)
532 {
533 for(G4int phi = 0; phi < fNSegment[IPHI]; phi++)
534 {
535 if(zphimin > zphicell[z][phi])
536 zphimin = zphicell[z][phi];
537 if(zphimax < zphicell[z][phi])
538 zphimax = zphicell[z][phi];
539 }
540 }
541
542 G4VisAttributes att;
543 att.SetForceSolid(true);
544 att.SetForceAuxEdgeVisible(true);
545
546 pVisManager->BeginDraw();
547
548 G4Scale3D scale;
549 // z-phi plane
550 if(projAxis == IR)
551 {
552 if(colorMap->IfFloatMinMax())
553 {
554 colorMap->SetMinMax(zphimin, zphimax);
555 }
556
557 G4double zhalf = fSize[2] / fNSegment[IZ];
558 G4double rsize[2] = {
559 fSize[0] + (fSize[1] - fSize[0]) / fNSegment[IR] * idxColumn,
560 fSize[0] + (fSize[1] - fSize[0]) / fNSegment[IR] * (idxColumn + 1)
561 };
562 for(int phi = 0; phi < fNSegment[IPHI]; phi++)
563 {
564 for(int z = 0; z < fNSegment[IZ]; z++)
565 {
566 G4double angle = fAngle[0] + fAngle[1] / fNSegment[IPHI] * phi;
567 G4double dphi = fAngle[1] / fNSegment[IPHI];
568 G4Tubs cylinder("z-phi", rsize[0], rsize[1], zhalf, angle,
569 dphi * 0.99999);
570
571 G4ThreeVector zpos(
572 0., 0., -fSize[2] + fSize[2] / fNSegment[IZ] * (1 + 2. * z));
573 G4Transform3D trans;
574 if(fRotationMatrix != nullptr)
575 {
576 trans =
577 G4Rotate3D(*fRotationMatrix).inverse() * G4Translate3D(zpos);
578 trans = G4Translate3D(fCenterPosition) * trans;
579 }
580 else
581 {
583 }
584 G4double c[4];
585 colorMap->GetMapColor(zphicell[z][phi], c);
586 att.SetColour(c[0], c[1], c[2]); //, c[3]);
587
588 G4Polyhedron* poly = cylinder.GetPolyhedron();
589 poly->Transform(trans);
590 poly->SetVisAttributes(att);
591 pVisManager->Draw(*poly);
592 }
593 }
594
595 // r-phi plane
596 }
597 else if(projAxis == IZ)
598 {
599 if(colorMap->IfFloatMinMax())
600 {
601 colorMap->SetMinMax(rphimin, rphimax);
602 }
603
604 G4double rsize = (fSize[1] - fSize[0]) / fNSegment[IR];
605 for(G4int phi = 0; phi < fNSegment[IPHI]; phi++)
606 {
607 for(G4int r = 0; r < fNSegment[IR]; r++)
608 {
609 G4double rs[2] = { fSize[0] + rsize * r, fSize[0] + rsize * (r + 1) };
610 G4double angle = fAngle[0] + fAngle[1] / fNSegment[IPHI] * phi;
611 G4double dz = fSize[2] / fNSegment[IZ];
612 G4double dphi = fAngle[1] / fNSegment[IPHI];
613 G4Tubs cylinder("r-phi", rs[0], rs[1], dz, angle, dphi * 0.99999);
614 G4ThreeVector zpos(
615 0., 0., -fSize[2] + fSize[2] / fNSegment[IZ] * (idxColumn * 2 + 1));
616 G4Transform3D trans;
617 if(fRotationMatrix != nullptr)
618 {
619 trans =
620 G4Rotate3D(*fRotationMatrix).inverse() * G4Translate3D(zpos);
621 trans = G4Translate3D(fCenterPosition) * trans;
622 }
623 else
624 {
626 }
627 G4double c[4];
628 colorMap->GetMapColor(rphicell[r][phi], c);
629 att.SetColour(c[0], c[1], c[2]); //, c[3]);
630
631 G4Polyhedron* poly = cylinder.GetPolyhedron();
632 poly->Transform(trans);
633 poly->SetVisAttributes(att);
634 pVisManager->Draw(*poly);
635 }
636 }
637
638 // r-z plane
639 }
640 else if(projAxis == IPHI)
641 {
642 if(colorMap->IfFloatMinMax())
643 {
644 colorMap->SetMinMax(rzmin, rzmax);
645 }
646
647 G4double rsize = (fSize[1] - fSize[0]) / fNSegment[IR];
648 G4double zhalf = fSize[2] / fNSegment[IZ];
649 G4double angle = fAngle[0] + fAngle[1] / fNSegment[IPHI] * idxColumn;
650 G4double dphi = fAngle[1] / fNSegment[IPHI];
651 for(G4int z = 0; z < fNSegment[IZ]; z++)
652 {
653 for(G4int r = 0; r < fNSegment[IR]; r++)
654 {
655 G4double rs[2] = { fSize[0] + rsize * r, fSize[0] + rsize * (r + 1) };
656 G4Tubs cylinder("z-phi", rs[0], rs[1], zhalf, angle, dphi);
657
658 G4ThreeVector zpos(
659 0., 0., -fSize[2] + fSize[2] / fNSegment[IZ] * (2. * z + 1));
660 G4Transform3D trans;
661 if(fRotationMatrix != nullptr)
662 {
663 trans =
664 G4Rotate3D(*fRotationMatrix).inverse() * G4Translate3D(zpos);
665 trans = G4Translate3D(fCenterPosition) * trans;
666 }
667 else
668 {
670 }
671 G4double c[4];
672 colorMap->GetMapColor(rzcell[r][z], c);
673 att.SetColour(c[0], c[1], c[2]); //, c[3]);
674
675 G4Polyhedron* poly = cylinder.GetPolyhedron();
676 poly->Transform(trans);
677 poly->SetVisAttributes(att);
678 pVisManager->Draw(*poly);
679 }
680 }
681 }
682 pVisManager->EndDraw();
683 }
684
685 colorMap->SetPSUnit(fDrawUnit);
686 colorMap->SetPSName(fDrawPSName);
687 colorMap->DrawColorChart();
688}
G4GLOB_DLL std::ostream G4cerr

◆ GetRZPhi()

void G4ScoringCylinder::GetRZPhi ( G4int index,
G4int q[3] ) const

Definition at line 690 of file G4ScoringCylinder.cc.

691{
692 // index = k + j * k-size + i * jk-plane-size
693
694 // nested : z -> phi -> r
695 G4int i = IZ;
696 G4int j = IPHI;
697 G4int k = IR;
698 G4int jk = fNSegment[j] * fNSegment[k];
699 q[i] = index / jk;
700 q[j] = (index - q[i] * jk) / fNSegment[k];
701 q[k] = index - q[j] * fNSegment[k] - q[i] * jk;
702}

Referenced by Draw(), and DrawColumn().

◆ List()

void G4ScoringCylinder::List ( ) const
overridevirtual

Reimplemented from G4VScoringMesh.

Definition at line 240 of file G4ScoringCylinder.cc.

241{
242 G4cout << "G4ScoringCylinder : " << fWorldName
243 << " --- Shape: Cylindrical mesh" << G4endl;
244
245 G4cout << " Size (Rmin, Rmax, Dz): (" << fSize[0] / cm << ", "
246 << fSize[1] / cm << ", " << fSize[2] / cm << ") [cm]" << G4endl;
247 G4cout << " Angle (start, span): (" << fAngle[0] / deg << ", "
248 << fAngle[1] / deg << ") [deg]" << G4endl;
249
251}
virtual void List() const

◆ RegisterPrimitives()

void G4ScoringCylinder::RegisterPrimitives ( std::vector< G4VPrimitiveScorer * > & vps)

◆ SetRMax()

void G4ScoringCylinder::SetRMax ( G4double rMax)
inline

Definition at line 56 of file G4ScoringCylinder.hh.

56{ fSize[1] = rMax; }

◆ SetRMin()

void G4ScoringCylinder::SetRMin ( G4double rMin)
inline

Definition at line 55 of file G4ScoringCylinder.hh.

55{ fSize[0] = rMin; }

◆ SetupGeometry()

void G4ScoringCylinder::SetupGeometry ( G4VPhysicalVolume * fWorldPhys)
overrideprotectedvirtual

Implements G4VScoringMesh.

Definition at line 67 of file G4ScoringCylinder.cc.

68{
69 if(verboseLevel > 9)
70 G4cout << "G4ScoringCylinder::SetupGeometry() ..." << G4endl;
71
72 // World
73 G4VPhysicalVolume* scoringWorld = fWorldPhys;
74 G4LogicalVolume* worldLogical = scoringWorld->GetLogicalVolume();
75
76 // Scoring Mesh
77 if(verboseLevel > 9)
79 G4String tubsName = fWorldName + "_mesh";
80
81 if(verboseLevel > 9)
82 {
83 G4cout << "R min, R max., Dz =: " << fSize[0] << ", " << fSize[1]
84 << ", " << fSize[2] << G4endl;
85 }
86 G4VSolid* tubsSolid = new G4Tubs(tubsName + "0", // name
87 fSize[0], // R min
88 fSize[1], // R max
89 fSize[2], // Dz
90 fAngle[0], // starting phi
91 fAngle[1]); // segment phi
92 auto tubsLogical = new G4LogicalVolume(tubsSolid, nullptr, tubsName);
93 new G4PVPlacement(fRotationMatrix, fCenterPosition, tubsLogical,
94 tubsName + "0", worldLogical, false, 0);
95
96 if(verboseLevel > 9)
97 G4cout << " # of segments : r, phi, z =: " << fNSegment[IR] << ", "
98 << fNSegment[IPHI] << ", " << fNSegment[IZ] << G4endl;
99
100 G4String layerName[2] = { tubsName + "1", tubsName + "2" };
101 G4VSolid* layerSolid[2];
102 G4LogicalVolume* layerLogical[2];
103
104 //-- fisrt nested layer (replicated along z direction)
105 if(verboseLevel > 9)
106 G4cout << "layer 1 :" << G4endl;
107 layerSolid[0] = new G4Tubs(layerName[0], // name
108 fSize[0], // inner radius
109 fSize[1], // outer radius
110 fSize[2] / fNSegment[IZ], // half len. in z
111 fAngle[0], // starting phi angle
112 fAngle[1]); // delta angle of the segment
113 layerLogical[0] = new G4LogicalVolume(layerSolid[0], nullptr, layerName[0]);
114 if(fNSegment[IZ] > 1)
115 {
116 if(verboseLevel > 9)
117 G4cout << "G4ScoringCylinder::Construct() : Replicate along z direction"
118 << G4endl;
120 {
121 if(verboseLevel > 9)
122 G4cout << "G4ScoringCylinder::Construct() : Replica" << G4endl;
123 new G4PVReplica(layerName[0], layerLogical[0], tubsLogical, kZAxis,
124 fNSegment[IZ], 2. * fSize[2] / fNSegment[IZ]);
125 }
126 else
127 {
128 if(verboseLevel > 9)
129 G4cout << "G4ScoringCylinder::Construct() : Division" << G4endl;
130 new G4PVDivision(layerName[0], layerLogical[0], tubsLogical, kZAxis,
131 fNSegment[IZ], 0.);
132 }
133 }
134 else if(fNSegment[IZ] == 1)
135 {
136 if(verboseLevel > 9)
137 G4cout << "G4ScoringCylinder::Construct() : Placement" << G4endl;
138 new G4PVPlacement(nullptr, G4ThreeVector(0., 0., 0.), layerLogical[0],
139 layerName[0], tubsLogical, false, 0);
140 }
141 else
142 {
143 G4cerr << "G4ScoringCylinder::SetupGeometry() : invalid parameter ("
144 << fNSegment[IZ] << ") "
145 << "in placement of the first nested layer." << G4endl;
146 }
147
148 // second nested layer (replicated along phi direction)
149 if(verboseLevel > 9)
150 G4cout << "layer 2 :" << G4endl;
151 layerSolid[1] =
152 new G4Tubs(layerName[1], fSize[0], fSize[1], fSize[2] / fNSegment[IZ],
153 fAngle[0], fAngle[1] / fNSegment[IPHI]);
154 layerLogical[1] = new G4LogicalVolume(layerSolid[1], nullptr, layerName[1]);
155 if(fNSegment[IPHI] > 1)
156 {
157 if(verboseLevel > 9)
158 G4cout << "G4ScoringCylinder::Construct() : Replicate along phi direction"
159 << G4endl;
161 {
162 if(verboseLevel > 9)
163 G4cout << "G4ScoringCylinder::Construct() : Replica" << G4endl;
164 new G4PVReplica(layerName[1], layerLogical[1], layerLogical[0], kPhi,
165 fNSegment[IPHI], fAngle[1] / fNSegment[IPHI], fAngle[0]);
166 }
167 else
168 {
169 if(verboseLevel > 9)
170 G4cout << "G4ScoringCylinder::Construct() : Division" << G4endl;
171 new G4PVDivision(layerName[1], layerLogical[1], layerLogical[0], kPhi,
172 fNSegment[IPHI], 0.);
173 }
174 }
175 else if(fNSegment[IPHI] == 1)
176 {
177 if(verboseLevel > 9)
178 G4cout << "G4ScoringCylinder::Construct() : Placement" << G4endl;
179 new G4PVPlacement(nullptr, G4ThreeVector(0., 0., 0.), layerLogical[1],
180 layerName[1], layerLogical[0], false, 0);
181 }
182 else
183 G4cerr << "ERROR : G4ScoringCylinder::SetupGeometry() : invalid parameter ("
184 << fNSegment[IPHI] << ") "
185 << "in placement of the second nested layer." << G4endl;
186
187 // mesh elements
188 if(verboseLevel > 9)
189 G4cout << "mesh elements :" << G4endl;
190 G4String elementName = tubsName + "3";
191 G4VSolid* elementSolid = new G4Tubs(
192 elementName, fSize[0], (fSize[1] - fSize[0]) / fNSegment[IR] + fSize[0],
193 fSize[2] / fNSegment[IZ], fAngle[0], fAngle[1] / fNSegment[IPHI]);
194 fMeshElementLogical = new G4LogicalVolume(elementSolid, nullptr, elementName);
195 if(fNSegment[IR] >= 1)
196 {
197 if(verboseLevel > 9)
198 G4cout << "G4ScoringCylinder::Construct() : Replicate along r direction"
199 << G4endl;
200
202 {
203 if(verboseLevel > 9)
204 G4cout << "G4ScoringCylinder::Construct() : Replica" << G4endl;
205 new G4PVReplica(elementName, fMeshElementLogical, layerLogical[1], kRho,
206 fNSegment[IR], (fSize[1] - fSize[0]) / fNSegment[IR],
207 fSize[0]);
208 }
209 else
210 {
211 if(verboseLevel > 9)
212 G4cout << "G4ScoringCylinder::Construct() : Division" << G4endl;
213 new G4PVDivision(elementName, fMeshElementLogical, layerLogical[1], kRho,
214 fNSegment[IR], 0.);
215 }
216 }
217 else
218 {
219 G4cerr << "G4ScoringCylinder::SetupGeometry() : "
220 << "invalid parameter (" << fNSegment[IR] << ") "
221 << "in mesh element placement." << G4endl;
222 }
223
224 // set the sensitive detector
225 fMeshElementLogical->SetSensitiveDetector(fMFD);
226
227 // vis. attributes
228 auto visatt = new G4VisAttributes(G4Colour(.5, .5, .5));
229 visatt->SetVisibility(true);
230 layerLogical[0]->SetVisAttributes(visatt);
231 layerLogical[1]->SetVisAttributes(visatt);
232 visatt = new G4VisAttributes(G4Colour(.5, .5, .5, 0.01));
233 // visatt->SetForceSolid(true);
234 fMeshElementLogical->SetVisAttributes(visatt);
235
236 if(verboseLevel > 9)
237 DumpVolumes();
238}
void SetVisAttributes(const G4VisAttributes *pVA)
static G4int GetReplicaLevel()
G4LogicalVolume * GetLogicalVolume() const
G4MultiFunctionalDetector * fMFD
G4LogicalVolume * fMeshElementLogical
@ kPhi
Definition geomdefs.hh:60
@ kZAxis
Definition geomdefs.hh:57
@ kRho
Definition geomdefs.hh:58

◆ SetZSize()

void G4ScoringCylinder::SetZSize ( G4double zSize)
inline

Definition at line 57 of file G4ScoringCylinder.hh.

57{ fSize[2] = zSize; } // half height

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