Geant4-11
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions
G4ScoringCylinder Class Reference

#include <G4ScoringCylinder.hh>

Inheritance diagram for G4ScoringCylinder:
G4VScoringMesh

Public Types

using EventScore = G4THitsMap< G4double >
 
enum  IDX { IZ , IPHI , IR }
 
using MeshScoreMap = std::map< G4String, RunScore * >
 
enum class  MeshShape {
  box , cylinder , sphere , realWorldLogVol ,
  probe , undefined = -1
}
 
using RunScore = G4THitsMap< G4StatDouble >
 

Public Member Functions

void Accumulate (G4THitsMap< G4double > *map)
 
void Accumulate (G4THitsMap< G4StatDouble > *map)
 
void Activate (G4bool vl=true)
 
virtual void Construct (G4VPhysicalVolume *fWorldPhys)
 
virtual void Draw (RunScore *map, G4VScoreColorMap *colorMap, G4int axflg=111)
 
virtual void DrawColumn (RunScore *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn)
 
void DrawMesh (const G4String &psName, G4int idxPlane, G4int iColumn, G4VScoreColorMap *colorMap)
 
void DrawMesh (const G4String &psName, G4VScoreColorMap *colorMap, G4int axflg=111)
 
void Dump ()
 
G4bool FindPrimitiveScorer (const G4String &psname)
 
 G4ScoringCylinder (G4String wName)
 
void GeometryHasBeenDestroyed ()
 
G4double GetAngleSpan () const
 
G4int GetCopyNumberLevel () const
 
G4String GetCurrentPSUnit ()
 
void GetDivisionAxisNames (G4String divisionAxisNames[3])
 
G4LogicalVolumeGetMeshElementLogical () const
 
void GetNumberOfSegments (G4int nSegment[3])
 
G4ParallelWorldProcessGetParallelWorldProcess () const
 
G4VPrimitiveScorerGetPrimitiveScorer (const G4String &name)
 
G4String GetPSUnit (const G4String &psname)
 
G4double GetPSUnitValue (const G4String &psname)
 
G4RotationMatrix GetRotationMatrix () const
 
void GetRZPhi (G4int index, G4int q[3]) const
 
MeshScoreMap GetScoreMap () const
 
MeshShape GetShape () const
 
G4ThreeVector GetSize () const
 
G4double GetStartAngle () const
 
G4ThreeVector GetTranslation () const
 
const G4StringGetWorldName () const
 
G4bool IsActive () const
 
G4bool IsCurrentPrimitiveScorerNull ()
 
G4bool LayeredMassFlg ()
 
virtual void List () const
 
void Merge (const G4VScoringMesh *scMesh)
 
G4bool ReadyForQuantity () const
 
void RegisterPrimitives (std::vector< G4VPrimitiveScorer * > &vps)
 
void ResetScore ()
 
void RotateX (G4double delta)
 
void RotateY (G4double delta)
 
void RotateZ (G4double delta)
 
void SetAngles (G4double, G4double)
 
void SetCenterPosition (G4double centerPosition[3])
 
void SetCopyNumberLevel (G4int val)
 
void SetCurrentPrimitiveScorer (const G4String &name)
 
void SetCurrentPSUnit (const G4String &unit)
 
void SetDrawPSName (const G4String &psname)
 
void SetFilter (G4VSDFilter *filter)
 
void SetMeshElementLogical (G4LogicalVolume *val)
 
void SetNullToCurrentPrimitiveScorer ()
 
void SetNumberOfSegments (G4int nSegment[3])
 
void SetParallelWorldProcess (G4ParallelWorldProcess *proc)
 
void SetPrimitiveScorer (G4VPrimitiveScorer *ps)
 
void SetRMax (G4double rMax)
 
void SetRMin (G4double rMin)
 
void SetSize (G4double size[3])
 
void SetVerboseLevel (G4int vl)
 
void SetZSize (G4double zSize)
 
virtual void WorkerConstruct (G4VPhysicalVolume *fWorldPhys)
 
 ~G4ScoringCylinder ()
 

Protected Member Functions

virtual void SetupGeometry (G4VPhysicalVolume *fWorldPhys)
 

Protected Attributes

G4int copyNumberLevel
 
G4bool fActive
 
G4double fAngle [2]
 
G4ThreeVector fCenterPosition
 
G4bool fConstructed
 
G4VPrimitiveScorerfCurrentPS
 
G4String fDivisionAxisNames [3]
 
G4String fDrawPSName
 
G4String fDrawUnit
 
G4double fDrawUnitValue
 
G4bool fGeometryHasBeenDestroyed
 
MeshScoreMap fMap
 
G4LogicalVolumefMeshElementLogical
 
G4MultiFunctionalDetectorfMFD
 
G4int fNSegment [3]
 
G4ParallelWorldProcessfParallelWorldProcess
 
G4RotationMatrixfRotationMatrix
 
MeshShape fShape
 
G4double fSize [3]
 
G4String fWorldName
 
G4bool layeredMassFlg
 
G4bool nMeshIsSet
 
G4bool sizeIsSet
 
G4int verboseLevel
 

Private Member Functions

void DumpLogVols (G4int)
 
void DumpPhysVols (G4int)
 
void DumpSolids (G4int)
 
void DumpVolumes ()
 

Detailed Description

Definition at line 40 of file G4ScoringCylinder.hh.

Member Typedef Documentation

◆ EventScore

Definition at line 66 of file G4VScoringMesh.hh.

◆ MeshScoreMap

using G4VScoringMesh::MeshScoreMap = std::map<G4String, RunScore*>
inherited

Definition at line 68 of file G4VScoringMesh.hh.

◆ RunScore

Definition at line 67 of file G4VScoringMesh.hh.

Member Enumeration Documentation

◆ IDX

Enumerator
IZ 
IPHI 
IR 

Definition at line 69 of file G4ScoringCylinder.hh.

◆ MeshShape

enum class G4VScoringMesh::MeshShape
stronginherited
Enumerator
box 
cylinder 
sphere 
realWorldLogVol 
probe 
undefined 

Definition at line 57 of file G4VScoringMesh.hh.

58 {
59 box,
60 cylinder,
61 sphere,
62 realWorldLogVol,
63 probe,
64 undefined = -1
65 };

Constructor & Destructor Documentation

◆ G4ScoringCylinder()

G4ScoringCylinder::G4ScoringCylinder ( G4String  wName)

Definition at line 55 of file G4ScoringCylinder.cc.

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

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

◆ ~G4ScoringCylinder()

G4ScoringCylinder::~G4ScoringCylinder ( )

Definition at line 65 of file G4ScoringCylinder.cc.

65{ ; }

Member Function Documentation

◆ Accumulate() [1/2]

void G4VScoringMesh::Accumulate ( G4THitsMap< G4double > *  map)
inherited

Definition at line 386 of file G4VScoringMesh.cc.

387{
388 G4String psName = map->GetName();
389 MeshScoreMap::const_iterator fMapItr = fMap.find(psName);
390 *(fMapItr->second) += *map;
391
392 if(verboseLevel > 9)
393 {
394 G4cout << G4endl;
395 G4cout << "G4VScoringMesh::Accumulate()" << G4endl;
396 G4cout << " PS name : " << psName << G4endl;
397 if(fMapItr == fMap.end())
398 {
399 G4cout << " " << psName << " was not found." << G4endl;
400 }
401 else
402 {
403 G4cout << " map size : " << map->GetSize() << G4endl;
404 map->PrintAllHits();
405 }
406 G4cout << G4endl;
407 }
408}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
MeshScoreMap fMap

References G4VScoringMesh::fMap, G4cout, G4endl, anonymous_namespace{G4QuasiElRatios.cc}::map, and G4VScoringMesh::verboseLevel.

◆ Accumulate() [2/2]

void G4VScoringMesh::Accumulate ( G4THitsMap< G4StatDouble > *  map)
inherited

Definition at line 410 of file G4VScoringMesh.cc.

411{
412 G4String psName = map->GetName();
413 MeshScoreMap::const_iterator fMapItr = fMap.find(psName);
414 *(fMapItr->second) += *map;
415
416 if(verboseLevel > 9)
417 {
418 G4cout << G4endl;
419 G4cout << "G4VScoringMesh::Accumulate()" << G4endl;
420 G4cout << " PS name : " << psName << G4endl;
421 if(fMapItr == fMap.end())
422 {
423 G4cout << " " << psName << " was not found." << G4endl;
424 }
425 else
426 {
427 G4cout << " map size : " << map->GetSize() << G4endl;
428 map->PrintAllHits();
429 }
430 G4cout << G4endl;
431 }
432}

References G4VScoringMesh::fMap, G4cout, G4endl, anonymous_namespace{G4QuasiElRatios.cc}::map, and G4VScoringMesh::verboseLevel.

◆ Activate()

void G4VScoringMesh::Activate ( G4bool  vl = true)
inlineinherited

Definition at line 95 of file G4VScoringMesh.hh.

95{ fActive = vl; }

References G4VScoringMesh::fActive.

◆ Construct()

void G4VScoringMesh::Construct ( G4VPhysicalVolume fWorldPhys)
virtualinherited

Definition at line 434 of file G4VScoringMesh.cc.

435{
436 if(fConstructed)
437 {
439 {
440 SetupGeometry(fWorldPhys);
442 }
443 if(verboseLevel > 0)
444 G4cout << fWorldName << " --- All quantities are reset." << G4endl;
445 ResetScore();
446 }
447 else
448 {
449 fConstructed = true;
450 SetupGeometry(fWorldPhys);
451 }
452}
virtual void SetupGeometry(G4VPhysicalVolume *fWorldPhys)=0
G4bool fGeometryHasBeenDestroyed

References G4VScoringMesh::fConstructed, G4VScoringMesh::fGeometryHasBeenDestroyed, G4VScoringMesh::fWorldName, G4cout, G4endl, G4VScoringMesh::ResetScore(), G4VScoringMesh::SetupGeometry(), and G4VScoringMesh::verboseLevel.

Referenced by G4RunManager::ConstructScoringWorlds().

◆ Draw()

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

Implements G4VScoringMesh.

Definition at line 259 of file G4ScoringCylinder.cc.

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

References angle, G4VVisManager::BeginDraw(), DBL_MAX, G4VVisManager::Draw(), G4VScoreColorMap::DrawColorChart(), G4VVisManager::EndDraw(), G4VScoringMesh::fAngle, G4VScoringMesh::fCenterPosition, G4VScoringMesh::fDrawPSName, G4VScoringMesh::fDrawUnit, G4VScoringMesh::fDrawUnitValue, G4VScoringMesh::fNSegment, G4VScoringMesh::fRotationMatrix, G4VScoringMesh::fSize, G4cout, G4endl, G4VVisManager::GetConcreteInstance(), G4VScoreColorMap::GetMapColor(), G4CSGSolid::GetPolyhedron(), GetRZPhi(), G4VScoreColorMap::IfFloatMinMax(), HepGeom::Transform3D::inverse(), IPHI, IR, IZ, anonymous_namespace{G4QuasiElRatios.cc}::map, G4VisAttributes::SetColour(), G4VisAttributes::SetForceAuxEdgeVisible(), G4VisAttributes::SetForceSolid(), G4VScoreColorMap::SetMinMax(), G4VScoreColorMap::SetPSName(), G4VScoreColorMap::SetPSUnit(), G4Visible::SetVisAttributes(), and HepPolyhedron::Transform().

◆ DrawColumn()

void G4ScoringCylinder::DrawColumn ( RunScore map,
G4VScoreColorMap colorMap,
G4int  idxProj,
G4int  idxColumn 
)
virtual

Implements G4VScoringMesh.

Definition at line 438 of file G4ScoringCylinder.cc.

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

References angle, G4VVisManager::BeginDraw(), DBL_MAX, G4VVisManager::Draw(), G4VScoreColorMap::DrawColorChart(), G4VVisManager::EndDraw(), G4VScoringMesh::fAngle, G4VScoringMesh::fCenterPosition, G4VScoringMesh::fDrawPSName, G4VScoringMesh::fDrawUnit, G4VScoringMesh::fDrawUnitValue, G4VScoringMesh::fNSegment, G4VScoringMesh::fRotationMatrix, G4VScoringMesh::fSize, G4cerr, G4cout, G4endl, G4VVisManager::GetConcreteInstance(), G4VScoreColorMap::GetMapColor(), G4CSGSolid::GetPolyhedron(), GetRZPhi(), G4VScoreColorMap::IfFloatMinMax(), HepGeom::Transform3D::inverse(), IPHI, IR, IZ, anonymous_namespace{G4QuasiElRatios.cc}::map, G4VisAttributes::SetColour(), G4VisAttributes::SetForceAuxEdgeVisible(), G4VisAttributes::SetForceSolid(), G4VScoreColorMap::SetMinMax(), G4VScoreColorMap::SetPSName(), G4VScoreColorMap::SetPSUnit(), G4Visible::SetVisAttributes(), and HepPolyhedron::Transform().

◆ DrawMesh() [1/2]

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

Definition at line 368 of file G4VScoringMesh.cc.

370{
371 fDrawPSName = psName;
372 MeshScoreMap::const_iterator fMapItr = fMap.find(psName);
373 if(fMapItr != fMap.end())
374 {
375 fDrawUnit = GetPSUnit(psName);
377 DrawColumn(fMapItr->second, colorMap, idxPlane, iColumn);
378 }
379 else
380 {
381 G4cerr << "Scorer <" << psName << "> is not defined. Method ignored."
382 << G4endl;
383 }
384}
G4double GetPSUnitValue(const G4String &psname)
G4String GetPSUnit(const G4String &psname)
virtual void DrawColumn(RunScore *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn)=0

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

◆ DrawMesh() [2/2]

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

Definition at line 350 of file G4VScoringMesh.cc.

352{
353 fDrawPSName = psName;
354 MeshScoreMap::const_iterator fMapItr = fMap.find(psName);
355 if(fMapItr != fMap.end())
356 {
357 fDrawUnit = GetPSUnit(psName);
359 Draw(fMapItr->second, colorMap, axflg);
360 }
361 else
362 {
363 G4cerr << "Scorer <" << psName << "> is not defined. Method ignored."
364 << G4endl;
365 }
366}
virtual void Draw(RunScore *map, G4VScoreColorMap *colorMap, G4int axflg=111)=0

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

Referenced by G4VSceneHandler::AddCompound(), and G4ScoringManager::DrawMesh().

◆ Dump()

void G4VScoringMesh::Dump ( )
inherited

Definition at line 338 of file G4VScoringMesh.cc.

339{
340 G4cout << "scoring mesh name: " << fWorldName << G4endl;
341 G4cout << "# of G4THitsMap : " << fMap.size() << G4endl;
342 for(auto mp : fMap)
343 {
344 G4cout << "[" << mp.first << "]" << G4endl;
345 mp.second->PrintAllHits();
346 }
347 G4cout << G4endl;
348}

References G4VScoringMesh::fMap, G4VScoringMesh::fWorldName, G4cout, and G4endl.

◆ DumpLogVols()

void G4ScoringCylinder::DumpLogVols ( G4int  lvl)
private

Definition at line 753 of file G4ScoringCylinder.cc.

754{
755 G4cout << "*********** List of registered logical volumes *************"
756 << G4endl;
758 auto itr = store->begin();
759 for(; itr != store->end(); itr++)
760 {
761 G4cout << (*itr)->GetName()
762 << "\t Solid = " << (*itr)->GetSolid()->GetName();
763 if((*itr)->GetMaterial())
764 {
765 G4cout << "\t Material = " << (*itr)->GetMaterial()->GetName() << G4endl;
766 }
767 else
768 {
769 G4cout << "\t Material : not defined " << G4endl;
770 }
771 if(lvl < 1)
772 continue;
773 G4cout << "\t region = ";
774 if((*itr)->GetRegion())
775 {
776 G4cout << (*itr)->GetRegion()->GetName();
777 }
778 else
779 {
780 G4cout << "not defined";
781 }
782 G4cout << "\t sensitive detector = ";
783 if((*itr)->GetSensitiveDetector())
784 {
785 G4cout << (*itr)->GetSensitiveDetector()->GetName();
786 }
787 else
788 {
789 G4cout << "not defined";
790 }
791 G4cout << G4endl;
792 G4cout << "\t daughters = " << (*itr)->GetNoDaughters();
793 if((*itr)->GetNoDaughters() > 0)
794 {
795 switch((*itr)->CharacteriseDaughters())
796 {
797 case kNormal:
798 G4cout << " (placement)";
799 break;
800 case kReplica:
801 G4cout << " (replica : " << (*itr)->GetDaughter(0)->GetMultiplicity()
802 << ")";
803 break;
804 case kParameterised:
805 G4cout << " (parameterized : "
806 << (*itr)->GetDaughter(0)->GetMultiplicity() << ")";
807 break;
808 default:;
809 }
810 }
811 G4cout << G4endl;
812 if(lvl < 2)
813 continue;
814 if((*itr)->GetMaterial())
815 {
816 G4cout << "\t weight = " << G4BestUnit((*itr)->GetMass(), "Mass")
817 << G4endl;
818 }
819 else
820 {
821 G4cout << "\t weight : not available" << G4endl;
822 }
823 }
824}
#define G4BestUnit(a, b)
static G4LogicalVolumeStore * GetInstance()
@ kNormal
Definition: geomdefs.hh:84
@ kParameterised
Definition: geomdefs.hh:86
@ kReplica
Definition: geomdefs.hh:85

References G4BestUnit, G4cout, G4endl, G4LogicalVolumeStore::GetInstance(), kNormal, kParameterised, and kReplica.

Referenced by DumpVolumes().

◆ DumpPhysVols()

void G4ScoringCylinder::DumpPhysVols ( G4int  lvl)
private

Definition at line 826 of file G4ScoringCylinder.cc.

827{
828 G4cout << "*********** List of registered physical volumes *************"
829 << G4endl;
831 auto itr = store->begin();
832 for(; itr != store->end(); itr++)
833 {
834 switch(lvl)
835 {
836 case 0:
837 G4cout << (*itr)->GetName() << G4endl;
838 break;
839 case 1:
840 G4cout << (*itr)->GetName() << "\t logical volume = "
841 << (*itr)->GetLogicalVolume()->GetName()
842 << "\t mother logical = ";
843 if((*itr)->GetMotherLogical())
844 {
845 G4cout << (*itr)->GetMotherLogical()->GetName();
846 }
847 else
848 {
849 G4cout << "not defined";
850 }
851 G4cout << G4endl;
852 break;
853 default:
854 G4cout << (*itr)->GetName() << "\t logical volume = "
855 << (*itr)->GetLogicalVolume()->GetName()
856 << "\t mother logical = ";
857 if((*itr)->GetMotherLogical())
858 {
859 G4cout << (*itr)->GetMotherLogical()->GetName();
860 }
861 else
862 {
863 G4cout << "not defined";
864 }
865 G4cout << "\t type = ";
866 switch((*itr)->VolumeType())
867 {
868 case kNormal:
869 G4cout << "placement";
870 break;
871 case kReplica:
872 G4cout << "replica";
873 break;
874 case kParameterised:
875 G4cout << "parameterized";
876 break;
877 default:;
878 }
879 G4cout << G4endl;
880 }
881 }
882}
static G4PhysicalVolumeStore * GetInstance()

References G4cout, G4endl, G4PhysicalVolumeStore::GetInstance(), kNormal, kParameterised, and kReplica.

Referenced by DumpVolumes().

◆ DumpSolids()

void G4ScoringCylinder::DumpSolids ( G4int  lvl)
private

Definition at line 728 of file G4ScoringCylinder.cc.

729{
730 G4cout << "*********** List of registered solids *************" << G4endl;
731 auto store = G4SolidStore::GetInstance();
732 auto itr = store->begin();
733 for(; itr != store->end(); itr++)
734 {
735 switch(lvl)
736 {
737 case 0:
738 G4cout << (*itr)->GetName() << G4endl;
739 break;
740 case 1:
741 G4cout << (*itr)->GetName() << "\t volume = "
742 << G4BestUnit((*itr)->GetCubicVolume(), "Volume")
743 << "\t surface = "
744 << G4BestUnit((*itr)->GetSurfaceArea(), "Surface") << G4endl;
745 break;
746 default:
747 (*itr)->DumpInfo();
748 break;
749 }
750 }
751}
static G4SolidStore * GetInstance()

References G4BestUnit, G4cout, G4endl, and G4SolidStore::GetInstance().

Referenced by DumpVolumes().

◆ DumpVolumes()

void G4ScoringCylinder::DumpVolumes ( )
private

Definition at line 720 of file G4ScoringCylinder.cc.

721{
722 G4int lvl = 2;
723 DumpSolids(lvl);
724 DumpLogVols(lvl);
725 DumpPhysVols(lvl);
726}

References DumpLogVols(), DumpPhysVols(), and DumpSolids().

Referenced by SetupGeometry().

◆ FindPrimitiveScorer()

G4bool G4VScoringMesh::FindPrimitiveScorer ( const G4String psname)
inherited

Definition at line 223 of file G4VScoringMesh.cc.

224{
225 MeshScoreMap::iterator itr = fMap.find(psname);
226 if(itr == fMap.end())
227 return false;
228 return true;
229}

References G4VScoringMesh::fMap.

Referenced by G4ScoreQuantityMessenger::CheckMeshPS().

◆ GeometryHasBeenDestroyed()

void G4VScoringMesh::GeometryHasBeenDestroyed ( )
inlineinherited

◆ GetAngleSpan()

G4double G4VScoringMesh::GetAngleSpan ( ) const
inlineinherited

Definition at line 129 of file G4VScoringMesh.hh.

129{ return fAngle[1]; }

References G4VScoringMesh::fAngle.

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ GetCopyNumberLevel()

G4int G4VScoringMesh::GetCopyNumberLevel ( ) const
inlineinherited

Definition at line 263 of file G4VScoringMesh.hh.

263{ return copyNumberLevel; }

References G4VScoringMesh::copyNumberLevel.

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ GetCurrentPSUnit()

G4String G4VScoringMesh::GetCurrentPSUnit ( )
inherited

Definition at line 244 of file G4VScoringMesh.cc.

245{
246 G4String unit = "";
247 if(!fCurrentPS)
248 {
249 G4String msg = "ERROR : G4VScoringMesh::GetCurrentPSUnit() : ";
250 msg += " Current primitive scorer is null.";
251 G4cerr << msg << G4endl;
252 }
253 else
254 {
255 unit = fCurrentPS->GetUnit();
256 }
257 return unit;
258}
const G4String & GetUnit() const
G4VPrimitiveScorer * fCurrentPS

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

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ GetDivisionAxisNames()

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

Definition at line 287 of file G4VScoringMesh.cc.

288{
289 for(int i = 0; i < 3; i++)
290 divisionAxisNames[i] = fDivisionAxisNames[i];
291}

References G4VScoringMesh::fDivisionAxisNames.

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

◆ GetMeshElementLogical()

G4LogicalVolume * G4VScoringMesh::GetMeshElementLogical ( ) const
inlineinherited

Definition at line 232 of file G4VScoringMesh.hh.

233 {
234 return fMeshElementLogical;
235 }

References G4VScoringMesh::fMeshElementLogical.

Referenced by G4WorkerRunManager::ConstructScoringWorlds().

◆ GetNumberOfSegments()

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

Definition at line 140 of file G4VScoringMesh.cc.

141{
142 for(int i = 0; i < 3; i++)
143 nSegment[i] = fNSegment[i];
144}

References G4VScoringMesh::fNSegment.

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

◆ GetParallelWorldProcess()

G4ParallelWorldProcess * G4VScoringMesh::GetParallelWorldProcess ( ) const
inlineinherited

Definition at line 246 of file G4VScoringMesh.hh.

247 {
249 }
G4ParallelWorldProcess * fParallelWorldProcess

References G4VScoringMesh::fParallelWorldProcess.

Referenced by G4RunManager::ConstructScoringWorlds(), and G4WorkerRunManager::ConstructScoringWorlds().

◆ GetPrimitiveScorer()

G4VPrimitiveScorer * G4VScoringMesh::GetPrimitiveScorer ( const G4String name)
inherited

◆ GetPSUnit()

G4String G4VScoringMesh::GetPSUnit ( const G4String psname)
inherited

Definition at line 231 of file G4VScoringMesh.cc.

232{
233 MeshScoreMap::iterator itr = fMap.find(psname);
234 if(itr == fMap.end())
235 {
236 return G4String("");
237 }
238 else
239 {
240 return GetPrimitiveScorer(psname)->GetUnit();
241 }
242}
G4VPrimitiveScorer * GetPrimitiveScorer(const G4String &name)

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

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

◆ GetPSUnitValue()

G4double G4VScoringMesh::GetPSUnitValue ( const G4String psname)
inherited

Definition at line 274 of file G4VScoringMesh.cc.

275{
276 MeshScoreMap::iterator itr = fMap.find(psname);
277 if(itr == fMap.end())
278 {
279 return 1.;
280 }
281 else
282 {
283 return GetPrimitiveScorer(psname)->GetUnitValue();
284 }
285}
G4double GetUnitValue() const

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

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

◆ GetRotationMatrix()

G4RotationMatrix G4VScoringMesh::GetRotationMatrix ( ) const
inlineinherited

Definition at line 141 of file G4VScoringMesh.hh.

142 {
144 return *fRotationMatrix;
145 else
147 }
static DLL_API const HepRotation IDENTITY
Definition: Rotation.h:366

References G4VScoringMesh::fRotationMatrix, and CLHEP::HepRotation::IDENTITY.

Referenced by G4GMocrenFileSceneHandler::AddSolid().

◆ GetRZPhi()

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

Definition at line 696 of file G4ScoringCylinder.cc.

697{
698 // index = k + j * k-size + i * jk-plane-size
699
700 // nested : z -> phi -> r
701 G4int i = IZ;
702 G4int j = IPHI;
703 G4int k = IR;
704 G4int jk = fNSegment[j] * fNSegment[k];
705 q[i] = index / jk;
706 q[j] = (index - q[i] * jk) / fNSegment[k];
707 q[k] = index - q[j] * fNSegment[k] - q[i] * jk;
708}

References G4VScoringMesh::fNSegment, IPHI, IR, and IZ.

Referenced by Draw(), and DrawColumn().

◆ GetScoreMap()

MeshScoreMap G4VScoringMesh::GetScoreMap ( ) const
inlineinherited

◆ GetShape()

MeshShape G4VScoringMesh::GetShape ( ) const
inlineinherited

◆ GetSize()

G4ThreeVector G4VScoringMesh::GetSize ( ) const
inherited

Definition at line 105 of file G4VScoringMesh.cc.

106{
107 if(sizeIsSet)
108 return G4ThreeVector(fSize[0], fSize[1], fSize[2]);
109 else
110 return G4ThreeVector(0., 0., 0.);
111}
CLHEP::Hep3Vector G4ThreeVector

References G4VScoringMesh::fSize, and G4VScoringMesh::sizeIsSet.

Referenced by G4GMocrenFileSceneHandler::AddSolid(), G4ScoreQuantityMessenger::SetNewValue(), and G4ScoringMessenger::SetNewValue().

◆ GetStartAngle()

G4double G4VScoringMesh::GetStartAngle ( ) const
inlineinherited

Definition at line 128 of file G4VScoringMesh.hh.

128{ return fAngle[0]; }

References G4VScoringMesh::fAngle.

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ GetTranslation()

G4ThreeVector G4VScoringMesh::GetTranslation ( ) const
inlineinherited

Definition at line 133 of file G4VScoringMesh.hh.

133{ return fCenterPosition; }

References G4VScoringMesh::fCenterPosition.

Referenced by G4GMocrenFileSceneHandler::AddSolid().

◆ GetWorldName()

const G4String & G4VScoringMesh::GetWorldName ( ) const
inlineinherited

◆ IsActive()

G4bool G4VScoringMesh::IsActive ( ) const
inlineinherited

Definition at line 93 of file G4VScoringMesh.hh.

93{ return fActive; }

References G4VScoringMesh::fActive.

Referenced by G4VSceneHandler::AddCompound(), and G4PSHitsModel::DescribeYourselfTo().

◆ IsCurrentPrimitiveScorerNull()

G4bool G4VScoringMesh::IsCurrentPrimitiveScorerNull ( )
inlineinherited

Definition at line 164 of file G4VScoringMesh.hh.

165 {
166 if(fCurrentPS == nullptr)
167 return true;
168 else
169 return false;
170 }

References G4VScoringMesh::fCurrentPS.

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ LayeredMassFlg()

G4bool G4VScoringMesh::LayeredMassFlg ( )
inlineinherited

◆ List()

void G4ScoringCylinder::List ( ) const
virtual

Reimplemented from G4VScoringMesh.

Definition at line 246 of file G4ScoringCylinder.cc.

247{
248 G4cout << "G4ScoringCylinder : " << fWorldName
249 << " --- Shape: Cylindrical mesh" << G4endl;
250
251 G4cout << " Size (Rmin, Rmax, Dz): (" << fSize[0] / cm << ", "
252 << fSize[1] / cm << ", " << fSize[2] / cm << ") [cm]" << G4endl;
253 G4cout << " Angle (start, span): (" << fAngle[0] / deg << ", "
254 << fAngle[1] / deg << ") [deg]" << G4endl;
255
257}
static constexpr double cm
Definition: G4SIunits.hh:99
static constexpr double deg
Definition: G4SIunits.hh:132
virtual void List() const

References cm, deg, G4VScoringMesh::fAngle, G4VScoringMesh::fSize, G4VScoringMesh::fWorldName, G4cout, G4endl, and G4VScoringMesh::List().

◆ Merge()

void G4VScoringMesh::Merge ( const G4VScoringMesh scMesh)
inherited

Definition at line 476 of file G4VScoringMesh.cc.

477{
478 const MeshScoreMap scMap = scMesh->GetScoreMap();
479
480 MeshScoreMap::const_iterator fMapItr = fMap.begin();
481 MeshScoreMap::const_iterator mapItr = scMap.begin();
482 for(; fMapItr != fMap.end(); fMapItr++)
483 {
484 if(verboseLevel > 9)
485 G4cout << "G4VScoringMesh::Merge()" << fMapItr->first << G4endl;
486 *(fMapItr->second) += *(mapItr->second);
487 mapItr++;
488 }
489}
std::map< G4String, RunScore * > MeshScoreMap
MeshScoreMap GetScoreMap() const

References G4VScoringMesh::fMap, G4cout, G4endl, G4VScoringMesh::GetScoreMap(), and G4VScoringMesh::verboseLevel.

Referenced by G4ScoringManager::Merge().

◆ ReadyForQuantity()

G4bool G4VScoringMesh::ReadyForQuantity ( ) const
inlineinherited

Definition at line 192 of file G4VScoringMesh.hh.

192{ return (sizeIsSet && nMeshIsSet); }

References G4VScoringMesh::nMeshIsSet, and G4VScoringMesh::sizeIsSet.

Referenced by G4VScoringMesh::SetPrimitiveScorer().

◆ RegisterPrimitives()

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

◆ ResetScore()

void G4VScoringMesh::ResetScore ( )
inherited

Definition at line 75 of file G4VScoringMesh.cc.

76{
77 if(verboseLevel > 9)
78 G4cout << "G4VScoringMesh::ResetScore() is called." << G4endl;
79 for(auto mp : fMap)
80 {
81 if(verboseLevel > 9)
82 G4cout << "G4VScoringMesh::ResetScore()" << mp.first << G4endl;
83 mp.second->clear();
84 }
85}

References G4VScoringMesh::fMap, G4cout, G4endl, and G4VScoringMesh::verboseLevel.

Referenced by G4VScoringMesh::Construct(), and G4VScoringMesh::WorkerConstruct().

◆ RotateX()

void G4VScoringMesh::RotateX ( G4double  delta)
inherited

Definition at line 145 of file G4VScoringMesh.cc.

146{
147 if(!fRotationMatrix)
149 fRotationMatrix->rotateX(delta);
150}
CLHEP::HepRotation G4RotationMatrix
HepRotation & rotateX(double delta)
Definition: Rotation.cc:61

References G4VScoringMesh::fRotationMatrix, and CLHEP::HepRotation::rotateX().

Referenced by G4ScoringMessenger::SetNewValue().

◆ RotateY()

void G4VScoringMesh::RotateY ( G4double  delta)
inherited

Definition at line 152 of file G4VScoringMesh.cc.

153{
154 if(!fRotationMatrix)
156 fRotationMatrix->rotateY(delta);
157}
HepRotation & rotateY(double delta)
Definition: Rotation.cc:74

References G4VScoringMesh::fRotationMatrix, and CLHEP::HepRotation::rotateY().

Referenced by G4ScoringMessenger::SetNewValue().

◆ RotateZ()

void G4VScoringMesh::RotateZ ( G4double  delta)
inherited

Definition at line 159 of file G4VScoringMesh.cc.

160{
161 if(!fRotationMatrix)
163 fRotationMatrix->rotateZ(delta);
164}
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:87

References G4VScoringMesh::fRotationMatrix, and CLHEP::HepRotation::rotateZ().

Referenced by G4ScoringMessenger::SetNewValue().

◆ SetAngles()

void G4VScoringMesh::SetAngles ( G4double  startAngle,
G4double  spanAngle 
)
inherited

Definition at line 112 of file G4VScoringMesh.cc.

113{
114 fAngle[0] = startAngle;
115 fAngle[1] = spanAngle;
116}

References G4VScoringMesh::fAngle.

Referenced by G4ScoringMessenger::SetNewValue().

◆ SetCenterPosition()

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

Definition at line 118 of file G4VScoringMesh.cc.

119{
121 G4ThreeVector(centerPosition[0], centerPosition[1], centerPosition[2]);
122}

References G4VScoringMesh::fCenterPosition.

Referenced by G4ScoringMessenger::SetNewValue().

◆ SetCopyNumberLevel()

void G4VScoringMesh::SetCopyNumberLevel ( G4int  val)
inlineinherited

Definition at line 262 of file G4VScoringMesh.hh.

262{ copyNumberLevel = val; }

References G4VScoringMesh::copyNumberLevel.

◆ SetCurrentPrimitiveScorer()

void G4VScoringMesh::SetCurrentPrimitiveScorer ( const G4String name)
inherited

Definition at line 212 of file G4VScoringMesh.cc.

213{
215 if(!fCurrentPS)
216 {
217 G4cerr << "ERROR : G4VScoringMesh::SetCurrentPrimitiveScorer() : The "
218 "primitive scorer <"
219 << name << "> does not found." << G4endl;
220 }
221}

References G4VScoringMesh::fCurrentPS, G4cerr, G4endl, G4VScoringMesh::GetPrimitiveScorer(), and G4InuclParticleNames::name().

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetCurrentPSUnit()

void G4VScoringMesh::SetCurrentPSUnit ( const G4String unit)
inherited

Definition at line 260 of file G4VScoringMesh.cc.

261{
262 if(!fCurrentPS)
263 {
264 G4String msg = "ERROR : G4VScoringMesh::GetCurrentPSUnit() : ";
265 msg += " Current primitive scorer is null.";
266 G4cerr << msg << G4endl;
267 }
268 else
269 {
270 fCurrentPS->SetUnit(unit);
271 }
272}
void SetUnit(const G4String &unit)

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

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetDrawPSName()

void G4VScoringMesh::SetDrawPSName ( const G4String psname)
inlineinherited

Definition at line 180 of file G4VScoringMesh.hh.

180{ fDrawPSName = psname; }

References G4VScoringMesh::fDrawPSName.

◆ SetFilter()

void G4VScoringMesh::SetFilter ( G4VSDFilter filter)
inherited

Definition at line 190 of file G4VScoringMesh.cc.

191{
192 if(!fCurrentPS)
193 {
194 G4cerr << "ERROR : G4VScoringMesh::SetSDFilter() : a quantity must be "
195 "defined first. This method is ignored."
196 << G4endl;
197 return;
198 }
199 if(verboseLevel > 0)
200 G4cout << "G4VScoringMesh::SetFilter() : " << filter->GetName()
201 << " is set to " << fCurrentPS->GetName() << G4endl;
202
203 G4VSDFilter* oldFilter = fCurrentPS->GetFilter();
204 if(oldFilter)
205 {
206 G4cout << "WARNING : G4VScoringMesh::SetFilter() : " << oldFilter->GetName()
207 << " is overwritten by " << filter->GetName() << G4endl;
208 }
209 fCurrentPS->SetFilter(filter);
210}
G4VSDFilter * GetFilter() const
void SetFilter(G4VSDFilter *f)
G4String GetName() const
Definition: G4VSDFilter.hh:55

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

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

◆ SetMeshElementLogical()

void G4VScoringMesh::SetMeshElementLogical ( G4LogicalVolume val)
inlineinherited

Definition at line 228 of file G4VScoringMesh.hh.

229 {
231 }

References G4VScoringMesh::fMeshElementLogical.

Referenced by G4WorkerRunManager::ConstructScoringWorlds().

◆ SetNullToCurrentPrimitiveScorer()

void G4VScoringMesh::SetNullToCurrentPrimitiveScorer ( )
inlineinherited

Definition at line 186 of file G4VScoringMesh.hh.

186{ fCurrentPS = nullptr; }

References G4VScoringMesh::fCurrentPS.

Referenced by G4ScoreQuantityMessenger::CheckMeshPS().

◆ SetNumberOfSegments()

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

Definition at line 123 of file G4VScoringMesh.cc.

124{
127 {
128 for(int i = 0; i < 3; i++)
129 fNSegment[i] = nSegment[i];
130 nMeshIsSet = true;
131 }
132 else
133 {
134 G4String message = " Number of bins has already been set and it cannot be changed.\n";
135 message += " This method is ignored.";
136 G4Exception("G4VScoringMesh::SetNumberOfSegments()",
137 "DigiHitsUtilsScoreVScoringMesh000", JustWarning, message);
138 }
139}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35

References G4VScoringMesh::fNSegment, G4VScoringMesh::fShape, G4Exception(), JustWarning, G4VScoringMesh::nMeshIsSet, G4VScoringMesh::probe, and G4VScoringMesh::realWorldLogVol.

Referenced by G4ScoringProbe::G4ScoringProbe(), G4ScoringRealWorld::G4ScoringRealWorld(), G4ScoringProbe::LocateProbe(), G4ScoringMessenger::MeshBinCommand(), and G4ScoringRealWorld::SetupGeometry().

◆ SetParallelWorldProcess()

void G4VScoringMesh::SetParallelWorldProcess ( G4ParallelWorldProcess proc)
inlineinherited

◆ SetPrimitiveScorer()

void G4VScoringMesh::SetPrimitiveScorer ( G4VPrimitiveScorer ps)
inherited

Definition at line 166 of file G4VScoringMesh.cc.

167{
168 if(!ReadyForQuantity())
169 {
170 G4cerr << "ERROR : G4VScoringMesh::SetPrimitiveScorer() : "
171 << prs->GetName()
172 << " does not yet have mesh size or number of bins. Set them first."
173 << G4endl << "This Method is ignored." << G4endl;
174 return;
175 }
176 if(verboseLevel > 0)
177 G4cout << "G4VScoringMesh::SetPrimitiveScorer() : " << prs->GetName()
178 << " is registered."
179 << " 3D size: (" << fNSegment[0] << ", " << fNSegment[1] << ", "
180 << fNSegment[2] << ")" << G4endl;
181
182 prs->SetNijk(fNSegment[0], fNSegment[1], fNSegment[2]);
183 fCurrentPS = prs;
186 new G4THitsMap<G4StatDouble>(fWorldName, prs->GetName());
187 fMap[prs->GetName()] = map;
188}
G4bool RegisterPrimitive(G4VPrimitiveScorer *)
G4bool ReadyForQuantity() const

References G4VScoringMesh::fCurrentPS, G4VScoringMesh::fMap, G4VScoringMesh::fMFD, G4VScoringMesh::fNSegment, G4VScoringMesh::fWorldName, G4cerr, G4cout, G4endl, G4VPrimitiveScorer::GetName(), anonymous_namespace{G4QuasiElRatios.cc}::map, G4VScoringMesh::ReadyForQuantity(), G4MultiFunctionalDetector::RegisterPrimitive(), G4VPrimitiveScorer::SetNijk(), and G4VScoringMesh::verboseLevel.

Referenced by G4ScoreQuantityMessenger::SetNewValue().

◆ SetRMax()

void G4ScoringCylinder::SetRMax ( G4double  rMax)
inline

Definition at line 57 of file G4ScoringCylinder.hh.

57{ fSize[1] = rMax; }

References G4VScoringMesh::fSize.

◆ SetRMin()

void G4ScoringCylinder::SetRMin ( G4double  rMin)
inline

Definition at line 56 of file G4ScoringCylinder.hh.

56{ fSize[0] = rMin; }

References G4VScoringMesh::fSize.

◆ SetSize()

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

Definition at line 87 of file G4VScoringMesh.cc.

88{
89 if(!sizeIsSet)
90 {
91 sizeIsSet = true;
92 for(int i = 0; i < 3; i++)
93 {
94 fSize[i] = size[i];
95 }
96 }
97 else
98 {
99 G4String message = " Mesh size has already been set and it cannot be changed.\n";
100 message += " This method is ignored.";
101 G4Exception("G4VScoringMesh::SetSize()",
102 "DigiHitsUtilsScoreVScoringMesh000", JustWarning, message);
103 }
104}

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

Referenced by G4ScoringProbe::G4ScoringProbe(), G4ScoringRealWorld::G4ScoringRealWorld(), and G4ScoringMessenger::SetNewValue().

◆ SetupGeometry()

void G4ScoringCylinder::SetupGeometry ( G4VPhysicalVolume fWorldPhys)
protectedvirtual

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 G4LogicalVolume* tubsLogical = new G4LogicalVolume(tubsSolid, 0, tubsName);
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], 0, 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(0, 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], 0, 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,
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(0, 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, 0, 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 // } else if(fNSegment[IR] == 1) {
217 // if(verboseLevel > 9) G4cout << "G4ScoringCylinder::Construct() :
218 // Placement" << G4endl;
219 // G4cout<<"#@#@#@#@#@ G4ScoringCylinder::Construct() : Placement "<<G4endl;
220 // new G4PVPlacement(0, G4ThreeVector(0.,0.,0.), fMeshElementLogical,
221 // elementName, layerLogical[1], false, 0);
222 }
223 else
224 {
225 G4cerr << "G4ScoringCylinder::SetupGeometry() : "
226 << "invalid parameter (" << fNSegment[IR] << ") "
227 << "in mesh element placement." << G4endl;
228 }
229
230 // set the sensitive detector
232
233 // vis. attributes
234 G4VisAttributes* visatt = new G4VisAttributes(G4Colour(.5, .5, .5));
235 visatt->SetVisibility(true);
236 layerLogical[0]->SetVisAttributes(visatt);
237 layerLogical[1]->SetVisAttributes(visatt);
238 visatt = new G4VisAttributes(G4Colour(.5, .5, .5, 0.01));
239 // visatt->SetForceSolid(true);
241
242 if(verboseLevel > 9)
243 DumpVolumes();
244}
void SetVisAttributes(const G4VisAttributes *pVA)
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
static G4int GetReplicaLevel()
G4LogicalVolume * GetLogicalVolume() const
void SetVisibility(G4bool=true)
@ kPhi
Definition: geomdefs.hh:60
@ kZAxis
Definition: geomdefs.hh:57
@ kRho
Definition: geomdefs.hh:58

References DumpVolumes(), G4VScoringMesh::fAngle, G4VScoringMesh::fCenterPosition, G4VScoringMesh::fMeshElementLogical, G4VScoringMesh::fMFD, G4VScoringMesh::fNSegment, G4VScoringMesh::fRotationMatrix, G4VScoringMesh::fSize, G4VScoringMesh::fWorldName, G4cerr, G4cout, G4endl, G4VPhysicalVolume::GetLogicalVolume(), G4ScoringManager::GetReplicaLevel(), IPHI, IR, IZ, kPhi, kRho, kZAxis, G4LogicalVolume::SetSensitiveDetector(), G4LogicalVolume::SetVisAttributes(), G4VisAttributes::SetVisibility(), and G4VScoringMesh::verboseLevel.

◆ SetVerboseLevel()

void G4VScoringMesh::SetVerboseLevel ( G4int  vl)
inlineinherited

Definition at line 188 of file G4VScoringMesh.hh.

188{ verboseLevel = vl; }

References G4VScoringMesh::verboseLevel.

Referenced by G4ScoringManager::RegisterScoringMesh().

◆ SetZSize()

void G4ScoringCylinder::SetZSize ( G4double  zSize)
inline

Definition at line 58 of file G4ScoringCylinder.hh.

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

References G4VScoringMesh::fSize.

◆ WorkerConstruct()

void G4VScoringMesh::WorkerConstruct ( G4VPhysicalVolume fWorldPhys)
virtualinherited

Field Documentation

◆ copyNumberLevel

G4int G4VScoringMesh::copyNumberLevel
protectedinherited

◆ fActive

G4bool G4VScoringMesh::fActive
protectedinherited

Definition at line 202 of file G4VScoringMesh.hh.

Referenced by G4VScoringMesh::Activate(), and G4VScoringMesh::IsActive().

◆ fAngle

G4double G4VScoringMesh::fAngle[2]
protectedinherited

◆ fCenterPosition

G4ThreeVector G4VScoringMesh::fCenterPosition
protectedinherited

◆ fConstructed

G4bool G4VScoringMesh::fConstructed
protectedinherited

◆ fCurrentPS

G4VPrimitiveScorer* G4VScoringMesh::fCurrentPS
protectedinherited

◆ fDivisionAxisNames

G4String G4VScoringMesh::fDivisionAxisNames[3]
protectedinherited

◆ fDrawPSName

G4String G4VScoringMesh::fDrawPSName
protectedinherited

◆ fDrawUnit

G4String G4VScoringMesh::fDrawUnit
protectedinherited

◆ fDrawUnitValue

G4double G4VScoringMesh::fDrawUnitValue
protectedinherited

◆ fGeometryHasBeenDestroyed

G4bool G4VScoringMesh::fGeometryHasBeenDestroyed
protectedinherited

◆ fMap

MeshScoreMap G4VScoringMesh::fMap
protectedinherited

◆ fMeshElementLogical

G4LogicalVolume* G4VScoringMesh::fMeshElementLogical
protectedinherited

◆ fMFD

G4MultiFunctionalDetector* G4VScoringMesh::fMFD
protectedinherited

◆ fNSegment

G4int G4VScoringMesh::fNSegment[3]
protectedinherited

◆ fParallelWorldProcess

G4ParallelWorldProcess* G4VScoringMesh::fParallelWorldProcess
protectedinherited

◆ fRotationMatrix

G4RotationMatrix* G4VScoringMesh::fRotationMatrix
protectedinherited

◆ fShape

MeshShape G4VScoringMesh::fShape
protectedinherited

◆ fSize

G4double G4VScoringMesh::fSize[3]
protectedinherited

◆ fWorldName

G4String G4VScoringMesh::fWorldName
protectedinherited

◆ layeredMassFlg

G4bool G4VScoringMesh::layeredMassFlg
protectedinherited

◆ nMeshIsSet

G4bool G4VScoringMesh::nMeshIsSet
protectedinherited

◆ sizeIsSet

G4bool G4VScoringMesh::sizeIsSet
protectedinherited

◆ verboseLevel

G4int G4VScoringMesh::verboseLevel
protectedinherited

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