source: trunk/source/digits_hits/utils/src/G4ScoringCylinder.cc@ 1331

Last change on this file since 1331 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

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