Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ScoringCylinder.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//
28
29#include "G4ScoringCylinder.hh"
30
31#include "G4LogicalVolume.hh"
33#include "G4PVDivision.hh"
34#include "G4PVPlacement.hh"
35#include "G4PVReplica.hh"
37#include "G4SDManager.hh"
38#include "G4ScoringManager.hh"
39#include "G4StatDouble.hh"
40#include "G4SystemOfUnits.hh"
41#include "G4Tubs.hh"
42#include "G4VPhysicalVolume.hh"
43#include "G4VPrimitiveScorer.hh"
44#include "G4VScoreColorMap.hh"
45#include "G4VVisManager.hh"
46#include "G4VisAttributes.hh"
47
57
59{
60 if(verboseLevel > 9)
61 G4cout << "G4ScoringCylinder::SetupGeometry() ..." << G4endl;
62
63 // World
64 G4VPhysicalVolume* scoringWorld = fWorldPhys;
65 G4LogicalVolume* worldLogical = scoringWorld->GetLogicalVolume();
66
67 // Scoring Mesh
68 if(verboseLevel > 9)
70 G4String tubsName = fWorldName + "_mesh";
71
72 if(verboseLevel > 9)
73 {
74 G4cout << "R min, R max., Dz =: " << fSize[0] << ", " << fSize[1]
75 << ", " << fSize[2] << G4endl;
76 }
77 G4VSolid* tubsSolid = new G4Tubs(tubsName + "0", // name
78 fSize[0], // R min
79 fSize[1], // R max
80 fSize[2], // Dz
81 fAngle[0], // starting phi
82 fAngle[1]); // segment phi
83 auto tubsLogical = new G4LogicalVolume(tubsSolid, nullptr, tubsName);
85 tubsName + "0", worldLogical, false, 0);
86
87 if(verboseLevel > 9)
88 G4cout << " # of segments : r, phi, z =: " << fNSegment[IR] << ", "
89 << fNSegment[IPHI] << ", " << fNSegment[IZ] << G4endl;
90
91 G4String layerName[2] = { tubsName + "1", tubsName + "2" };
92 G4VSolid* layerSolid[2];
93 G4LogicalVolume* layerLogical[2];
94
95 //-- fisrt nested layer (replicated along z direction)
96 if(verboseLevel > 9)
97 G4cout << "layer 1 :" << G4endl;
98 layerSolid[0] = new G4Tubs(layerName[0], // name
99 fSize[0], // inner radius
100 fSize[1], // outer radius
101 fSize[2] / fNSegment[IZ], // half len. in z
102 fAngle[0], // starting phi angle
103 fAngle[1]); // delta angle of the segment
104 layerLogical[0] = new G4LogicalVolume(layerSolid[0], nullptr, layerName[0]);
105 if(fNSegment[IZ] > 1)
106 {
107 if(verboseLevel > 9)
108 G4cout << "G4ScoringCylinder::Construct() : Replicate along z direction"
109 << G4endl;
111 {
112 if(verboseLevel > 9)
113 G4cout << "G4ScoringCylinder::Construct() : Replica" << G4endl;
114 new G4PVReplica(layerName[0], layerLogical[0], tubsLogical, kZAxis,
115 fNSegment[IZ], 2. * fSize[2] / fNSegment[IZ]);
116 }
117 else
118 {
119 if(verboseLevel > 9)
120 G4cout << "G4ScoringCylinder::Construct() : Division" << G4endl;
121 new G4PVDivision(layerName[0], layerLogical[0], tubsLogical, kZAxis,
122 fNSegment[IZ], 0.);
123 }
124 }
125 else if(fNSegment[IZ] == 1)
126 {
127 if(verboseLevel > 9)
128 G4cout << "G4ScoringCylinder::Construct() : Placement" << G4endl;
129 new G4PVPlacement(nullptr, G4ThreeVector(0., 0., 0.), layerLogical[0],
130 layerName[0], tubsLogical, false, 0);
131 }
132 else
133 {
134 G4cerr << "G4ScoringCylinder::SetupGeometry() : invalid parameter ("
135 << fNSegment[IZ] << ") "
136 << "in placement of the first nested layer." << G4endl;
137 }
138
139 // second nested layer (replicated along phi direction)
140 if(verboseLevel > 9)
141 G4cout << "layer 2 :" << G4endl;
142 layerSolid[1] =
143 new G4Tubs(layerName[1], fSize[0], fSize[1], fSize[2] / fNSegment[IZ],
144 fAngle[0], fAngle[1] / fNSegment[IPHI]);
145 layerLogical[1] = new G4LogicalVolume(layerSolid[1], nullptr, layerName[1]);
146 if(fNSegment[IPHI] > 1)
147 {
148 if(verboseLevel > 9)
149 G4cout << "G4ScoringCylinder::Construct() : Replicate along phi direction"
150 << G4endl;
152 {
153 if(verboseLevel > 9)
154 G4cout << "G4ScoringCylinder::Construct() : Replica" << G4endl;
155 new G4PVReplica(layerName[1], layerLogical[1], layerLogical[0], kPhi,
156 fNSegment[IPHI], fAngle[1] / fNSegment[IPHI], fAngle[0]);
157 }
158 else
159 {
160 if(verboseLevel > 9)
161 G4cout << "G4ScoringCylinder::Construct() : Division" << G4endl;
162 new G4PVDivision(layerName[1], layerLogical[1], layerLogical[0], kPhi,
163 fNSegment[IPHI], 0.);
164 }
165 }
166 else if(fNSegment[IPHI] == 1)
167 {
168 if(verboseLevel > 9)
169 G4cout << "G4ScoringCylinder::Construct() : Placement" << G4endl;
170 new G4PVPlacement(nullptr, G4ThreeVector(0., 0., 0.), layerLogical[1],
171 layerName[1], layerLogical[0], false, 0);
172 }
173 else
174 G4cerr << "ERROR : G4ScoringCylinder::SetupGeometry() : invalid parameter ("
175 << fNSegment[IPHI] << ") "
176 << "in placement of the second nested layer." << G4endl;
177
178 // mesh elements
179 if(verboseLevel > 9)
180 G4cout << "mesh elements :" << G4endl;
181 G4String elementName = tubsName + "3";
182 G4VSolid* elementSolid = new G4Tubs(
183 elementName, fSize[0], (fSize[1] - fSize[0]) / fNSegment[IR] + fSize[0],
184 fSize[2] / fNSegment[IZ], fAngle[0], fAngle[1] / fNSegment[IPHI]);
185 fMeshElementLogical = new G4LogicalVolume(elementSolid, nullptr, elementName);
186 if(fNSegment[IR] >= 1)
187 {
188 if(verboseLevel > 9)
189 G4cout << "G4ScoringCylinder::Construct() : Replicate along r direction"
190 << G4endl;
191
193 {
194 if(verboseLevel > 9)
195 G4cout << "G4ScoringCylinder::Construct() : Replica" << G4endl;
196 new G4PVReplica(elementName, fMeshElementLogical, layerLogical[1], kRho,
197 fNSegment[IR], (fSize[1] - fSize[0]) / fNSegment[IR],
198 fSize[0]);
199 }
200 else
201 {
202 if(verboseLevel > 9)
203 G4cout << "G4ScoringCylinder::Construct() : Division" << G4endl;
204 new G4PVDivision(elementName, fMeshElementLogical, layerLogical[1], kRho,
205 fNSegment[IR], 0.);
206 }
207 }
208 else
209 {
210 G4cerr << "G4ScoringCylinder::SetupGeometry() : "
211 << "invalid parameter (" << fNSegment[IR] << ") "
212 << "in mesh element placement." << G4endl;
213 }
214
215 // set the sensitive detector
217
218 // vis. attributes
219 auto visatt = new G4VisAttributes(G4Colour(.5, .5, .5));
220 visatt->SetVisibility(true);
221 layerLogical[0]->SetVisAttributes(visatt);
222 layerLogical[1]->SetVisAttributes(visatt);
223 visatt = new G4VisAttributes(G4Colour(.5, .5, .5, 0.01));
224 // visatt->SetForceSolid(true);
226
227 if(verboseLevel > 9)
228 DumpVolumes();
229}
230
232{
233 G4cout << "G4ScoringCylinder : " << fWorldName
234 << " --- Shape: Cylindrical mesh" << G4endl;
235
236 G4cout << " Size (Rmin, Rmax, Dz): (" << fSize[0] / cm << ", "
237 << fSize[1] / cm << ", " << fSize[2] / cm << ") [cm]" << G4endl;
238 G4cout << " Angle (start, span): (" << fAngle[0] / deg << ", "
239 << fAngle[1] / deg << ") [deg]" << G4endl;
240
242}
243
245 G4int axflg)
246{
248 if(pVisManager != nullptr)
249 {
250 // cell vectors
251 std::vector<double> ephi;
252 for(int phi = 0; phi < fNSegment[IPHI]; phi++)
253 ephi.push_back(0.);
254 //-
255 std::vector<std::vector<double>> zphicell; // zphicell[Z][PHI]
256 for(int z = 0; z < fNSegment[IZ]; z++)
257 zphicell.push_back(ephi);
258 //-
259 std::vector<std::vector<double>> rphicell; // rphicell[R][PHI]
260 for(int r = 0; r < fNSegment[IR]; r++)
261 rphicell.push_back(ephi);
262
263 // projections
264 G4int q[3];
265 auto itr = map->GetMap()->begin();
266 for(; itr != map->GetMap()->end(); itr++)
267 {
268 if(itr->first < 0)
269 {
270 G4cout << itr->first << G4endl;
271 continue;
272 }
273 GetRZPhi(itr->first, q);
274
275 zphicell[q[IZ]][q[IPHI]] += (itr->second->sum_wx()) / fDrawUnitValue;
276 rphicell[q[IR]][q[IPHI]] += (itr->second->sum_wx()) / fDrawUnitValue;
277 }
278
279 // search min./max. values
280 G4double zphimin = DBL_MAX, rphimin = DBL_MAX;
281 G4double zphimax = 0., rphimax = 0.;
282 for(int iphi = 0; iphi < fNSegment[IPHI]; iphi++)
283 {
284 for(int iz = 0; iz < fNSegment[IZ]; iz++)
285 {
286 if(zphimin > zphicell[iz][iphi])
287 zphimin = zphicell[iz][iphi];
288 if(zphimax < zphicell[iz][iphi])
289 zphimax = zphicell[iz][iphi];
290 }
291 for(int ir = 0; ir < fNSegment[IR]; ir++)
292 {
293 if(rphimin > rphicell[ir][iphi])
294 rphimin = rphicell[ir][iphi];
295 if(rphimax < rphicell[ir][iphi])
296 rphimax = rphicell[ir][iphi];
297 }
298 }
299
300 G4VisAttributes att;
301 att.SetForceSolid(true);
302 att.SetForceAuxEdgeVisible(true);
303
304 G4Scale3D scale;
305 if(axflg / 100 == 1)
306 {
307 // rz plane
308 }
309 axflg = axflg % 100;
310 if(axflg / 10 == 1)
311 {
312 pVisManager->BeginDraw();
313
314 // z-phi plane
315 if(colorMap->IfFloatMinMax())
316 {
317 colorMap->SetMinMax(zphimin, zphimax);
318 }
319
320 G4double zhalf = fSize[2] / fNSegment[IZ];
321 for(int phi = 0; phi < fNSegment[IPHI]; phi++)
322 {
323 for(int z = 0; z < fNSegment[IZ]; z++)
324 {
325 //-
326 G4double angle = fAngle[0] + fAngle[1] / fNSegment[IPHI] * phi;
327 G4double dphi = fAngle[1] / fNSegment[IPHI];
329 "z-phi", // name
330 fSize[1] * 0.99, fSize[1], // inner radius, outer radius
331 zhalf, // half length in z
332 angle, dphi * 0.99999); // starting phi angle, delta angle
333 //-
334 G4ThreeVector zpos(
335 0., 0., -fSize[2] + fSize[2] / fNSegment[IZ] * (1 + 2. * z));
336 G4Transform3D trans;
337 if(fRotationMatrix != nullptr)
338 {
339 trans =
340 G4Rotate3D(*fRotationMatrix).inverse() * G4Translate3D(zpos);
341 trans = G4Translate3D(fCenterPosition) * trans;
342 }
343 else
344 {
346 }
347 G4double c[4];
348 colorMap->GetMapColor(zphicell[z][phi], c);
349 att.SetColour(c[0], c[1], c[2]); //, c[3]);
350 //-
351 G4Polyhedron* poly = cylinder.GetPolyhedron();
352 poly->Transform(trans);
353 poly->SetVisAttributes(att);
354 pVisManager->Draw(*poly);
355 }
356 }
357 pVisManager->EndDraw();
358 }
359 axflg = axflg % 10;
360 if(axflg == 1)
361 {
362 pVisManager->BeginDraw();
363
364 // r-phi plane
365 if(colorMap->IfFloatMinMax())
366 {
367 colorMap->SetMinMax(rphimin, rphimax);
368 }
369
370 G4double rsize = (fSize[1] - fSize[0]) / fNSegment[IR];
371 for(int phi = 0; phi < fNSegment[IPHI]; phi++)
372 {
373 for(int r = 0; r < fNSegment[IR]; r++)
374 {
375 G4double rs[2] = { fSize[0] + rsize * r, fSize[0] + rsize * (r + 1) };
376 G4double angle = fAngle[0] + fAngle[1] / fNSegment[IPHI] * phi;
377 G4double dphi = fAngle[1] / fNSegment[IPHI];
378 G4Tubs cylindern("z-phi", rs[0], rs[1], 0.001, angle, dphi * 0.99999);
379 G4Tubs cylinderp = cylindern;
380
381 G4ThreeVector zposn(0., 0., -fSize[2]);
382 G4ThreeVector zposp(0., 0., fSize[2]);
383 G4Transform3D transn, transp;
384 if(fRotationMatrix != nullptr)
385 {
386 transn =
387 G4Rotate3D(*fRotationMatrix).inverse() * G4Translate3D(zposn);
388 transn = G4Translate3D(fCenterPosition) * transn;
389 transp =
390 G4Rotate3D(*fRotationMatrix).inverse() * G4Translate3D(zposp);
391 transp = G4Translate3D(fCenterPosition) * transp;
392 }
393 else
394 {
397 }
398 G4double c[4];
399 colorMap->GetMapColor(rphicell[r][phi], c);
400 att.SetColour(c[0], c[1], c[2]); //, c[3]);
401
402 G4Polyhedron* polyn = cylindern.GetPolyhedron();
403 polyn->Transform(transn);
404 polyn->SetVisAttributes(att);
405 pVisManager->Draw(*polyn);
406
407 G4Polyhedron* polyp = cylinderp.GetPolyhedron();
408 polyp->Transform(transp);
409 polyp->SetVisAttributes(att);
410 pVisManager->Draw(*polyp);
411 }
412 }
413
414 pVisManager->EndDraw();
415 }
416
417 colorMap->SetPSUnit(fDrawUnit);
418 colorMap->SetPSName(fDrawPSName);
419 colorMap->DrawColorChart();
420 }
421}
422
424 G4int idxProj, G4int idxColumn)
425{
426 G4int projAxis = 0;
427 switch(idxProj)
428 {
429 case 0:
430 projAxis = IR;
431 break;
432 case 1:
433 projAxis = IZ;
434 break;
435 case 2:
436 projAxis = IPHI;
437 break;
438 }
439
440 if(idxColumn < 0 || idxColumn >= fNSegment[projAxis])
441 {
442 G4cerr << "Warning : Column number " << idxColumn
443 << " is out of scoring mesh [0," << fNSegment[projAxis] - 1
444 << "]. Method ignored." << G4endl;
445 return;
446 }
448 if(pVisManager != nullptr)
449 {
450 // cell vectors
451 std::vector<std::vector<std::vector<double>>> cell; // cell[R][Z][PHI]
452 std::vector<double> ephi;
453 for(int phi = 0; phi < fNSegment[IPHI]; phi++)
454 ephi.push_back(0.);
455 std::vector<std::vector<double>> ezphi;
456 for(int z = 0; z < fNSegment[IZ]; z++)
457 ezphi.push_back(ephi);
458 for(int r = 0; r < fNSegment[IR]; r++)
459 cell.push_back(ezphi);
460
461 std::vector<std::vector<double>> rzcell; // rzcell[R][Z]
462 std::vector<double> ez;
463 for(int z = 0; z < fNSegment[IZ]; z++)
464 ez.push_back(0.);
465 for(int r = 0; r < fNSegment[IR]; r++)
466 rzcell.push_back(ez);
467
468 std::vector<std::vector<double>> zphicell; // zphicell[Z][PHI]
469 for(int z = 0; z < fNSegment[IZ]; z++)
470 zphicell.push_back(ephi);
471
472 std::vector<std::vector<double>> rphicell; // rphicell[R][PHI]
473 for(int r = 0; r < fNSegment[IR]; r++)
474 rphicell.push_back(ephi);
475
476 // projections
477 G4int q[3];
478 auto itr = map->GetMap()->begin();
479 for(; itr != map->GetMap()->end(); itr++)
480 {
481 if(itr->first < 0)
482 {
483 G4cout << itr->first << G4endl;
484 continue;
485 }
486 GetRZPhi(itr->first, q);
487
488 if(projAxis == IR && q[IR] == idxColumn)
489 { // zphi plane
490 zphicell[q[IZ]][q[IPHI]] += (itr->second->sum_wx()) / fDrawUnitValue;
491 }
492 if(projAxis == IZ && q[IZ] == idxColumn)
493 { // rphi plane
494 rphicell[q[IR]][q[IPHI]] += (itr->second->sum_wx()) / fDrawUnitValue;
495 }
496 if(projAxis == IPHI && q[IPHI] == idxColumn)
497 { // rz plane
498 rzcell[q[IR]][q[IZ]] += (itr->second->sum_wx()) / fDrawUnitValue;
499 }
500 }
501
502 // search min./max. values
503 G4double rzmin = DBL_MAX, zphimin = DBL_MAX, rphimin = DBL_MAX;
504 G4double rzmax = 0., zphimax = 0., rphimax = 0.;
505 for(int r = 0; r < fNSegment[IR]; r++)
506 {
507 for(int phi = 0; phi < fNSegment[IPHI]; phi++)
508 {
509 if(rphimin > rphicell[r][phi])
510 rphimin = rphicell[r][phi];
511 if(rphimax < rphicell[r][phi])
512 rphimax = rphicell[r][phi];
513 }
514 for(int z = 0; z < fNSegment[IZ]; z++)
515 {
516 if(rzmin > rzcell[r][z])
517 rzmin = rzcell[r][z];
518 if(rzmax < rzcell[r][z])
519 rzmax = rzcell[r][z];
520 }
521 }
522 for(int z = 0; z < fNSegment[IZ]; z++)
523 {
524 for(int phi = 0; phi < fNSegment[IPHI]; phi++)
525 {
526 if(zphimin > zphicell[z][phi])
527 zphimin = zphicell[z][phi];
528 if(zphimax < zphicell[z][phi])
529 zphimax = zphicell[z][phi];
530 }
531 }
532
533 G4VisAttributes att;
534 att.SetForceSolid(true);
535 att.SetForceAuxEdgeVisible(true);
536
537 pVisManager->BeginDraw();
538
539 G4Scale3D scale;
540 // z-phi plane
541 if(projAxis == IR)
542 {
543 if(colorMap->IfFloatMinMax())
544 {
545 colorMap->SetMinMax(zphimin, zphimax);
546 }
547
548 G4double zhalf = fSize[2] / fNSegment[IZ];
549 G4double rsize[2] = {
550 fSize[0] + (fSize[1] - fSize[0]) / fNSegment[IR] * idxColumn,
551 fSize[0] + (fSize[1] - fSize[0]) / fNSegment[IR] * (idxColumn + 1)
552 };
553 for(int phi = 0; phi < fNSegment[IPHI]; phi++)
554 {
555 for(int z = 0; z < fNSegment[IZ]; z++)
556 {
557 G4double angle = fAngle[0] + fAngle[1] / fNSegment[IPHI] * phi;
558 G4double dphi = fAngle[1] / fNSegment[IPHI];
559 G4Tubs cylinder("z-phi", rsize[0], rsize[1], zhalf, angle,
560 dphi * 0.99999);
561
562 G4ThreeVector zpos(
563 0., 0., -fSize[2] + fSize[2] / fNSegment[IZ] * (1 + 2. * z));
564 G4Transform3D trans;
565 if(fRotationMatrix != nullptr)
566 {
567 trans =
568 G4Rotate3D(*fRotationMatrix).inverse() * G4Translate3D(zpos);
569 trans = G4Translate3D(fCenterPosition) * trans;
570 }
571 else
572 {
574 }
575 G4double c[4];
576 colorMap->GetMapColor(zphicell[z][phi], c);
577 att.SetColour(c[0], c[1], c[2]); //, c[3]);
578
579 G4Polyhedron* poly = cylinder.GetPolyhedron();
580 poly->Transform(trans);
581 poly->SetVisAttributes(att);
582 pVisManager->Draw(*poly);
583 }
584 }
585
586 // r-phi plane
587 }
588 else if(projAxis == IZ)
589 {
590 if(colorMap->IfFloatMinMax())
591 {
592 colorMap->SetMinMax(rphimin, rphimax);
593 }
594
595 G4double rsize = (fSize[1] - fSize[0]) / fNSegment[IR];
596 for(int phi = 0; phi < fNSegment[IPHI]; phi++)
597 {
598 for(int r = 0; r < fNSegment[IR]; r++)
599 {
600 G4double rs[2] = { fSize[0] + rsize * r, fSize[0] + rsize * (r + 1) };
601 G4double angle = fAngle[0] + fAngle[1] / fNSegment[IPHI] * phi;
602 G4double dz = fSize[2] / fNSegment[IZ];
603 G4double dphi = fAngle[1] / fNSegment[IPHI];
604 G4Tubs cylinder("r-phi", rs[0], rs[1], dz, angle, dphi * 0.99999);
605 G4ThreeVector zpos(
606 0., 0., -fSize[2] + fSize[2] / fNSegment[IZ] * (idxColumn * 2 + 1));
607 G4Transform3D trans;
608 if(fRotationMatrix != nullptr)
609 {
610 trans =
611 G4Rotate3D(*fRotationMatrix).inverse() * G4Translate3D(zpos);
612 trans = G4Translate3D(fCenterPosition) * trans;
613 }
614 else
615 {
617 }
618 G4double c[4];
619 colorMap->GetMapColor(rphicell[r][phi], c);
620 att.SetColour(c[0], c[1], c[2]); //, c[3]);
621
622 G4Polyhedron* poly = cylinder.GetPolyhedron();
623 poly->Transform(trans);
624 poly->SetVisAttributes(att);
625 pVisManager->Draw(*poly);
626 }
627 }
628
629 // r-z plane
630 }
631 else if(projAxis == IPHI)
632 {
633 if(colorMap->IfFloatMinMax())
634 {
635 colorMap->SetMinMax(rzmin, rzmax);
636 }
637
638 G4double rsize = (fSize[1] - fSize[0]) / fNSegment[IR];
639 G4double zhalf = fSize[2] / fNSegment[IZ];
640 G4double angle = fAngle[0] + fAngle[1] / fNSegment[IPHI] * idxColumn;
641 G4double dphi = fAngle[1] / fNSegment[IPHI];
642 for(int z = 0; z < fNSegment[IZ]; z++)
643 {
644 for(int r = 0; r < fNSegment[IR]; r++)
645 {
646 G4double rs[2] = { fSize[0] + rsize * r, fSize[0] + rsize * (r + 1) };
647 G4Tubs cylinder("z-phi", rs[0], rs[1], zhalf, angle, dphi);
648
649 G4ThreeVector zpos(
650 0., 0., -fSize[2] + fSize[2] / fNSegment[IZ] * (2. * z + 1));
651 G4Transform3D trans;
652 if(fRotationMatrix != nullptr)
653 {
654 trans =
655 G4Rotate3D(*fRotationMatrix).inverse() * G4Translate3D(zpos);
656 trans = G4Translate3D(fCenterPosition) * trans;
657 }
658 else
659 {
661 }
662 G4double c[4];
663 colorMap->GetMapColor(rzcell[r][z], c);
664 att.SetColour(c[0], c[1], c[2]); //, c[3]);
665
666 G4Polyhedron* poly = cylinder.GetPolyhedron();
667 poly->Transform(trans);
668 poly->SetVisAttributes(att);
669 pVisManager->Draw(*poly);
670 }
671 }
672 }
673 pVisManager->EndDraw();
674 }
675
676 colorMap->SetPSUnit(fDrawUnit);
677 colorMap->SetPSName(fDrawPSName);
678 colorMap->DrawColorChart();
679}
680
682{
683 // index = k + j * k-size + i * jk-plane-size
684
685 // nested : z -> phi -> r
686 G4int i = IZ;
687 G4int j = IPHI;
688 G4int k = IR;
689 G4int jk = fNSegment[j] * fNSegment[k];
690 q[i] = index / jk;
691 q[j] = (index - q[i] * jk) / fNSegment[k];
692 q[k] = index - q[j] * fNSegment[k] - q[i] * jk;
693}
694
696#include "G4Material.hh"
698#include "G4SolidStore.hh"
699#include "G4UnitsTable.hh"
701#include "G4VSolid.hh"
702
703void G4ScoringCylinder::DumpVolumes()
704{
705 G4int lvl = 2;
706 DumpSolids(lvl);
707 DumpLogVols(lvl);
708 DumpPhysVols(lvl);
709}
710
711void G4ScoringCylinder::DumpSolids(G4int lvl)
712{
713 G4cout << "*********** List of registered solids *************" << G4endl;
714 auto store = G4SolidStore::GetInstance();
715 auto itr = store->begin();
716 for(; itr != store->end(); itr++)
717 {
718 switch(lvl)
719 {
720 case 0:
721 G4cout << (*itr)->GetName() << G4endl;
722 break;
723 case 1:
724 G4cout << (*itr)->GetName() << "\t volume = "
725 << G4BestUnit((*itr)->GetCubicVolume(), "Volume")
726 << "\t surface = "
727 << G4BestUnit((*itr)->GetSurfaceArea(), "Surface") << G4endl;
728 break;
729 default:
730 (*itr)->DumpInfo();
731 break;
732 }
733 }
734}
735
736void G4ScoringCylinder::DumpLogVols(G4int lvl)
737{
738 G4cout << "*********** List of registered logical volumes *************"
739 << G4endl;
741 auto itr = store->begin();
742 for(; itr != store->end(); itr++)
743 {
744 G4cout << (*itr)->GetName()
745 << "\t Solid = " << (*itr)->GetSolid()->GetName();
746 if((*itr)->GetMaterial() != nullptr)
747 {
748 G4cout << "\t Material = " << (*itr)->GetMaterial()->GetName() << G4endl;
749 }
750 else
751 {
752 G4cout << "\t Material : not defined " << G4endl;
753 }
754 if(lvl < 1)
755 continue;
756 G4cout << "\t region = ";
757 if((*itr)->GetRegion() != nullptr)
758 {
759 G4cout << (*itr)->GetRegion()->GetName();
760 }
761 else
762 {
763 G4cout << "not defined";
764 }
765 G4cout << "\t sensitive detector = ";
766 if((*itr)->GetSensitiveDetector() != nullptr)
767 {
768 G4cout << (*itr)->GetSensitiveDetector()->GetName();
769 }
770 else
771 {
772 G4cout << "not defined";
773 }
774 G4cout << G4endl;
775 G4cout << "\t daughters = " << (*itr)->GetNoDaughters();
776 if((*itr)->GetNoDaughters() > 0)
777 {
778 switch((*itr)->CharacteriseDaughters())
779 {
780 case kNormal:
781 G4cout << " (placement)";
782 break;
783 case kReplica:
784 G4cout << " (replica : " << (*itr)->GetDaughter(0)->GetMultiplicity()
785 << ")";
786 break;
787 case kParameterised:
788 G4cout << " (parameterized : "
789 << (*itr)->GetDaughter(0)->GetMultiplicity() << ")";
790 break;
791 default:;
792 }
793 }
794 G4cout << G4endl;
795 if(lvl < 2)
796 continue;
797 if((*itr)->GetMaterial() != nullptr)
798 {
799 G4cout << "\t weight = " << G4BestUnit((*itr)->GetMass(), "Mass")
800 << G4endl;
801 }
802 else
803 {
804 G4cout << "\t weight : not available" << G4endl;
805 }
806 }
807}
808
809void G4ScoringCylinder::DumpPhysVols(G4int lvl)
810{
811 G4cout << "*********** List of registered physical volumes *************"
812 << G4endl;
814 auto itr = store->begin();
815 for(; itr != store->end(); itr++)
816 {
817 switch(lvl)
818 {
819 case 0:
820 G4cout << (*itr)->GetName() << G4endl;
821 break;
822 case 1:
823 G4cout << (*itr)->GetName() << "\t logical volume = "
824 << (*itr)->GetLogicalVolume()->GetName()
825 << "\t mother logical = ";
826 if((*itr)->GetMotherLogical() != nullptr)
827 {
828 G4cout << (*itr)->GetMotherLogical()->GetName();
829 }
830 else
831 {
832 G4cout << "not defined";
833 }
834 G4cout << G4endl;
835 break;
836 default:
837 G4cout << (*itr)->GetName() << "\t logical volume = "
838 << (*itr)->GetLogicalVolume()->GetName()
839 << "\t mother logical = ";
840 if((*itr)->GetMotherLogical() != nullptr)
841 {
842 G4cout << (*itr)->GetMotherLogical()->GetName();
843 }
844 else
845 {
846 G4cout << "not defined";
847 }
848 G4cout << "\t type = ";
849 switch((*itr)->VolumeType())
850 {
851 case kNormal:
852 G4cout << "placement";
853 break;
854 case kReplica:
855 G4cout << "replica";
856 break;
857 case kParameterised:
858 G4cout << "parameterized";
859 break;
860 default:;
861 }
862 G4cout << G4endl;
863 }
864 }
865}
#define G4BestUnit(a, b)
CLHEP::Hep3Vector G4ThreeVector
HepGeom::Translate3D G4Translate3D
HepGeom::Rotate3D G4Rotate3D
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4Polyhedron * GetPolyhedron() const override
static G4LogicalVolumeStore * GetInstance()
void SetVisAttributes(const G4VisAttributes *pVA)
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
static G4PhysicalVolumeStore * GetInstance()
void GetRZPhi(G4int index, G4int q[3]) const
void SetupGeometry(G4VPhysicalVolume *fWorldPhys) override
void Draw(RunScore *map, G4VScoreColorMap *colorMap, G4int axflg=111) override
void List() const override
G4ScoringCylinder(G4String wName)
void DrawColumn(RunScore *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn) override
static G4int GetReplicaLevel()
static G4SolidStore * GetInstance()
G4LogicalVolume * GetLogicalVolume() const
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
virtual void List() const
G4double fDrawUnitValue
G4MultiFunctionalDetector * fMFD
G4String fDivisionAxisNames[3]
G4LogicalVolume * fMeshElementLogical
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)
@ kPhi
Definition geomdefs.hh:60
@ kZAxis
Definition geomdefs.hh:57
@ kRho
Definition geomdefs.hh:58
@ kNormal
Definition geomdefs.hh:84
@ kParameterised
Definition geomdefs.hh:86
@ kReplica
Definition geomdefs.hh:85
#define DBL_MAX
Definition templates.hh:62