Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 // $Id: G4ScoringCylinder.cc 99154 2016-09-07 08:06:30Z gcosmo $
28 //
29 
30 #include "G4ScoringCylinder.hh"
31 
32 #include "G4PhysicalConstants.hh"
33 #include "G4SystemOfUnits.hh"
34 #include "G4VPhysicalVolume.hh"
35 #include "G4Tubs.hh"
36 #include "G4LogicalVolume.hh"
37 #include "G4VPhysicalVolume.hh"
38 #include "G4PVPlacement.hh"
39 #include "G4PVReplica.hh"
40 #include "G4PVDivision.hh"
41 #include "G4VisAttributes.hh"
42 #include "G4VVisManager.hh"
43 #include "G4VScoreColorMap.hh"
44 
45 #include "G4SDManager.hh"
47 #include "G4SDParticleFilter.hh"
48 #include "G4VPrimitiveScorer.hh"
49 #include "G4PSEnergyDeposit.hh"
50 #include "G4PSTrackLength.hh"
51 #include "G4PSNofStep.hh"
52 #include "G4ScoringManager.hh"
53 #include "G4StatDouble.hh"
54 
55 
57  :G4VScoringMesh(wName)
58 {
60 
61  fDivisionAxisNames[0] = "Z";
62  fDivisionAxisNames[1] = "PHI";
63  fDivisionAxisNames[2] = "R";
64 }
65 
67 {;}
68 
70 
71  if(verboseLevel > 9) G4cout << "G4ScoringCylinder::SetupGeometry() ..." << G4endl;
72 
73  // World
74  G4VPhysicalVolume * scoringWorld = fWorldPhys;
75  G4LogicalVolume * worldLogical = scoringWorld->GetLogicalVolume();
76 
77  // Scoring Mesh
78  if(verboseLevel > 9) G4cout << fWorldName << G4endl;
79  G4String tubsName = fWorldName+"_mesh";
80 
81  if(verboseLevel > 9) G4cout << "R max., Dz =: " << fSize[0] << ", " << fSize[1] << G4endl;
82  G4VSolid * tubsSolid = new G4Tubs(tubsName+"0", // name
83  0., // R min
84  fSize[0], // R max
85  fSize[1], // Dz
86  0., // starting phi
87  twopi*rad); // segment phi
88  G4LogicalVolume * tubsLogical = new G4LogicalVolume(tubsSolid, 0, tubsName);
90  tubsLogical, tubsName+"0", worldLogical, false, 0);
91 
92  if(verboseLevel > 9) G4cout << " # of segments : r, phi, z =: "
93  << fNSegment[IR] << ", " << fNSegment[IPHI] << ", " << fNSegment[IZ] << G4endl;
94 
95  G4String layerName[2] = {tubsName + "1", tubsName + "2"};
96  G4VSolid * layerSolid[2];
97  G4LogicalVolume * layerLogical[2];
98 
99  //-- fisrt nested layer (replicated along z direction)
100  if(verboseLevel > 9) G4cout << "layer 1 :" << G4endl;
101  layerSolid[0] = new G4Tubs(layerName[0], // name
102  0., // inner radius
103  fSize[0], // outer radius
104  fSize[1]/fNSegment[IZ], // half len. in z
105  0., // starting phi angle
106  twopi*rad); // delta angle of the segment
107  layerLogical[0] = new G4LogicalVolume(layerSolid[0], 0, layerName[0]);
108  if(fNSegment[IZ] > 1) {
109  if(verboseLevel > 9) G4cout << "G4ScoringCylinder::Construct() : Replicate along z direction" << G4endl;
111  if(verboseLevel > 9) G4cout << "G4ScoringCylinder::Construct() : Replica" << G4endl;
112  new G4PVReplica(layerName[0], layerLogical[0], tubsLogical, kZAxis, fNSegment[IZ], 2.*fSize[1]/fNSegment[IZ]);
113  } else {
114  if(verboseLevel > 9) G4cout << "G4ScoringCylinder::Construct() : Division" << G4endl;
115  new G4PVDivision(layerName[0], layerLogical[0], tubsLogical, kZAxis, fNSegment[IZ], 0.);
116  }
117  } else if(fNSegment[IZ] == 1) {
118  if(verboseLevel > 9) G4cout << "G4ScoringCylinder::Construct() : Placement" << G4endl;
119  new G4PVPlacement(0, G4ThreeVector(0.,0.,0.), layerLogical[0], layerName[0], tubsLogical, false, 0);
120  } else {
121  G4cerr << "G4ScoringCylinder::SetupGeometry() : invalid parameter ("
122  << fNSegment[IZ] << ") "
123  << "in placement of the first nested layer." << G4endl;
124  }
125 
126  // second nested layer (replicated along phi direction)
127  if(verboseLevel > 9) G4cout << "layer 2 :" << G4endl;
128  layerSolid[1] = new G4Tubs(layerName[1],
129  0.,
130  fSize[0],
131  fSize[1]/fNSegment[IZ],
132  0.,
133  twopi*rad/fNSegment[IPHI]);
134  layerLogical[1] = new G4LogicalVolume(layerSolid[1], 0, layerName[1]);
135  if(fNSegment[IPHI] > 1) {
136  if(verboseLevel > 9) G4cout << "G4ScoringCylinder::Construct() : Replicate along phi direction" << G4endl;
138  if(verboseLevel > 9) G4cout << "G4ScoringCylinder::Construct() : Replica" << G4endl;
139  new G4PVReplica(layerName[1], layerLogical[1], layerLogical[0], kPhi,
140  fNSegment[IPHI], twopi*rad/fNSegment[IPHI]);
141  } else {
142  if(verboseLevel > 9) G4cout << "G4ScoringCylinder::Construct() : Division" << G4endl;
143  new G4PVDivision(layerName[1], layerLogical[1], layerLogical[0], kPhi, fNSegment[IPHI], 0.);
144  }
145  } else if(fNSegment[IPHI] == 1) {
146  if(verboseLevel > 9) G4cout << "G4ScoringCylinder::Construct() : Placement" << G4endl;
147  new G4PVPlacement(0, G4ThreeVector(0.,0.,0.), layerLogical[1], layerName[1], layerLogical[0], false, 0);
148  } else
149  G4cerr << "ERROR : G4ScoringCylinder::SetupGeometry() : invalid parameter ("
150  << fNSegment[IPHI] << ") "
151  << "in placement of the second nested layer." << G4endl;
152 
153  // mesh elements
154  if(verboseLevel > 9) G4cout << "mesh elements :" << G4endl;
155  G4String elementName = tubsName +"3";
156  G4VSolid * elementSolid = new G4Tubs(elementName,
157  0.,
158  fSize[0]/fNSegment[IR],
159  fSize[1]/fNSegment[IZ],
160  0.,
161  twopi*rad/fNSegment[IPHI]);
162  fMeshElementLogical = new G4LogicalVolume(elementSolid, 0, elementName);
163  if(fNSegment[IR] > 1) {
164 
165  if(verboseLevel > 9) G4cout << "G4ScoringCylinder::Construct() : Replicate along r direction" << G4endl;
166 
168  if(verboseLevel > 9) G4cout << "G4ScoringCylinder::Construct() : Replica" << G4endl;
169  new G4PVReplica(elementName, fMeshElementLogical, layerLogical[1], kRho,
170  fNSegment[IR], fSize[0]/fNSegment[IR]);
171  } else {
172  if(verboseLevel > 9) G4cout << "G4ScoringCylinder::Construct() : Division" << G4endl;
173  new G4PVDivision(elementName, fMeshElementLogical, layerLogical[1], kRho, fNSegment[IR], 0.);
174  }
175  } else if(fNSegment[IR] == 1) {
176  if(verboseLevel > 9) G4cout << "G4ScoringCylinder::Construct() : Placement" << G4endl;
177  new G4PVPlacement(0, G4ThreeVector(0.,0.,0.), fMeshElementLogical, elementName, layerLogical[1], false, 0);
178  } else {
179  G4cerr << "G4ScoringCylinder::SetupGeometry() : "
180  << "invalid parameter (" << fNSegment[IR] << ") "
181  << "in mesh element placement." << G4endl;
182  }
183 
184  // set the sensitive detector
186 
187 
188  // vis. attributes
189  G4VisAttributes * visatt = new G4VisAttributes(G4Colour(.5,.5,.5));
190  visatt->SetVisibility(true);
191  layerLogical[0]->SetVisAttributes(visatt);
192  layerLogical[1]->SetVisAttributes(visatt);
193  visatt = new G4VisAttributes(G4Colour(.5,.5,.5,0.01));
194  //visatt->SetForceSolid(true);
196 }
197 
199  G4cout << "G4ScoringCylinder : " << fWorldName << " --- Shape: Cylindrical mesh" << G4endl;
200 
201  G4cout << " Size (R, Dz): ("
202  << fSize[0]/cm << ", "
203  << fSize[1]/cm << ") [cm]"
204  << G4endl;
205 
207 }
208 
209 
211  G4VScoreColorMap* colorMap, G4int axflg) {
212 
214  if(pVisManager) {
215 
216  // cell vectors
217  std::vector<double> ephi;
218  for(int phi = 0; phi < fNSegment[IPHI]; phi++) ephi.push_back(0.);
219  //-
220  std::vector<std::vector<double> > zphicell; // zphicell[Z][PHI]
221  for(int z = 0; z < fNSegment[IZ]; z++) zphicell.push_back(ephi);
222  //-
223  std::vector<std::vector<double> > rphicell; // rphicell[R][PHI]
224  for(int r = 0; r < fNSegment[IR]; r++) rphicell.push_back(ephi);
225 
226  // projections
227  G4int q[3];
228  std::map<G4int, G4StatDouble*>::iterator itr = map->GetMap()->begin();
229  for(; itr != map->GetMap()->end(); itr++) {
230  if(itr->first < 0) {
231  G4cout << itr->first << G4endl;
232  continue;
233  }
234  GetRZPhi(itr->first, q);
235 
236  zphicell[q[IZ]][q[IPHI]] += (itr->second->sum_wx())/fDrawUnitValue;
237  rphicell[q[IR]][q[IPHI]] += (itr->second->sum_wx())/fDrawUnitValue;
238  }
239 
240  // search min./max. values
241  G4double zphimin = DBL_MAX, rphimin = DBL_MAX;
242  G4double zphimax = 0., rphimax = 0.;
243  for(int iphi = 0; iphi < fNSegment[IPHI]; iphi++) {
244  for(int iz = 0; iz < fNSegment[IZ]; iz++) {
245  if(zphimin > zphicell[iz][iphi]) zphimin = zphicell[iz][iphi];
246  if(zphimax < zphicell[iz][iphi]) zphimax = zphicell[iz][iphi];
247  }
248  for(int ir = 0; ir < fNSegment[IR]; ir++) {
249  if(rphimin > rphicell[ir][iphi]) rphimin = rphicell[ir][iphi];
250  if(rphimax < rphicell[ir][iphi]) rphimax = rphicell[ir][iphi];
251  }
252  }
253 
254  G4VisAttributes att;
255  att.SetForceSolid(true);
256  att.SetForceAuxEdgeVisible(true);
257 
258 
259  G4Scale3D scale;
260  if(axflg/100==1) {
261  // rz plane
262  }
263  axflg = axflg%100;
264  if(axflg/10==1) {
265  pVisManager->BeginDraw();
266 
267  // z-phi plane
268  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(zphimin, zphimax); }
269 
270  G4double zhalf = fSize[1]/fNSegment[IZ];
271  for(int phi = 0; phi < fNSegment[IPHI]; phi++) {
272  for(int z = 0; z < fNSegment[IZ]; z++) {
273  //-
274  G4double angle = twopi/fNSegment[IPHI]*phi;
275  G4double dphi = twopi/fNSegment[IPHI];
276  G4Tubs cylinder("z-phi", // name
277  fSize[0]*0.99, fSize[0], // inner radius, outer radius
278  zhalf, // half length in z
279  angle, dphi*0.99999); // starting phi angle, delta angle
280  //-
281  G4ThreeVector zpos(0., 0., -fSize[1] + fSize[1]/fNSegment[IZ]*(1 + 2.*z));
282  G4Transform3D trans;
283  if(fRotationMatrix) {
284  trans = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(zpos);
285  trans = G4Translate3D(fCenterPosition)*trans;
286  } else {
288  }
289  G4double c[4];
290  colorMap->GetMapColor(zphicell[z][phi], c);
291  att.SetColour(c[0], c[1], c[2]);//, c[3]);
292  //-
293  G4Polyhedron * poly = cylinder.GetPolyhedron();
294  poly->Transform(trans);
295  poly->SetVisAttributes(att);
296  pVisManager->Draw(*poly);
297  }
298  }
299  pVisManager->EndDraw();
300  }
301  axflg = axflg%10;
302  if(axflg==1) {
303  pVisManager->BeginDraw();
304 
305  // r-phi plane
306  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(rphimin, rphimax); }
307 
308  G4double rsize = fSize[0]/fNSegment[IR];
309  for(int phi = 0; phi < fNSegment[IPHI]; phi++) {
310  for(int r = 0; r < fNSegment[IR]; r++) {
311 
312  G4double rs[2] = {rsize*r, rsize*(r+1)};
313  G4double angle = twopi/fNSegment[IPHI]*phi;
314  G4double dphi = twopi/fNSegment[IPHI];
315  G4Tubs cylindern("z-phi", rs[0], rs[1], 0.001,
316  angle, dphi*0.99999);
317  G4Tubs cylinderp = cylindern;
318 
319  G4ThreeVector zposn(0., 0., -fSize[1]);
320  G4ThreeVector zposp(0., 0., fSize[1]);
321  G4Transform3D transn, transp;
322  if(fRotationMatrix) {
323  transn = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(zposn);
324  transn = G4Translate3D(fCenterPosition)*transn;
325  transp = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(zposp);
326  transp = G4Translate3D(fCenterPosition)*transp;
327  } else {
330  }
331  G4double c[4];
332  colorMap->GetMapColor(rphicell[r][phi], c);
333  att.SetColour(c[0], c[1], c[2]);//, c[3]);
334 
335  G4Polyhedron * polyn = cylindern.GetPolyhedron();
336  polyn->Transform(transn);
337  polyn->SetVisAttributes(att);
338  pVisManager->Draw(*polyn);
339 
340  G4Polyhedron * polyp = cylinderp.GetPolyhedron();
341  polyp->Transform(transp);
342  polyp->SetVisAttributes(att);
343  pVisManager->Draw(*polyp);
344  }
345  }
346 
347  pVisManager->EndDraw();
348  }
349 
350  colorMap->SetPSUnit(fDrawUnit);
351  colorMap->SetPSName(fDrawPSName);
352  colorMap->DrawColorChart();
353 
354  }
355 }
356 
358  G4int idxProj, G4int idxColumn)
359 {
360  G4int projAxis = 0;
361  switch(idxProj) {
362  case 0:
363  projAxis = IR;
364  break;
365  case 1:
366  projAxis = IZ;
367  break;
368  case 2:
369  projAxis = IPHI;
370  break;
371  }
372 
373  if(idxColumn<0 || idxColumn>=fNSegment[projAxis])
374  {
375  G4cerr << "Warning : Column number " << idxColumn << " is out of scoring mesh [0," << fNSegment[projAxis]-1 <<
376  "]. Method ignored." << G4endl;
377  return;
378  }
380  if(pVisManager) {
381 
382  // cell vectors
383  std::vector<std::vector<std::vector<double> > > cell; // cell[R][Z][PHI]
384  std::vector<double> ephi;
385  for(int phi = 0; phi < fNSegment[IPHI]; phi++) ephi.push_back(0.);
386  std::vector<std::vector<double> > ezphi;
387  for(int z = 0; z < fNSegment[IZ]; z++) ezphi.push_back(ephi);
388  for(int r = 0; r < fNSegment[IR]; r++) cell.push_back(ezphi);
389 
390  std::vector<std::vector<double> > rzcell; // rzcell[R][Z]
391  std::vector<double> ez;
392  for(int z = 0; z < fNSegment[IZ]; z++) ez.push_back(0.);
393  for(int r = 0; r < fNSegment[IR]; r++) rzcell.push_back(ez);
394 
395  std::vector<std::vector<double> > zphicell; // zphicell[Z][PHI]
396  for(int z = 0; z < fNSegment[IZ]; z++) zphicell.push_back(ephi);
397 
398  std::vector<std::vector<double> > rphicell; // rphicell[R][PHI]
399  for(int r = 0; r < fNSegment[IR]; r++) rphicell.push_back(ephi);
400 
401  // projections
402  G4int q[3];
403  std::map<G4int,G4StatDouble*>::iterator itr = map->GetMap()->begin();
404  for(; itr != map->GetMap()->end(); itr++) {
405  if(itr->first < 0) {
406  G4cout << itr->first << G4endl;
407  continue;
408  }
409  GetRZPhi(itr->first, q);
410 
411  if(projAxis == IR && q[IR] == idxColumn) { // zphi plane
412  zphicell[q[IZ]][q[IPHI]] += (itr->second->sum_wx())/fDrawUnitValue;
413  }
414  if(projAxis == IZ && q[IZ] == idxColumn) { // rphi plane
415  rphicell[q[IR]][q[IPHI]] += (itr->second->sum_wx())/fDrawUnitValue;
416  }
417  if(projAxis == IPHI && q[IPHI] == idxColumn) { // rz plane
418  rzcell[q[IR]][q[IZ]] += (itr->second->sum_wx())/fDrawUnitValue;
419  }
420  }
421 
422  // search min./max. values
423  G4double rzmin = DBL_MAX, zphimin = DBL_MAX, rphimin = DBL_MAX;
424  G4double rzmax = 0., zphimax = 0., rphimax = 0.;
425  for(int r = 0; r < fNSegment[IR]; r++) {
426  for(int phi = 0; phi < fNSegment[IPHI]; phi++) {
427  if(rphimin > rphicell[r][phi]) rphimin = rphicell[r][phi];
428  if(rphimax < rphicell[r][phi]) rphimax = rphicell[r][phi];
429  }
430  for(int z = 0; z < fNSegment[IZ]; z++) {
431  if(rzmin > rzcell[r][z]) rzmin = rzcell[r][z];
432  if(rzmax < rzcell[r][z]) rzmax = rzcell[r][z];
433  }
434  }
435  for(int z = 0; z < fNSegment[IZ]; z++) {
436  for(int phi = 0; phi < fNSegment[IPHI]; phi++) {
437  if(zphimin > zphicell[z][phi]) zphimin = zphicell[z][phi];
438  if(zphimax < zphicell[z][phi]) zphimax = zphicell[z][phi];
439  }
440  }
441 
442 
443  G4VisAttributes att;
444  att.SetForceSolid(true);
445  att.SetForceAuxEdgeVisible(true);
446 
447  pVisManager->BeginDraw();
448 
449  G4Scale3D scale;
450  // z-phi plane
451  if(projAxis == IR) {
452 
453  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(zphimin,zphimax); }
454 
455  G4double zhalf = fSize[1]/fNSegment[IZ];
456  G4double rsize[2] = {fSize[0]/fNSegment[IR]*idxColumn,
457  fSize[0]/fNSegment[IR]*(idxColumn+1)};
458  for(int phi = 0; phi < fNSegment[IPHI]; phi++) {
459  for(int z = 0; z < fNSegment[IZ]; z++) {
460 
461  G4double angle = twopi/fNSegment[IPHI]*phi*radian;
462  G4double dphi = twopi/fNSegment[IPHI]*radian;
463  G4Tubs cylinder("z-phi", rsize[0], rsize[1], zhalf,
464  angle, dphi*0.99999);
465 
466  G4ThreeVector zpos(0., 0., -fSize[1] + fSize[1]/fNSegment[IZ]*(1 + 2.*z));
467  G4Transform3D trans;
468  if(fRotationMatrix) {
469  trans = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(zpos);
470  trans = G4Translate3D(fCenterPosition)*trans;
471  } else {
473  }
474  G4double c[4];
475  colorMap->GetMapColor(zphicell[z][phi], c);
476  att.SetColour(c[0], c[1], c[2]);//, c[3]);
477 
478  G4Polyhedron * poly = cylinder.GetPolyhedron();
479  poly->Transform(trans);
480  poly->SetVisAttributes(att);
481  pVisManager->Draw(*poly);
482  }
483  }
484 
485  // r-phi plane
486  } else if(projAxis == IZ) {
487  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(rphimin,rphimax); }
488 
489  G4double rsize = fSize[0]/fNSegment[IR];
490  for(int phi = 0; phi < fNSegment[IPHI]; phi++) {
491  for(int r = 0; r < fNSegment[IR]; r++) {
492 
493  G4double rs[2] = {rsize*r, rsize*(r+1)};
494  G4double angle = twopi/fNSegment[IPHI]*phi*radian;
495  G4double dz = fSize[1]/fNSegment[IZ];
496  G4double dphi = twopi/fNSegment[IPHI]*radian;
497  G4Tubs cylinder("r-phi", rs[0], rs[1], dz,
498  angle, dphi*0.99999);
499  G4ThreeVector zpos(0., 0.,
500  -fSize[1]+fSize[1]/fNSegment[IZ]*(idxColumn*2+1));
501  G4Transform3D trans;
502  if(fRotationMatrix) {
503  trans = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(zpos);
504  trans = G4Translate3D(fCenterPosition)*trans;
505  } else {
507  }
508  G4double c[4];
509  colorMap->GetMapColor(rphicell[r][phi], c);
510  att.SetColour(c[0], c[1], c[2]);//, c[3]);
511 
512  G4Polyhedron * poly = cylinder.GetPolyhedron();
513  poly->Transform(trans);
514  poly->SetVisAttributes(att);
515  pVisManager->Draw(*poly);
516  }
517  }
518 
519  // r-z plane
520  } else if(projAxis == IPHI) {
521  if(colorMap->IfFloatMinMax()) { colorMap->SetMinMax(rzmin,rzmax); }
522 
523  G4double rsize = fSize[0]/fNSegment[IR];
524  G4double zhalf = fSize[1]/fNSegment[IZ];
525  G4double angle = twopi/fNSegment[IPHI]*idxColumn*radian;
526  G4double dphi = twopi/fNSegment[IPHI]*radian;
527  for(int z = 0; z < fNSegment[IZ]; z++) {
528  for(int r = 0; r < fNSegment[IR]; r++) {
529 
530  G4double rs[2] = {rsize*r, rsize*(r+1)};
531  G4Tubs cylinder("z-phi", rs[0], rs[1], zhalf,
532  angle, dphi);
533 
534  G4ThreeVector zpos(0., 0.,
535  -fSize[1]+fSize[1]/fNSegment[IZ]*(2.*z+1));
536  G4Transform3D trans;
537  if(fRotationMatrix) {
538  trans = G4Rotate3D(*fRotationMatrix).inverse()*G4Translate3D(zpos);
539  trans = G4Translate3D(fCenterPosition)*trans;
540  } else {
542  }
543  G4double c[4];
544  colorMap->GetMapColor(rzcell[r][z], c);
545  att.SetColour(c[0], c[1], c[2]);//, c[3]);
546 
547  G4Polyhedron * poly = cylinder.GetPolyhedron();
548  poly->Transform(trans);
549  poly->SetVisAttributes(att);
550  pVisManager->Draw(*poly);
551  }
552  }
553  }
554  pVisManager->EndDraw();
555  }
556 
557  colorMap->SetPSUnit(fDrawUnit);
558  colorMap->SetPSName(fDrawPSName);
559  colorMap->DrawColorChart();
560 
561 }
562 
563 void G4ScoringCylinder::GetRZPhi(G4int index, G4int q[3]) const {
564  // index = k + j * k-size + i * jk-plane-size
565 
566  // nested : z -> phi -> r
567  G4int i = IZ;
568  G4int j = IPHI;
569  G4int k = IR;
570  G4int jk = fNSegment[j]*fNSegment[k];
571  q[i] = index/jk;
572  q[j] = (index - q[i]*jk)/fNSegment[k];
573  q[k] = index - q[j]*fNSegment[k] - q[i]*jk;
574 }
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
Definition: geomdefs.hh:54
void SetColour(const G4Colour &)
void SetPSName(G4String &psName)
G4ThreeVector fCenterPosition
CLHEP::Hep3Vector G4ThreeVector
static G4VVisManager * GetConcreteInstance()
G4MultiFunctionalDetector * fMFD
Definition: G4Tubs.hh:85
void SetPSUnit(G4String &unit)
static G4double angle[DIM]
virtual void GetMapColor(G4double val, G4double color[4])=0
Transform3D inverse() const
Definition: Transform3D.cc:142
G4double fSize[3]
static constexpr double rad
Definition: G4SIunits.hh:149
HepPolyhedron & Transform(const G4Transform3D &t)
int G4int
Definition: G4Types.hh:78
void SetForceSolid(G4bool=true)
static constexpr double twopi
Definition: G4SIunits.hh:76
void SetMinMax(G4double minVal, G4double maxVal)
virtual void List() const
G4double fDrawUnitValue
G4GLOB_DLL std::ostream G4cout
static G4int GetReplicaLevel()
void SetVisibility(G4bool=true)
virtual void List() const
G4RotationMatrix * fRotationMatrix
virtual void EndDraw()=0
static constexpr double cm
Definition: G4SIunits.hh:119
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4CSGSolid.cc:134
virtual void BeginDraw(const G4Transform3D &objectTransformation=G4Transform3D())=0
void SetVisAttributes(const G4VisAttributes *)
Definition: G4Visible.cc:80
G4LogicalVolume * fMeshElementLogical
static constexpr double radian
Definition: G4SIunits.hh:142
G4ScoringCylinder(G4String wName)
HepGeom::Rotate3D G4Rotate3D
G4LogicalVolume * GetLogicalVolume() const
virtual void SetupGeometry(G4VPhysicalVolume *fWorldPhys)
tuple z
Definition: test.py:28
HepGeom::Translate3D G4Translate3D
std::map< G4int, T * > * GetMap() const
Definition: G4THitsMap.hh:99
G4bool IfFloatMinMax() const
#define G4endl
Definition: G4ios.hh:61
void SetForceAuxEdgeVisible(G4bool=true)
G4String fDivisionAxisNames[3]
double G4double
Definition: G4Types.hh:76
Definition: geomdefs.hh:54
tuple c
Definition: test.py:13
#define DBL_MAX
Definition: templates.hh:83
void SetVisAttributes(const G4VisAttributes *pVA)
void SetSensitiveDetector(G4VSensitiveDetector *pSDetector)
virtual void DrawColumn(RunScore *map, G4VScoreColorMap *colorMap, G4int idxProj, G4int idxColumn)
G4String fDrawPSName
G4GLOB_DLL std::ostream G4cerr
virtual void Draw(RunScore *map, G4VScoreColorMap *colorMap, G4int axflg=111)
virtual void DrawColorChart(G4int nPoint=5)
void GetRZPhi(G4int index, G4int q[3]) const