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