Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4GMocrenFileSceneHandler.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// Created: Mar. 31, 2009 Akinori Kimura
30// Sep. 22, 2009 Akinori Kimura : modify and fix code to support
31// PrimitiveScorers and clean up
32//
33// GMocrenFile scene handler
34
35
36//----- header files
37#include <fstream>
38#include <cstdlib>
39#include <cstring>
40#include <sstream>
41#include <sstream>
42#include <iomanip>
43
44#include "globals.hh"
45#include "G4VisManager.hh"
46
47#include "G4GMocrenFile.hh"
50#include "G4Point3D.hh"
51#include "G4VisAttributes.hh"
52#include "G4Scene.hh"
53#include "G4Transform3D.hh"
54#include "G4Polyhedron.hh"
55#include "G4Box.hh"
56#include "G4Cons.hh"
57#include "G4Polyline.hh"
58#include "G4Trd.hh"
59#include "G4Tubs.hh"
60#include "G4Trap.hh"
61#include "G4Torus.hh"
62#include "G4Sphere.hh"
63#include "G4Para.hh"
64#include "G4Text.hh"
65#include "G4Circle.hh"
66#include "G4Square.hh"
67#include "G4VPhysicalVolume.hh"
69#include "G4LogicalVolume.hh"
70#include "G4Material.hh"
71
74#include "G4VisTrajContext.hh"
76#include "G4VTrajectoryModel.hh"
78#include "G4HitsModel.hh"
79#include "G4GMocrenMessenger.hh"
80//#include "G4PSHitsModel.hh"
81#include "G4GMocrenIO.hh"
83#include "G4GMocrenTouchable.hh"
86
87#include "G4ScoringManager.hh"
88#include "G4ScoringBox.hh"
89
91#include "G4SystemOfUnits.hh"
92
93//----- constants
94const char GDD_FILE_HEADER [] = "g4_";
95const char DEFAULT_GDD_FILE_NAME[] = "g4_00.gdd";
96
99
100//-- for a debugging
101const G4bool GFDEBUG = false;
102const G4bool GFDEBUG_TRK = false;//true;
103const G4bool GFDEBUG_HIT = false;//true;
104const G4bool GFDEBUG_DIGI = false;//true;
105const G4int GFDEBUG_DET = 0; // 0: false
106
107//////////////////////
108// static variables //
109//////////////////////
110
111//----- static variables
112G4int G4GMocrenFileSceneHandler::kSceneIdCount = 0;
113
114///////////////////////////
115// Driver-dependent part //
116///////////////////////////
117
118
119//----- G4GMocrenFileSceneHandler, constructor
121 G4GMocrenMessenger & messenger,
122 const G4String& name)
123 : G4VSceneHandler(system, kSceneIdCount++, name),
124 kSystem(system),
125 kMessenger(messenger),
126 kgMocrenIO(new G4GMocrenIO()),
127 kbSetModalityVoxelSize(false),
128 kbModelingTrajectory(false),
129// kGddDest(0),
130 kFlagInModeling(false),
131 kFlagSaving_g4_gdd(false),
132 kFlagParameterization(0),
133 kFlagProcessedInteractiveScorer(false) {
134
135 // g4.gdd filename and its directory
136 if(std::getenv("G4GMocrenFile_DEST_DIR") == NULL) {
137 kGddDestDir[0] = '\0';
138 //std::strcpy(kGddDestDir , ""); // output dir
139 //std::strcpy(kGddFileName, DEFAULT_GDD_FILE_NAME); // filename
140 std::strncpy(kGddFileName, DEFAULT_GDD_FILE_NAME,
141 std::strlen(DEFAULT_GDD_FILE_NAME)+1); // filename
142 } else {
143 const char * env = std::getenv("G4GMocrenFile_DEST_DIR");
144 G4int len = (G4int)std::strlen(env);
145 if(len > 256) {
146 G4Exception("G4GMocrenFileSceneHandler::G4GMocrenFileSceneHandler(*)",
147 "gMocren1000", FatalException,
148 "Invalid length of string set in G4GMocrenFile_DEST_DIR");
149 }
150 std::strncpy(kGddDestDir, env, len+1); // output dir
151 std::strncpy(kGddFileName, DEFAULT_GDD_FILE_NAME,
152 std::strlen(DEFAULT_GDD_FILE_NAME)+1); // filename
153 }
154
155 // maximum number of g4.gdd files in the dest directory
156 kMaxFileNum = FR_MAX_FILE_NUM ; // initialization
157 if ( std::getenv( "G4GMocrenFile_MAX_FILE_NUM" ) != NULL ) {
158 char * pcFileNum = std::getenv("G4GMocrenFile_MAX_FILE_NUM");
159 char c10FileNum[10];
160 std::strncpy(c10FileNum, pcFileNum, 9);
161 c10FileNum[9] = '\0';
162 kMaxFileNum = std::atoi(c10FileNum);
163
164 } else {
165 kMaxFileNum = FR_MAX_FILE_NUM ;
166 }
167 if( kMaxFileNum < 1 ) { kMaxFileNum = 1 ; }
168
169 InitializeParameters();
170
171}
172
173
174//----- G4GMocrenFileSceneHandler, destructor
176{
178 G4cout << "***** ~G4GMocrenFileSceneHandler" << G4endl;
179
180 if(kGddDest) {
181 //----- End of modeling
182 // close g4.gdd
184 }
185 if(kgMocrenIO != NULL) delete kgMocrenIO;
186
187}
188
189//----- initialize all parameters
190void G4GMocrenFileSceneHandler::InitializeParameters() {
191
192 kbSetModalityVoxelSize = false;
193
194 for(G4int i = 0; i < 3; i++) {
195 kModalitySize[i] = 0;
196 kNestedVolumeDimension[i] = 0;
197 kNestedVolumeDirAxis[i] = -1;
198 }
199
200 // delete kgMocrenIO;
201
202}
203
204//-----
206{
207 // g4_00.gdd, g4_01.gdd, ..., g4_MAX_FILE_INDEX.gdd
208 const G4int MAX_FILE_INDEX = kMaxFileNum - 1 ;
209
210 // dest directory (null if no environmental variables is set)
211 std::strncpy(kGddFileName, kGddDestDir, sizeof(kGddFileName)-1);
212 kGddFileName[sizeof(kGddFileName)-1] = '\0';
213
214 // create full path name (default)
215 std::strncat ( kGddFileName, DEFAULT_GDD_FILE_NAME,
216 sizeof(kGddFileName) - std::strlen(kGddFileName) - 1);
217
218 // Automatic updation of file names
219 static G4int currentNumber = 0;
220 for( G4int i = currentNumber ; i < kMaxFileNum ; i++) {
221
222 // Message in the final execution
223 if( i == MAX_FILE_INDEX )
224 {
226 G4cout << "===========================================" << G4endl;
227 G4cout << "WARNING MESSAGE from GMocrenFile driver: " << G4endl;
228 G4cout << " This file name is the final one in the " << G4endl;
229 G4cout << " automatic updation of the output file name." << G4endl;
230 G4cout << " You may overwrite existing files, i.e. " << G4endl;
231 G4cout << " g4_XX.gdd." << G4endl;
232 G4cout << "===========================================" << G4endl;
233 }
234 }
235
236 // re-determine file name as G4GMocrenFile_DEST_DIR/g4_XX.gdd
237 std::ostringstream filename;
238 filename
239 << kGddDestDir << GDD_FILE_HEADER
240 << std::setw(2) << std::setfill('0') << i << ".wrl";
241 strncpy(kGddFileName,filename.str().c_str(),sizeof(kGddFileName)-1);
242 kGddFileName[sizeof(kGddFileName)-1] = '\0';
243
244 // check validity of the file name
245 std::ifstream fin(kGddFileName);
246 if(GFDEBUG)
247 G4cout << "FILEOPEN: " << i << " : " << kGddFileName << fin.fail()
248 << G4endl;
249 if(!fin) {
250 // new file
251 fin.close();
252 currentNumber = i+1;
253 break;
254 } else {
255 // already exists (try next)
256 fin.close();
257 }
258
259 } // for
260
261 G4cout << "======================================================================" << G4endl;
262 G4cout << "Output file: " << kGddFileName << G4endl;
263 G4cout << "Destination directory (current dir if NULL): " << kGddDestDir << G4endl;
264 G4cout << "Maximum number of files in the destination directory: " << kMaxFileNum << G4endl;
265 G4cout << "Note:" << G4endl;
266 G4cout << " * The maximum number is customizable as: " << G4endl;
267 G4cout << " % setenv G4GMocrenFile_MAX_FILE_NUM number " << G4endl;
268 G4cout << " * The destination directory is customizable as:" << G4endl;
269 G4cout << " % setenv G4GMocrenFile_DEST_DIR dir_name/ " << G4endl;
270 G4cout << " ** Do not forget \"/\" at the end of the dir_name, e.g. \"./tmp/\"." << G4endl;
271 //G4cout << " dir_name, e.g. \"./tmp/\"." << G4endl;
272 G4cout << G4endl;
273 G4cout << "Maximum number of trajectories is set to " << MAX_NUM_TRAJECTORIES << "."<< G4endl;
274 G4cout << "======================================================================" << G4endl;
275
276} // G4GMocrenFileSceneHandler::SetGddFileName()
277
278
279//-----
281{
283 G4cout << "***** BeginSavingGdd (called)" << G4endl;
284
285 if( !IsSavingGdd() ) {
286
288 G4cout << "***** (started) " ;
289 G4cout << "(open g4.gdd, ##)" << G4endl;
290 }
291
292 SetGddFileName() ; // result set to kGddFileName
293 kFlagSaving_g4_gdd = true;
294
295
297 short minmax[2];
298 minmax[0] = ctdens.GetMinCT();
299 minmax[1] = ctdens.GetMaxCT();
300 kgMocrenIO->setModalityImageMinMax(minmax);
301 std::vector<G4float> map;
302 G4float dens;
303 for(G4int i = minmax[0]; i <= minmax[1]; i++) {
304 dens = ctdens.GetDensity(i);
305 map.push_back(dens);
306 }
307 kgMocrenIO->setModalityImageDensityMap(map);
308
309 /*
310 G4String fname = "modality-map.dat";
311 std::ifstream ifile(fname);
312 if(ifile) {
313 short minmax[2];
314 ifile >> minmax[0] >> minmax[1];
315 kgMocrenIO->setModalityImageMinMax(minmax);
316 std::vector<G4float> map;
317 G4float dens;
318 for(G4int i = minmax[0]; i <= minmax[1]; i++) {
319 ifile >> dens;
320 map.push_back(dens);
321 }
322 kgMocrenIO->setModalityImageDensityMap(map);
323
324 } else {
325 G4cout << "cann't open the file : " << fname << G4endl;
326 }
327 */
328
329 // mesh size
330 //kMessenger.getNoVoxels(kModalitySize[0], kModalitySize[1], kModalitySize[2]);
331 //kgMocrenIO->setModalityImageSize(kModalitySize);
332
333 // initializations
334 //kgMocrenIO->clearModalityImage();
335 kgMocrenIO->clearDoseDistAll();
336 kgMocrenIO->clearROIAll();
337 kgMocrenIO->clearTracks();
338 kgMocrenIO->clearDetector();
339 std::vector<Detector>::iterator itr = kDetectors.begin();
340 for(; itr != kDetectors.end(); itr++) {
341 itr->clear();
342 }
343 kDetectors.clear();
344
345 kNestedHitsList.clear();
346 kNestedVolumeNames.clear();
347
348 }
349}
350
352{
354 G4cout << "***** EndSavingGdd (called)" << G4endl;
355
356 if(IsSavingGdd()) {
358 G4cout << "***** (started) (close "
359 << kGddFileName << ")" << G4endl;
360
361 if(kGddDest) kGddDest.close();
362 kFlagSaving_g4_gdd = false;
363
364 std::map<Index3D, G4float>::iterator itr = kNestedModality.begin();
365 G4int xmax=0, ymax=0, zmax=0;
366 for(; itr != kNestedModality.end(); itr++) {
367 if(itr->first.x > xmax) xmax = itr->first.x;
368 if(itr->first.y > ymax) ymax = itr->first.y;
369 if(itr->first.z > zmax) zmax = itr->first.z;
370 }
371 // mesh size
372 kModalitySize[0] = xmax+1;
373 kModalitySize[1] = ymax+1;
374 kModalitySize[2] = zmax+1;
375 kgMocrenIO->setModalityImageSize(kModalitySize);
376 if(GFDEBUG) G4cout << "gMocren-file driver : modality size : "
377 << kModalitySize[0] << " x "
378 << kModalitySize[1] << " x "
379 << kModalitySize[2] << G4endl;
380
381 G4int nxy = kModalitySize[0]*kModalitySize[1];
382 //std::map<G4int, G4float>::iterator itr;
383 for(G4int z = 0; z < kModalitySize[2]; z++) {
384 short * modality = new short[nxy];
385 for(G4int y = 0; y < kModalitySize[1]; y++) {
386 for(G4int x = 0; x < kModalitySize[0]; x++) {
387 //for(G4int x = kModalitySize[0]-1; x >= 0 ; x--) {
388 //G4int ixy = x + (kModalitySize[1]-y-1)*kModalitySize[0];
389
390 G4int ixy = x + y*kModalitySize[0];
391 Index3D idx(x,y,z);
392 itr = kNestedModality.find(idx);
393 if(itr != kNestedModality.end()) {
394
395 modality[ixy] = kgMocrenIO->convertDensityToHU(itr->second);
396 } else {
397 modality[ixy] = -1024;
398 }
399
400 }
401 }
402 kgMocrenIO->setModalityImage(modality);
403 }
404
405 //-- dose
406 size_t nhits = kNestedHitsList.size();
407 if(GFDEBUG) G4cout << "gMocren-file driver : # hits = " << nhits << G4endl;
408
409 std::map<Index3D, G4double>::iterator hitsItr;
410 std::map<G4String, std::map<Index3D, G4double> >::iterator hitsListItr = kNestedHitsList.begin();
411
412 for(G4int n = 0; hitsListItr != kNestedHitsList.end(); hitsListItr++, n++) {
413
414 kgMocrenIO->newDoseDist();
415 kgMocrenIO->setDoseDistName(hitsListItr->first, n);
416 kgMocrenIO->setDoseDistSize(kModalitySize, n);
417
418 G4double minmax[2] = {DBL_MAX, -DBL_MAX};
419 for(G4int z = 0 ; z < kModalitySize[2]; z++) {
420 G4double * values = new G4double[nxy];
421 for(G4int y = 0; y < kModalitySize[1]; y++) {
422 for(G4int x = 0; x < kModalitySize[0]; x++) {
423
424 G4int ixy = x + y*kModalitySize[0];
425 Index3D idx(x,y,z);
426 hitsItr = hitsListItr->second.find(idx);
427 if(hitsItr != hitsListItr->second.end()) {
428
429 values[ixy] = hitsItr->second;
430 } else {
431 values[ixy] = 0.;
432 }
433 if(values[ixy] < minmax[0]) minmax[0] = values[ixy];
434 if(values[ixy] > minmax[1]) minmax[1] = values[ixy];
435 }
436 }
437 kgMocrenIO->setDoseDist(values, n);
438 }
439 kgMocrenIO->setDoseDistMinMax(minmax, n);
440 G4double lower = 0.;
441 if(minmax[0] < 0) lower = minmax[0];
442 G4double scale = (minmax[1]-lower)/25000.;
443 kgMocrenIO->setDoseDistScale(scale, n);
444 G4String sunit("unit?"); //temporarily
445 kgMocrenIO->setDoseDistUnit(sunit, n);
446 }
447
448
449 //-- draw axes
450 if(false) {//true,false
451 G4ThreeVector trans;
453 trans = kVolumeTrans3D.getTranslation();
454 rot = kVolumeTrans3D.getRotation().inverse();
455 // x
456 std::vector<G4float *> tracks;
457 unsigned char colors[3];
458 G4float * trk = new G4float[6];
459 tracks.push_back(trk);
460
461 G4ThreeVector orig(0.,0.,0), xa(2000.,0.,0.), ya(0.,2000.,0.), za(0.,0.,2000.);
462 orig -= trans;
463 orig.transform(rot);
464 xa -= trans;
465 xa.transform(rot);
466 ya -= trans;
467 ya.transform(rot);
468 za -= trans;
469 za.transform(rot);
470 for(G4int i = 0; i < 3; i++) trk[i] = orig[i];
471 for(G4int i = 0; i < 3; i++) trk[i+3] = xa[i];
472 colors[0] = 255; colors[1] = 0; colors[2] = 0;
473 kgMocrenIO->addTrack(tracks, colors);
474 // y
475 for(G4int i = 0; i < 3; i++) trk[i+3] = ya[i];
476 colors[0] = 0; colors[1] = 255; colors[2] = 0;
477 kgMocrenIO->addTrack(tracks, colors);
478 // z
479 for(G4int i = 0; i < 3; i++) trk[i+3] = za[i];
480 colors[0] = 0; colors[1] = 0; colors[2] = 255;
481 kgMocrenIO->addTrack(tracks, colors);
482 }
483
484 //-- detector
485 ExtractDetector();
486
487
488 if(GFDEBUG_DET) G4cout << ">>>>>>>>>>>>>>>>>>>>>> (";
489 std::vector<G4float> transformObjects;
490 for(G4int i = 0; i < 3; i++) {
491 // need to check!!
492 transformObjects.push_back((kVolumeSize[i]/2. - kVoxelDimension[i]/2.));
493 if(GFDEBUG_DET) G4cout << transformObjects[i] << ", ";
494 }
495 if(GFDEBUG_DET) G4cout << ")" << G4endl;
496
497
498 kgMocrenIO->translateTracks(transformObjects);
499 kgMocrenIO->translateDetector(transformObjects);
500
501 // store
502 kgMocrenIO->storeData(kGddFileName);
503 }
504
505}
506
507
508//-----
510{
512
513 if( !GFIsInModeling() ) {
514
516 G4cout << "***** G4GMocrenFileSceneHandler::GFBeginModeling (called & started)" << G4endl;
517
518 //----- Send saving command and heading comment
520
521 kFlagInModeling = true ;
522
523 // These models are entrusted to user commands /vis/scene/add/psHits or hits
524 //GetScene()->AddEndOfEventModel(new G4PSHitsModel());
525 //GetScene()->AddEndOfRunModel(new G4PSHitsModel());
526 //scene->AddEndOfEventModel(new G4HitsModel());
527 if(GFDEBUG_HIT) {
528 G4Scene * scene = GetScene();
529 std::vector<G4Scene::Model> vmodel = scene->GetEndOfEventModelList();
530 std::vector<G4Scene::Model>::iterator itr = vmodel.begin();
531 for(; itr != vmodel.end(); itr++) {
532 G4cout << " IIIIII model name: " << itr->fpModel->GetGlobalTag() << G4endl;
533 }
534 }
535 }
536}
537
538
539//========== AddPrimitive() functions ==========//
540
541//----- Add polyline
543{
545 G4cout << "***** AddPrimitive" << G4endl;
546
547 if (fProcessing2D) {
548 static G4bool warned = false;
549 if (!warned) {
550 warned = true;
552 ("G4GMocrenFileSceneHandler::AddPrimitive (const G4Polyline&)",
553 "gMocren1001", JustWarning,
554 "2D polylines not implemented. Ignored.");
555 }
556 return;
557 }
558
559 //----- Initialize if necessary
561
562 static G4int numTrajectories = 0;
563 if(numTrajectories >= MAX_NUM_TRAJECTORIES) return;
564
565 // draw trajectories
566 if(kbModelingTrajectory) {
567
568 G4TrajectoriesModel * pTrModel = dynamic_cast<G4TrajectoriesModel*>(fpModel);
569 if (!pTrModel) {
571 ("G4VSceneHandler::AddCompound(const G4Polyline&)",
572 "gMocren0002", FatalException, "Not a G4TrajectoriesModel.");
573 }
574
575 G4ThreeVector trans;
577 trans = kVolumeTrans3D.getTranslation();
578 rot = kVolumeTrans3D.getRotation().inverse();
579
580 if(GFDEBUG_TRK) G4cout << " trajectory points : " << G4endl;
581 std::vector<G4float *> trajectory;
582 if(polyline.size() < 2) return;
583 G4Polyline::const_iterator preitr = polyline.begin();
584 G4Polyline::const_iterator postitr = preitr; postitr++;
585 for(; postitr != polyline.end(); preitr++, postitr++) {
586 G4ThreeVector prePts(preitr->x(), preitr->y(), preitr->z());
587 prePts -= trans;
588 prePts.transform(rot);
589 G4ThreeVector postPts(postitr->x(), postitr->y(), postitr->z());
590 postPts -= trans;
591 postPts.transform(rot);
592 G4float * stepPts = new G4float[6];
593 stepPts[0] = prePts.x();
594 stepPts[1] = prePts.y();
595 stepPts[2] = prePts.z();
596 stepPts[3] = postPts.x();
597 stepPts[4] = postPts.y();
598 stepPts[5] = postPts.z();
599 trajectory.push_back(stepPts);
600
601 if(GFDEBUG_TRK) {
602 G4cout << " ("
603 << stepPts[0] << ", "
604 << stepPts[1] << ", "
605 << stepPts[2] << ") - ("
606 << stepPts[3] << ", "
607 << stepPts[4] << ", "
608 << stepPts[5] << ")" << G4endl;
609 }
610 }
611
612 const G4VisAttributes * att = polyline.GetVisAttributes();
613 G4Color color = att->GetColor();
614 unsigned char trkcolor[3];
615 trkcolor[0] = (unsigned char)(color.GetRed()*255);
616 trkcolor[1] = (unsigned char)(color.GetGreen()*255);
617 trkcolor[2] = (unsigned char)(color.GetBlue()*255);
618 if(GFDEBUG_TRK) {
619 G4cout << " color : ["
620 << color.GetRed() << ", "
621 << color.GetGreen() << ", "
622 << color.GetBlue() << "]" << G4endl;
623 }
624
625 kgMocrenIO->addTrack(trajectory, trkcolor);
626
627 numTrajectories++;
628 }
629
630} // G4GMocrenFileSceneHandler::AddPrimitive (polyline)
631
632
633//----- Add text
635{
636 if (fProcessing2D) {
637 static G4bool warned = false;
638 if (!warned) {
639 warned = true;
641 ("G4GMocrenFileSceneHandler::AddPrimitive (const G4Text&)",
642 "gMocren1002", JustWarning,
643 "2D text not implemented. Ignored.");
644 }
645 return;
646 }
647
648 // to avoid a warning in the compile process
649 G4Text dummytext = text;
650
651 //-----
653 G4cout << "***** AddPrimitive( G4Text )" << G4endl;
654
655 //----- Initialize IF NECESSARY
657
658} // G4GMocrenFileSceneHandler::AddPrimitive ( text )
659
660
661//----- Add circle
663{
664 // to avoid a warning in the compile process
665 G4Circle dummycircle = mark_circle;
666
667 if (fProcessing2D) {
668 static G4bool warned = false;
669 if (!warned) {
670 warned = true;
672 ("G4GMocrenFileSceneHandler::AddPrimitive (const G4Circle&)",
673 "gMocren1003", JustWarning,
674 "2D circles not implemented. Ignored.");
675 }
676 return;
677 }
678
679 //-----
681 G4cout << "***** AddPrimitive( G4Circle )" << G4endl;
682
683 //----- Initialize IF NECESSARY
685
686
687} // G4GMocrenFileSceneHandler::AddPrimitive ( mark_circle )
688
689
690//----- Add square
692{
693 // to avoid a warning in the compile process
694 G4Square dummysquare = mark_square;
695
696 if (fProcessing2D) {
697 static G4bool warned = false;
698 if (!warned) {
699 warned = true;
701 ("G4GMocrenFileSceneHandler::AddPrimitive (const G4Square&)",
702 "gMocren1004", JustWarning,
703 "2D squares not implemented. Ignored.");
704 }
705 return;
706 }
707
708 //-----
710 G4cout << "***** AddPrimitive( G4Square )" << G4endl;
711
712 //----- Initialize if necessary
714
715} // G4GMocrenFileSceneHandler::AddPrimitive ( mark_square )
716
717
718//----- Add polyhedron
720{
721 //-----
723 G4cout << "***** AddPrimitive( G4Polyhedron )" << G4endl;
724
725
726 if (polyhedron.GetNoFacets() == 0) return;
727
728 if (fProcessing2D) {
729 static G4bool warned = false;
730 if (!warned) {
731 warned = true;
733 ("G4GMocrenFileSceneHandler::AddPrimitive (const G4Polyhedron&)",
734 "gMocren1005", JustWarning,
735 "2D polyhedra not implemented. Ignored.");
736 }
737 return;
738 }
739
740 //----- Initialize if necessary
742
743 //---------- (3) Facet block
744 for (G4int f = polyhedron.GetNoFacets(); f; f--){
745 G4bool notLastEdge = true;
746 G4int index = -1; // initialization
747 G4int edgeFlag = 1;
748 //G4int preedgeFlag = 1;
749 //G4int work[4], i = 0;
750 G4int i = 0;
751 do {
752 //preedgeFlag = edgeFlag;
753 notLastEdge = polyhedron.GetNextVertexIndex(index, edgeFlag);
754 //work[i++] = index;
755 i++;
756 }while (notLastEdge);
757 switch (i){
758 case 3:
759 //SendStrInt3(FR_FACET, work[0], work[1], work[2] );
760 break;
761 case 4:
762 //SendStrInt4(FR_FACET, work[0], work[1], work[2], work[3] );
763 break;
764 default:
766 G4cout <<
767 "ERROR G4GMocrenFileSceneHandler::AddPrimitive(G4Polyhedron)" << G4endl;
768 G4PhysicalVolumeModel* pPVModel =
769 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
770 if (pPVModel)
772 G4cout << "Volume " << pPVModel->GetCurrentPV()->GetName() <<
773 ", Solid " << pPVModel->GetCurrentLV()->GetSolid()->GetName() <<
774 " (" << pPVModel->GetCurrentLV()->GetSolid()->GetEntityType();
775
777 G4cout <<
778 "\nG4Polyhedron facet with " << i << " edges" << G4endl;
779 }
780 }
781
782} // G4GMocrenFileSceneHandler::AddPrimitive (polyhedron)
783
784
785//-----
787{
789
790 //-----
792 G4cout << "***** GFEndModeling (called)" << G4endl;
793
794 if( GFIsInModeling() ) {
795
797 G4cout << "***** GFEndModeling (started) " ;
798 G4cout << "(/EndModeling, /DrawAll, /CloseDevice)" << G4endl;
799 }
800
801 //----- End saving data to g4.gdd
802 EndSavingGdd() ;
803
804 //------ Reset flag
805 kFlagInModeling = false ;
806
807 }
808
809}
810
811
812//-----
814{
816 G4cout << "***** BeginPrimitives " << G4endl;
817
819
820 G4VSceneHandler::BeginPrimitives (objectTransformation);
821
822}
823
824
825//-----
827{
829 G4cout << "***** EndPrimitives " << G4endl;
830
832}
833
834
835//========== AddSolid() functions ==========//
836
837//----- Add box
839{
841 G4cout << "***** AddSolid ( box )" << G4endl;
842
843 if(GFDEBUG_DET > 0)
844 G4cout << "G4GMocrenFileSceneHandler::AddSolid(const G4Box&) : "
845 << box.GetName() << G4endl;
846
847 //----- skip drawing invisible primitive
848 if( !IsVisible() ) { return ; }
849
850 //----- Initialize if necessary
852
853
854 //--
855 if(GFDEBUG_DET > 1) {
856 G4cout << "-------" << G4endl;
857 G4cout << " " << box.GetName() << G4endl;
858 G4Polyhedron * poly = box.CreatePolyhedron();
860 //G4int nv = poly->GetNoVertices();
861 G4Point3D v1, v2;
862 G4int next;
863 //while(1) { // next flag isn't functional.
864 for(G4int i = 0; i < 12; i++) { // # of edges is 12.
865 poly->GetNextEdge(v1, v2, next);
866 if(next == 0) break;
867 G4cout << " (" << v1.x() << ", "
868 << v1.y() << ", "
869 << v1.z() << ") - ("
870 << v2.x() << ", "
871 << v2.y() << ", "
872 << v2.z() << ") [" << next << "]"
873 << G4endl;
874 }
875 delete poly;
876 }
877
878
879 // the volume name set by /vis/gMocren/setVolumeName
880 G4String volName = kMessenger.getVolumeName();
881
882
883 if(kFlagParameterization != 2) {
885 if(pScrMan) {
886 G4ScoringBox * pScBox = dynamic_cast<G4ScoringBox*>(pScrMan->FindMesh(volName));
887 G4bool bMesh = false;
888 if(pScBox != NULL) bMesh = true;
889 if(bMesh) kFlagParameterization = 2;
890 if(GFDEBUG_DET > 0) G4cout << " G4ScoringManager::FindMesh() : "
891 << volName << " - " << bMesh << G4endl;
892 }
893 }
894
895 const G4VModel* pv_model = GetModel();
896 if (!pv_model) { return ; }
897 G4PhysicalVolumeModel* pPVModel =
898 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
899 if (!pPVModel) { return ; }
900
901
902 //-- debug information
903 if(GFDEBUG_DET > 0) {
904 G4Material * mat = pPVModel->GetCurrentMaterial();
905 G4String name = mat->GetName();
906 G4double dens = mat->GetDensity()/(g/cm3);
907 G4int copyNo = pPVModel->GetCurrentPV()->GetCopyNo();
908 G4int depth = pPVModel->GetCurrentDepth();
909 G4cout << " copy no.: " << copyNo << G4endl;
910 G4cout << " depth : " << depth << G4endl;
911 G4cout << " density : " << dens << " [g/cm3]" << G4endl;
912 G4cout << " location: " << pPVModel->GetCurrentPV()->GetObjectTranslation() << G4endl;
913 G4cout << " Multiplicity : " << pPVModel->GetCurrentPV()->GetMultiplicity() << G4endl;
914 G4cout << " Is replicated? : " << pPVModel->GetCurrentPV()->IsReplicated() << G4endl;
915 G4cout << " Is parameterised? : " << pPVModel->GetCurrentPV()->IsParameterised() << G4endl;
916 G4cout << " top phys. vol. name : " << pPVModel->GetTopPhysicalVolume()->GetName() << G4endl;
917 }
918
919 //-- check the parameterised volume
920 if(box.GetName() == volName) {
921
922 kVolumeTrans3D = fObjectTransformation;
923 // coordination system correction for gMocren
924 G4ThreeVector raxis(1., 0., 0.), dummy(0.,0.,0.);
925 G4RotationMatrix rot(raxis, pi*rad);
926 G4Transform3D trot(rot, dummy);
927 if(GFDEBUG_DET) {
928 G4ThreeVector trans1 = kVolumeTrans3D.getTranslation();
929 G4RotationMatrix rot1 = kVolumeTrans3D.getRotation().inverse();
930 G4cout << "kVolumeTrans3D: " << trans1 << G4endl << rot1 << G4endl;
931 }
932 kVolumeTrans3D = kVolumeTrans3D*trot;
933 if(GFDEBUG_DET) G4cout << " Parameterised volume : " << box.GetName() << G4endl;
934
935
936
937 //
938 G4VPhysicalVolume * pv[3] = {0,0,0};
939 pv[0] = pPVModel->GetCurrentPV()->GetLogicalVolume()->GetDaughter(0);
940 if(!pv[0]) {
941 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
942 "gMocren0003", FatalException, "Unexpected volume.");
943 }
944 G4int dirAxis[3] = {-1,-1,-1};
945 G4int nDaughters[3] = {0,0,0};
946
947 EAxis axis; G4int nReplicas; G4double width; G4double offset; G4bool consuming;
948 pv[0]->GetReplicationData(axis, nReplicas, width, offset, consuming);
949 nDaughters[0] = nReplicas;
950 switch(axis) {
951 case kXAxis: dirAxis[0] = 0; break;
952 case kYAxis: dirAxis[0] = 1; break;
953 case kZAxis: dirAxis[0] = 2; break;
954 default:
955 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
956 "gMocren0004", FatalException, "Error.");
957 }
958 kNestedVolumeNames.push_back(pv[0]->GetName());
959 if(GFDEBUG_DET)
960 G4cout << " daughter name : " << pv[0]->GetName()
961 << " # : " << nDaughters[0] << G4endl;
962
963 //
964 if(GFDEBUG_DET) {
965 if(pv[0]->GetLogicalVolume()->GetNoDaughters()) {
966 G4cout << "# of daughters : "
967 << pv[0]->GetLogicalVolume()->GetNoDaughters() << G4endl;
968 } else {
969 //G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
970 // "gMocren0005", FatalException, "Error.");
971 }
972 }
973
974 // check whether nested or regular parameterization
975 if(GFDEBUG_DET) G4cout << "# of daughters : "
976 << pv[0]->GetLogicalVolume()->GetNoDaughters() << G4endl;
977 if(pv[0]->GetLogicalVolume()->GetNoDaughters() == 0) {
978 kFlagParameterization = 1;
979 //G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
980 // "gMocren0006", FatalException, "Error.");
981 }
982
983 if(kFlagParameterization == 0) {
984
985 pv[1] = pv[0]->GetLogicalVolume()->GetDaughter(0);
986 if(pv[1]) {
987 pv[1]->GetReplicationData(axis, nReplicas, width, offset, consuming);
988 nDaughters[1] = nReplicas;
989 switch(axis) {
990 case kXAxis: dirAxis[1] = 0; break;
991 case kYAxis: dirAxis[1] = 1; break;
992 case kZAxis: dirAxis[1] = 2; break;
993 default:
994 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
995 "gMocren0007", FatalException, "Error.");
996 }
997 kNestedVolumeNames.push_back(pv[1]->GetName());
998 if(GFDEBUG_DET)
999 G4cout << " sub-daughter name : " << pv[1]->GetName()
1000 << " # : " << nDaughters[1]<< G4endl;
1001
1002 //
1003 pv[2] = pv[1]->GetLogicalVolume()->GetDaughter(0);
1004 if(pv[2]) {
1005 nDaughters[2] = pv[2]->GetMultiplicity();
1006 kNestedVolumeNames.push_back(pv[2]->GetName());
1007 if(GFDEBUG_DET)
1008 G4cout << " sub-sub-daughter name : " << pv[2]->GetName()
1009 << " # : " << nDaughters[2] << G4endl;
1010
1011 if(nDaughters[2] > 1) {
1012 G4VNestedParameterisation * nestPara
1013 = dynamic_cast<G4VNestedParameterisation*>(pv[2]->GetParameterisation());
1014 if(nestPara == NULL)
1015 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
1016 "gMocren0008", FatalException, "Non-nested parameterisation");
1017
1018 nestPara->ComputeTransformation(0, pv[2]);
1019 G4ThreeVector trans0 = pv[2]->GetObjectTranslation();
1020 nestPara->ComputeTransformation(1, pv[2]);
1021 G4ThreeVector trans1 = pv[2]->GetObjectTranslation();
1022 G4ThreeVector diff(trans0 - trans1);
1023 if(GFDEBUG_DET)
1024 G4cout << trans0 << " - " << trans1 << " - " << diff << G4endl;
1025
1026 if(diff.x() != 0.) dirAxis[2] = 0;
1027 else if(diff.y() != 0.) dirAxis[2] = 1;
1028 else if(diff.z() != 0.) dirAxis[2] = 2;
1029 else
1030 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
1031 "gMocren0009", FatalException, "Unexpected nested parameterisation");
1032 }
1033 }
1034 }
1035
1036 for(G4int i = 0; i < 3; i++) {
1037 kNestedVolumeDimension[i] = nDaughters[i];
1038 //kNestedVolumeDimension[i] = nDaughters[dirAxis[i]];
1039 kNestedVolumeDirAxis[i] = dirAxis[i];
1040 }
1041 //G4cout << "@@@@@@@@@ "
1042 // << dirAxis[0] << ", " << dirAxis[1] << ", " << dirAxis[2] << G4endl;
1043
1044 // get densities
1045 G4VNestedParameterisation * nestPara
1046 = dynamic_cast<G4VNestedParameterisation*>(pv[2]->GetParameterisation());
1047 if(nestPara != NULL) {
1048 G4double prexyz[3] = {0.,0.,0.}, xyz[3] = {0.,0.,0.};
1049 for(G4int n0 = 0; n0 < nDaughters[0]; n0++) {
1050 for(G4int n1 = 0; n1 < nDaughters[1]; n1++) {
1051 for(G4int n2 = 0; n2 < nDaughters[2]; n2++) {
1052
1053 G4GMocrenTouchable * touch = new G4GMocrenTouchable(n1, n0);
1054 if(GFDEBUG_DET)
1055 G4cout << " retrieve volume : copy # : " << n0
1056 << ", " << n1 << ", " << n2 << G4endl;
1057 G4Material * mat = nestPara->ComputeMaterial(pv[2], n2, touch);
1058 delete touch;
1059 G4double dens = mat->GetDensity()/(g/cm3);
1060
1061 if(GFDEBUG_DET)
1062 G4cout << " density :" << dens << " [g/cm3]" << G4endl;
1063
1064 G4Box tbox(box);
1065 nestPara->ComputeDimensions(tbox, n2, pv[2]);
1066 xyz[0] = tbox.GetXHalfLength()/mm;
1067 xyz[1] = tbox.GetYHalfLength()/mm;
1068 xyz[2] = tbox.GetZHalfLength()/mm;
1069 if(n0 != 0 || n1 != 0 || n2 != 0) {
1070 for(G4int i = 0; i < 3; i++) {
1071 if(xyz[i] != prexyz[i])
1072 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
1073 "gMocren0010", FatalException, "Unsupported parameterisation");
1074 }
1075 }
1076 if(GFDEBUG_DET)
1077 G4cout << " size : " << tbox.GetXHalfLength()/mm << " x "
1078 << tbox.GetYHalfLength()/mm << " x "
1079 << tbox.GetZHalfLength()/mm << " [mm3]" << G4endl;
1080
1081 G4int idx[3];
1082 idx[dirAxis[0]] = n0;
1083 idx[dirAxis[1]] = n1;
1084 idx[dirAxis[2]] = n2;
1085 Index3D i3d(idx[0],idx[1],idx[2]);
1086 kNestedModality[i3d] = dens;
1087 if(GFDEBUG_DET)
1088 G4cout << " index: " << idx[0] << ", " << idx[1] << ", " << idx[2]
1089 << " density: " << dens << G4endl;
1090
1091 for(G4int i = 0; i < 3; i++) prexyz[i] = xyz[i];
1092 }
1093 }
1094 }
1095
1096 kVolumeSize.set(box.GetXHalfLength()*2/mm,
1097 box.GetYHalfLength()*2/mm,
1098 box.GetZHalfLength()*2/mm);
1099 // mesh size
1100 if(!kbSetModalityVoxelSize) {
1101 G4float spacing[3] = {static_cast<G4float>(2*xyz[0]),
1102 static_cast<G4float>(2*xyz[1]),
1103 static_cast<G4float>(2*xyz[2])};
1104 kgMocrenIO->setVoxelSpacing(spacing);
1105 kVoxelDimension.set(spacing[0], spacing[1], spacing[2]);
1106 kbSetModalityVoxelSize = true;
1107 }
1108
1109 } else {
1110 if(GFDEBUG_DET)
1111 G4cout << pv[2]->GetName() << G4endl;
1112 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
1113 "gMocren0011", FatalException, "Non-nested parameterisation");
1114 }
1115
1116
1117
1118 //-- debug
1119 if(GFDEBUG_DET > 1) {
1120 if(pPVModel->GetCurrentPV()->IsParameterised()) {
1122 G4cout << " Is nested parameterisation? : " << para->IsNested() << G4endl;
1123
1124
1125 G4int npvp = pPVModel->GetDrawnPVPath().size();
1126 G4cout << " physical volume node id : "
1127 << "size: " << npvp << ", PV name: ";
1128 for(G4int i = 0; i < npvp; i++) {
1129 G4cout << pPVModel->GetDrawnPVPath()[i].GetPhysicalVolume()->GetName()
1130 << " [param:"
1131 << pPVModel->GetDrawnPVPath()[i].GetPhysicalVolume()->IsParameterised()
1132 << ",rep:"
1133 << pPVModel->GetDrawnPVPath()[i].GetPhysicalVolume()->IsReplicated();
1134 if(pPVModel->GetDrawnPVPath()[i].GetPhysicalVolume()->GetParameterisation()) {
1135 G4cout << ",nest:"
1136 << pPVModel->GetDrawnPVPath()[i].GetPhysicalVolume()->GetParameterisation()->IsNested();
1137 }
1138 G4cout << ",copyno:"
1139 << pPVModel->GetDrawnPVPath()[i].GetPhysicalVolume()->GetCopyNo();
1140 G4cout << "] - ";
1141 }
1142 G4cout << G4endl;
1143
1144
1145 pPVModel->GetCurrentPV()->GetReplicationData(axis, nReplicas, width, offset, consuming);
1146 G4cout << " # replicas : " << nReplicas << G4endl;
1147 G4double pareDims[3] = {0.,0.,0.};
1148 G4Box * pbox = dynamic_cast<G4Box *>(pPVModel->GetDrawnPVPath()[npvp-2].GetPhysicalVolume()->GetLogicalVolume()->GetSolid());
1149 if(pbox) {
1150 pareDims[0] = 2.*pbox->GetXHalfLength()*mm;
1151 pareDims[1] = 2.*pbox->GetYHalfLength()*mm;
1152 pareDims[2] = 2.*pbox->GetZHalfLength()*mm;
1153 G4cout << " mother size ["
1154 << pPVModel->GetDrawnPVPath()[npvp-2].GetPhysicalVolume()->GetName()
1155 << "] : "
1156 << pareDims[0] << " x "
1157 << pareDims[1] << " x "
1158 << pareDims[2] << " [mm3]"
1159 << G4endl;
1160 }
1161 G4double paraDims[3];
1162 G4Box * boxP = dynamic_cast<G4Box *>(pPVModel->GetDrawnPVPath()[npvp-1].GetPhysicalVolume()->GetLogicalVolume()->GetSolid());
1163 if(boxP) {
1164 paraDims[0] = 2.*boxP->GetXHalfLength()*mm;
1165 paraDims[1] = 2.*boxP->GetYHalfLength()*mm;
1166 paraDims[2] = 2.*boxP->GetZHalfLength()*mm;
1167 G4cout << " parameterised volume? ["
1168 << pPVModel->GetDrawnPVPath()[npvp-1].GetPhysicalVolume()->GetName()
1169 << "] : "
1170 << paraDims[0] << " x "
1171 << paraDims[1] << " x "
1172 << paraDims[2] << " [mm3] : "
1173 << G4int(pareDims[0]/paraDims[0]) << " x "
1174 << G4int(pareDims[1]/paraDims[1]) << " x "
1175 << G4int(pareDims[2]/paraDims[2]) << G4endl;
1176 } else {
1177 G4cout << pPVModel->GetDrawnPVPath()[npvp-2].GetPhysicalVolume()->GetName()
1178 << " isn't a G4Box." << G4endl;
1179 }
1180 }
1181 }
1182
1183
1184 } else if(kFlagParameterization == 1) { // G4PhantomParameterisation based geom. construnction
1185
1186 // get the dimension of the parameterized patient geometry
1187 G4PhantomParameterisation * phantomPara
1188 = dynamic_cast<G4PhantomParameterisation*>(pv[0]->GetParameterisation());
1189 if(phantomPara == NULL) {
1190 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
1191 "gMocren0012", FatalException, "no G4PhantomParameterisation");
1192 } else {
1193 ;
1194 }
1195
1196 kNestedVolumeDimension[0] = (G4int)phantomPara->GetNoVoxelsX();
1197 kNestedVolumeDimension[1] = (G4int)phantomPara->GetNoVoxelsY();
1198 kNestedVolumeDimension[2] = (G4int)phantomPara->GetNoVoxelsZ();
1199 kNestedVolumeDirAxis[0] = 0;
1200 kNestedVolumeDirAxis[1] = 1;
1201 kNestedVolumeDirAxis[2] = 2;
1202
1203 // get densities of the parameterized patient geometry
1204 G4int nX = kNestedVolumeDimension[0];
1205 G4int nXY = kNestedVolumeDimension[0]*kNestedVolumeDimension[1];
1206
1207 for(G4int n0 = 0; n0 < kNestedVolumeDimension[0]; n0++) {
1208 for(G4int n1 = 0; n1 < kNestedVolumeDimension[1]; n1++) {
1209 for(G4int n2 = 0; n2 < kNestedVolumeDimension[2]; n2++) {
1210
1211 G4int repNo = n0 + n1*nX + n2*nXY;
1212 G4Material * mat = phantomPara->ComputeMaterial(repNo, pv[0]);
1213 G4double dens = mat->GetDensity()/(g/cm3);
1214
1215
1216 G4int idx[3];
1217 idx[kNestedVolumeDirAxis[0]] = n0;
1218 idx[kNestedVolumeDirAxis[1]] = n1;
1219 idx[kNestedVolumeDirAxis[2]] = n2;
1220 Index3D i3d(idx[0],idx[1],idx[2]);
1221 kNestedModality[i3d] = dens;
1222
1223 if(GFDEBUG_DET)
1224 G4cout << " index: " << idx[0] << ", " << idx[1] << ", " << idx[2]
1225 << " density: " << dens << G4endl;
1226
1227 }
1228 }
1229 }
1230
1231 kVolumeSize.set(box.GetXHalfLength()*2/mm,
1232 box.GetYHalfLength()*2/mm,
1233 box.GetZHalfLength()*2/mm);
1234
1235 // mesh size
1236 if(!kbSetModalityVoxelSize) {
1237 G4float spacing[3] = {static_cast<G4float>(2*phantomPara->GetVoxelHalfX()),
1238 static_cast<G4float>(2*phantomPara->GetVoxelHalfY()),
1239 static_cast<G4float>(2*phantomPara->GetVoxelHalfZ())};
1240 kgMocrenIO->setVoxelSpacing(spacing);
1241 kVoxelDimension.set(spacing[0], spacing[1], spacing[2]);
1242 kbSetModalityVoxelSize = true;
1243 }
1244 }
1245
1246 } // if(box.GetName() == volName)
1247
1248
1249 // processing geometry construction based on the interactive PS
1250 if(!kFlagProcessedInteractiveScorer) {
1251
1252
1253 // get the dimension of the geometry defined in G4VScoringMesh
1255 //if(!pScrMan) return;
1256 if(pScrMan) {
1257 G4ScoringBox * scoringBox
1258 = dynamic_cast<G4ScoringBox*>(pScrMan->FindMesh(volName));
1259 //if(scoringBox == NULL) return;
1260 if(scoringBox) {
1261
1262
1263
1264 G4int nVoxels[3];
1265 scoringBox->GetNumberOfSegments(nVoxels);
1266 // this order depends on the G4ScoringBox
1267 kNestedVolumeDimension[0] = nVoxels[2];
1268 kNestedVolumeDimension[1] = nVoxels[1];
1269 kNestedVolumeDimension[2] = nVoxels[0];
1270 kNestedVolumeDirAxis[0] = 2;
1271 kNestedVolumeDirAxis[1] = 1;
1272 kNestedVolumeDirAxis[2] = 0;
1273
1274 // get densities of the parameterized patient geometry
1275 for(G4int n0 = 0; n0 < kNestedVolumeDimension[0]; n0++) {
1276 for(G4int n1 = 0; n1 < kNestedVolumeDimension[1]; n1++) {
1277 for(G4int n2 = 0; n2 < kNestedVolumeDimension[2]; n2++) {
1278
1279 G4double dens = 0.*(g/cm3);
1280
1281 G4int idx[3];
1282 idx[kNestedVolumeDirAxis[0]] = n0;
1283 idx[kNestedVolumeDirAxis[1]] = n1;
1284 idx[kNestedVolumeDirAxis[2]] = n2;
1285 Index3D i3d(idx[0],idx[1],idx[2]);
1286 kNestedModality[i3d] = dens;
1287
1288 }
1289 }
1290 }
1291
1292 G4ThreeVector boxSize = scoringBox->GetSize();
1293 if(GFDEBUG_DET > 1) {
1294 G4cout << "Interactive Scorer : size - "
1295 << boxSize.x()/cm << " x "
1296 << boxSize.y()/cm << " x "
1297 << boxSize.z()/cm << " [cm3]" << G4endl;
1298 G4cout << "Interactive Scorer : # voxels - "
1299 << nVoxels[0] << " x "
1300 << nVoxels[1] << " x "
1301 << nVoxels[2] << G4endl;
1302 }
1303 kVolumeSize.set(boxSize.x()*2,
1304 boxSize.y()*2,
1305 boxSize.z()*2);
1306
1307 // mesh size
1308 if(!kbSetModalityVoxelSize) {
1309 G4float spacing[3] = {static_cast<G4float>(boxSize.x()*2/nVoxels[0]),
1310 static_cast<G4float>(boxSize.y()*2/nVoxels[1]),
1311 static_cast<G4float>(boxSize.z()*2/nVoxels[2])};
1312
1313 kgMocrenIO->setVoxelSpacing(spacing);
1314 kVoxelDimension.set(spacing[0], spacing[1], spacing[2]);
1315 kbSetModalityVoxelSize = true;
1316
1317 }
1318
1319
1320 kVolumeTrans3D = fObjectTransformation;
1321
1322 // translation for the scoring mesh
1323 G4ThreeVector sbth = scoringBox->GetTranslation();
1324 G4Translate3D sbtranslate(sbth);
1325 kVolumeTrans3D = kVolumeTrans3D*sbtranslate;
1326
1327 // rotation matrix for the scoring mesh
1328 G4RotationMatrix sbrm;
1329 sbrm = scoringBox->GetRotationMatrix();
1330 if(!sbrm.isIdentity()) {
1331 G4ThreeVector sbdummy(0.,0.,0.);
1332 G4Transform3D sbrotate(sbrm.inverse(), sbdummy);
1333 kVolumeTrans3D = kVolumeTrans3D*sbrotate;
1334 }
1335
1336
1337 // coordination system correction for gMocren
1338 G4ThreeVector raxisY(0., 1., 0.), dummyY(0.,0.,0.);
1339 G4RotationMatrix rotY(raxisY, pi*rad);
1340 G4Transform3D trotY(rotY, dummyY);
1341 G4ThreeVector raxisZ(0., 0., 1.), dummyZ(0.,0.,0.);
1342 G4RotationMatrix rotZ(raxisZ, pi*rad);
1343 G4Transform3D trotZ(rotZ, dummyZ);
1344
1345 kVolumeTrans3D = kVolumeTrans3D*trotY*trotZ;
1346
1347 }
1348 }
1349 //
1350 kFlagProcessedInteractiveScorer = true;
1351 }
1352
1353
1354 static G4VPhysicalVolume * volPV = NULL;
1355 if(pPVModel->GetCurrentPV()->GetName() == volName) {
1356 volPV = pPVModel->GetCurrentPV();
1357 }
1358
1359 //-- add detectors
1360 G4bool bAddDet = true;
1361 if(!kMessenger.getDrawVolumeGrid()) {
1362
1363 if(kFlagParameterization == 0) { // nested parameterisation
1364 /*
1365 G4String volDSolidName;
1366 if(volPV) {
1367 G4int nDaughter = volPV->GetLogicalVolume()->GetNoDaughters();
1368 G4VPhysicalVolume * volDPV = NULL;
1369 if(nDaughter > 0) volDPV = volPV->GetLogicalVolume()->GetDaughter(0);
1370 if(volDPV) {
1371 nDaughter = volDPV->GetLogicalVolume()->GetNoDaughters();
1372 if(nDaughter > 0)
1373 volDSolidName = volDPV->GetLogicalVolume()->GetDaughter(0)
1374 ->GetLogicalVolume()->GetSolid()->GetName();
1375 }
1376 }
1377 */
1378
1379 //std::cout << "Parameterization volume: " << volName << " - "
1380 // << box.GetName() << std::endl;
1381
1382 if(volName == box.GetName()) {
1383 bAddDet = false;
1384 }
1385
1386 std::vector<G4String>::iterator itr = kNestedVolumeNames.begin();
1387 for(; itr != kNestedVolumeNames.end(); itr++) {
1388 if(*itr == box.GetName()) {
1389 bAddDet = false;
1390 break;
1391 }
1392 }
1393 } else if(kFlagParameterization == 1) { // phantom paramemterisation
1394
1395 G4String volDSolidName;
1396 if(volPV) {
1397 volDSolidName = volPV->GetLogicalVolume()->GetDaughter(0)
1399 }
1400
1401 //std::cout << "Phantom Parameterization volume: " << volDSolidName
1402 // << " - " << box.GetName() << std::endl;
1403
1404 if(volDSolidName == box.GetName()) {
1405 bAddDet = false;
1406 }
1407
1408 } else if(kFlagParameterization == 2) { // interactive primitive scorer
1409 //std::cout << "Regular Parameterization volume: " << box.GetName() << std::endl;
1410 }
1411
1412 }
1413 if(bAddDet) AddDetector(box);
1414
1415
1416} // void G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )
1417
1418
1419//----- Add tubes
1420void
1422{
1424 G4cout << "***** AddSolid ( tubes )" << G4endl;
1425
1426 //----- skip drawing invisible primitive
1427 if( !IsVisible() ) { return ; }
1428
1429 //----- Initialize if necessary
1431
1432 //
1433 AddDetector(tubes);
1434
1435
1436 // for a debug
1437 if(GFDEBUG_DET > 0) {
1438 G4cout << "-------" << G4endl;
1439 G4cout << " " << tubes.GetName() << G4endl;
1440 G4Polyhedron * poly = tubes.CreatePolyhedron();
1441 G4int nv = poly->GetNoVertices();
1442 for(G4int i = 0; i < nv; i++) {
1443 G4cout << " (" << poly->GetVertex(i).x() << ", "
1444 << poly->GetVertex(i).y() << ", "
1445 << poly->GetVertex(i).z() << ")" << G4endl;
1446 }
1447 delete poly;
1448 }
1449
1450 const G4VModel* pv_model = GetModel();
1451 if (!pv_model) { return ; }
1452 G4PhysicalVolumeModel* pPVModel =
1453 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
1454 if (!pPVModel) { return ; }
1455 G4Material * mat = pPVModel->GetCurrentMaterial();
1456 G4String name = mat->GetName();
1457
1458} // void G4GMocrenFileSceneHandler::AddSolid( const G4Tubs& )
1459
1460
1461
1462//----- Add cons
1463void
1465{
1467 G4cout << "***** AddSolid ( cons )" << G4endl;
1468
1469 //----- skip drawing invisible primitive
1470 if( !IsVisible() ) { return ; }
1471
1472 //----- Initialize if necessary
1474
1475 //
1476 AddDetector(cons);
1477
1478}// G4GMocrenFileSceneHandler::AddSolid( cons )
1479
1480
1481//----- Add trd
1483{
1485 G4cout << "***** AddSolid ( trd )" << G4endl;
1486
1487
1488 //----- skip drawing invisible primitive
1489 if( !IsVisible() ) { return ; }
1490
1491 //----- Initialize if necessary
1493
1494 //
1495 AddDetector(trd);
1496
1497} // G4GMocrenFileSceneHandler::AddSolid ( trd )
1498
1499
1500//----- Add sphere
1502{
1504 G4cout << "***** AddSolid ( sphere )" << G4endl;
1505
1506 //----- skip drawing invisible primitive
1507 if( !IsVisible() ) { return ; }
1508
1509 //----- Initialize if necessary
1511
1512 //
1513 AddDetector(sphere);
1514
1515} // G4GMocrenFileSceneHandler::AddSolid ( sphere )
1516
1517
1518//----- Add para
1520{
1522 G4cout << "***** AddSolid ( para )" << G4endl;
1523
1524 //----- skip drawing invisible primitive
1525 if( !IsVisible() ) { return ; }
1526
1527 //----- Initialize if necessary
1529
1530 //
1531 AddDetector(para);
1532
1533} // G4GMocrenFileSceneHandler::AddSolid ( para )
1534
1535
1536//----- Add trap
1538{
1540 G4cout << "***** AddSolid ( trap )" << G4endl;
1541
1542 //----- skip drawing invisible primitive
1543 if( !IsVisible() ) { return ; }
1544
1545 //----- Initialize if necessary
1547
1548 //
1549 AddDetector(trap);
1550
1551} // G4GMocrenFileSceneHandler::AddSolid (const G4Trap& trap)
1552
1553
1554//----- Add torus
1555void
1557{
1559 G4cout << "***** AddSolid ( torus )" << G4endl;
1560
1561 //----- skip drawing invisible primitive
1562 if( !IsVisible() ) { return ; }
1563
1564 //----- Initialize if necessary
1566
1567 //
1568 AddDetector(torus);
1569
1570} // void G4GMocrenFileSceneHandler::AddSolid( const G4Torus& )
1571
1572
1573
1574//----- Add a shape which is not treated above
1576{
1577 //----- skip drawing invisible primitive
1578 if( !IsVisible() ) { return ; }
1579
1580 //----- Initialize if necessary
1582
1583 //
1584 AddDetector(solid);
1585
1586 //----- Send a primitive
1587 G4VSceneHandler::AddSolid( solid ) ;
1588
1589} //G4GMocrenFileSceneHandler::AddSolid ( const G4VSolid& )
1590
1591#include "G4TrajectoriesModel.hh"
1592#include "G4VTrajectory.hh"
1593#include "G4VTrajectoryPoint.hh"
1594
1595//----- Add a trajectory
1597
1598 kbModelingTrajectory = true;
1599
1601
1602 if(GFDEBUG_TRK) {
1603 G4cout << " ::AddCompound(const G4VTrajectory&) >>>>>>>>> " << G4endl;
1604 G4TrajectoriesModel * pTrModel = dynamic_cast<G4TrajectoriesModel*>(fpModel);
1605 if (!pTrModel) {
1607 ("G4VSceneHandler::AddCompound(const G4VTrajectory&)",
1608 "gMocren0013", FatalException, "Not a G4TrajectoriesModel.");
1609 } else {
1610 traj.DrawTrajectory();
1611
1612 const G4VTrajectory * trj = pTrModel->GetCurrentTrajectory();
1613 G4cout << "------ track" << G4endl;
1614 G4cout << " name: " << trj->GetParticleName() << G4endl;
1615 G4cout << " id: " << trj->GetTrackID() << G4endl;
1616 G4cout << " charge: " << trj->GetCharge() << G4endl;
1617 G4cout << " momentum: " << trj->GetInitialMomentum() << G4endl;
1618
1619 G4int nPnt = trj->GetPointEntries();
1620 G4cout << " point: ";
1621 for(G4int i = 0; i < nPnt; i++) {
1622 G4cout << trj->GetPoint(i)->GetPosition() << ", ";
1623 }
1624 G4cout << G4endl;
1625 }
1626 G4cout << G4endl;
1627 }
1628
1629 kbModelingTrajectory = false;
1630}
1631
1632#include <vector>
1633#include "G4VHit.hh"
1634#include "G4AttValue.hh"
1635//----- Add a hit
1637 if(GFDEBUG_HIT) G4cout << " ::AddCompound(const G4VHit&) >>>>>>>>> " << G4endl;
1638
1640
1641 /*
1642 const std::map<G4String, G4AttDef> * map = hit.GetAttDefs();
1643 if(!map) return;
1644 std::map<G4String, G4AttDef>::const_iterator itr = map->begin();
1645 for(; itr != map->end(); itr++) {
1646 G4cout << itr->first << " : " << itr->second.GetName()
1647 << " , " << itr->second.GetDesc() << G4endl;
1648 }
1649 */
1650
1651 std::vector<G4String> hitNames = kMessenger.getHitNames();
1652 if(GFDEBUG_HIT) {
1653 std::vector<G4String>::iterator itr = hitNames.begin();
1654 for(; itr != hitNames.end(); itr++)
1655 G4cout << " hit name : " << *itr << G4endl;
1656 }
1657
1658 std::vector<G4AttValue> * attval = hit.CreateAttValues();
1659 if(!attval) {G4cout << "0 empty " << G4endl;}
1660 else {
1661
1662 G4bool bid[3] = {false, false, false};
1663 Index3D id;
1664
1665 std::vector<G4AttValue>::iterator itr;
1666 // First, get IDs
1667 for(itr = attval->begin(); itr != attval->end(); itr++) {
1668 std::string stmp = itr->GetValue();
1669 std::istringstream sval(stmp.c_str());
1670
1671 if(itr->GetName() == G4String("XID")) {
1672 sval >> id.x;
1673 bid[0] = true;
1674 continue;
1675 }
1676 if(itr->GetName() == G4String("YID")) {
1677 sval >> id.y;
1678 bid[1] = true;
1679 continue;
1680 }
1681 if(itr->GetName() == G4String("ZID")) {
1682 sval >> id.z;
1683 bid[2] = true;
1684 continue;
1685 }
1686 }
1687
1688 G4int nhitname = (G4int)hitNames.size();
1689
1690 if(bid[0] && bid[1] && bid[2]) {
1691
1692 if(GFDEBUG_HIT)
1693 G4cout << " Hit : index(" << id.x << ", " << id.y << ", "
1694 << id.z << ")" << G4endl;
1695
1696 // Get attributes
1697 for(itr = attval->begin(); itr != attval->end(); itr++) {
1698 for(G4int i = 0; i < nhitname; i++) {
1699 if(itr->GetName() == hitNames[i]) {
1700
1701 std::string stmp = itr->GetValue();
1702 std::istringstream sval(stmp.c_str());
1703 G4double value;
1704 G4String unit;
1705 sval >> value >> unit;
1706
1707 std::map<G4String, std::map<Index3D, G4double> >::iterator kNestedHitsListItr;
1708 kNestedHitsListItr = kNestedHitsList.find(hitNames[i]);
1709 if(kNestedHitsListItr != kNestedHitsList.end()) {
1710 //fTempNestedHits = &kNestedHitsListItr->second;
1711 //(*fTempNestedHits)[id] = value;
1712 kNestedHitsListItr->second[id] = value;
1713 } else {
1714 std::map<Index3D, G4double> hits;
1715 hits.insert(std::map<Index3D, G4double>::value_type(id, value));
1716 kNestedHitsList[hitNames[i]] = hits;
1717 }
1718
1719
1720 if(GFDEBUG_HIT)
1721 G4cout << " : " << hitNames[i] << " -> " << value
1722 << " [" << unit << "]" << G4endl;
1723 }
1724 }
1725 }
1726 } else {
1727 G4Exception("G4GMocrenFileSceneHandler::AddCompound(const G4VHit &)",
1728 "gMocren0014", FatalException, "Error");
1729 }
1730
1731 delete attval;
1732 }
1733
1734}
1735
1737 if(GFDEBUG_DIGI) G4cout << " ::AddCompound(const G4VDigi&) >>>>>>>>> " << G4endl;
1739}
1740
1742 if(GFDEBUG_HIT)
1743 G4cout << " ::AddCompound(const std::map<G4int, G4double*> &) >>>>>>>>> " << G4endl;
1744
1745
1746 std::vector<G4String> hitScorerNames = kMessenger.getHitScorerNames();
1747 G4int nhitname = (G4int)hitScorerNames.size();
1748 G4String scorername = static_cast<G4VHitsCollection>(hits).GetName();
1749
1750 //-- --//
1751 /*
1752 std::vector<G4String> hitScorerNames = kMessenger.getHitScorerNames();
1753 if(GFDEBUG_HIT) {
1754 std::vector<G4String>::iterator itr = hitScorerNames.begin();
1755 for(; itr != hitScorerNames.end(); itr++)
1756 G4cout << " PS name : " << *itr << G4endl;
1757 }
1758 */
1759
1760 { // Scope bracket to avoid compiler messages about shadowing (JA).
1761 //for(G4int i = 0; i < nhitname; i++) { // this selection trusts
1762 //if(scorername == hitScorerNames[i]) { // thea command /vis/scene/add/psHits hit_name.
1763
1764 G4int idx[3];
1765 std::map<G4int, G4double*> * map = hits.GetMap();
1766 std::map<G4int, G4double*>::const_iterator itr = map->begin();
1767 for(; itr != map->end(); itr++) {
1768 GetNestedVolumeIndex(itr->first, idx);
1769 Index3D id(idx[0], idx[1], idx[2]);
1770
1771 std::map<G4String, std::map<Index3D, G4double> >::iterator nestedHitsListItr;
1772 nestedHitsListItr = kNestedHitsList.find(scorername);
1773 if(nestedHitsListItr != kNestedHitsList.end()) {
1774 nestedHitsListItr->second[id] = *(itr->second);
1775 } else {
1776 std::map<Index3D, G4double> hit;
1777 hit.insert(std::map<Index3D, G4double>::value_type(id, *(itr->second)));
1778 kNestedHitsList[scorername] = hit;
1779 }
1780 }
1781
1782 //break;
1783 //}
1784 //}
1785 }
1786
1787 if(GFDEBUG_HIT) {
1788 G4String meshname = static_cast<G4VHitsCollection>(hits).GetSDname();
1789 G4cout << " >>>>> " << meshname << " : " << scorername << G4endl;
1790
1791 for(G4int i = 0; i < nhitname; i++)
1792 if(scorername == hitScorerNames[i])
1793 G4cout << " !!!! Hit scorer !!!! " << scorername << G4endl;
1794
1795 G4cout << " dimension: "
1796 << kNestedVolumeDimension[0] << " x "
1797 << kNestedVolumeDimension[1] << " x "
1798 << kNestedVolumeDimension[2] << G4endl;
1799
1800 G4int id[3];
1801 std::map<G4int, G4double*> * map = hits.GetMap();
1802 std::map<G4int, G4double*>::const_iterator itr = map->begin();
1803 for(; itr != map->end(); itr++) {
1804 GetNestedVolumeIndex(itr->first, id);
1805 G4cout << "[" << itr->first << "] "
1806 << "("<< id[0] << "," << id[1] << "," << id[2] << ")"
1807 << *(itr->second) << ", ";
1808 }
1809 G4cout << G4endl;
1810 }
1811}
1812
1814 if(GFDEBUG_HIT)
1815 G4cout << " ::AddCompound(const std::map<G4int, G4StatDouble*> &) >>>>>>>>> " << G4endl;
1816
1817
1818 std::vector<G4String> hitScorerNames = kMessenger.getHitScorerNames();
1819 G4int nhitname = (G4int)hitScorerNames.size();
1820 G4String scorername = static_cast<G4VHitsCollection>(hits).GetName();
1821
1822 //-- --//
1823 /*
1824 std::vector<G4String> hitScorerNames = kMessenger.getHitScorerNames();
1825 if(GFDEBUG_HIT) {
1826 std::vector<G4String>::iterator itr = hitScorerNames.begin();
1827 for(; itr != hitScorerNames.end(); itr++)
1828 G4cout << " PS name : " << *itr << G4endl;
1829 }
1830 */
1831
1832 { // Scope bracket to avoid compiler messages about shadowing (JA).
1833 //for(G4int i = 0; i < nhitname; i++) { // this selection trusts
1834 //if(scorername == hitScorerNames[i]) { // thea command /vis/scene/add/psHits hit_name.
1835
1836 G4int idx[3];
1837 std::map<G4int, G4StatDouble*> * map = hits.GetMap();
1838 std::map<G4int, G4StatDouble*>::const_iterator itr = map->begin();
1839 for(; itr != map->end(); itr++) {
1840 GetNestedVolumeIndex(itr->first, idx);
1841 Index3D id(idx[0], idx[1], idx[2]);
1842
1843 std::map<G4String, std::map<Index3D, G4double> >::iterator nestedHitsListItr;
1844 nestedHitsListItr = kNestedHitsList.find(scorername);
1845 if(nestedHitsListItr != kNestedHitsList.end()) {
1846 nestedHitsListItr->second[id] = itr->second->sum_wx();
1847 } else {
1848 std::map<Index3D, G4double> hit;
1849 hit.insert(std::map<Index3D, G4double>::value_type(id, itr->second->sum_wx()));
1850 kNestedHitsList[scorername] = hit;
1851 }
1852 }
1853
1854 //break;
1855 //}
1856 //}
1857 }
1858
1859 if(GFDEBUG_HIT) {
1860 G4String meshname = static_cast<G4VHitsCollection>(hits).GetSDname();
1861 G4cout << " >>>>> " << meshname << " : " << scorername << G4endl;
1862
1863 for(G4int i = 0; i < nhitname; i++)
1864 if(scorername == hitScorerNames[i])
1865 G4cout << " !!!! Hit scorer !!!! " << scorername << G4endl;
1866
1867 G4cout << " dimension: "
1868 << kNestedVolumeDimension[0] << " x "
1869 << kNestedVolumeDimension[1] << " x "
1870 << kNestedVolumeDimension[2] << G4endl;
1871
1872 G4int id[3];
1873 std::map<G4int, G4StatDouble*> * map = hits.GetMap();
1874 std::map<G4int, G4StatDouble*>::const_iterator itr = map->begin();
1875 for(; itr != map->end(); itr++) {
1876 GetNestedVolumeIndex(itr->first, id);
1877 G4cout << "[" << itr->first << "] "
1878 << "("<< id[0] << "," << id[1] << "," << id[2] << ")"
1879 << itr->second->sum_wx() << ", ";
1880 }
1881 G4cout << G4endl;
1882 }
1883}
1884
1885//-----
1886G4bool G4GMocrenFileSceneHandler::IsVisible()
1887{
1888 //-----
1889 G4bool visibility = true ;
1890
1891 const G4VisAttributes* pVisAttribs =
1893
1894 if(pVisAttribs) {
1895 visibility = pVisAttribs->IsVisible();
1896 }
1897
1898 return visibility ;
1899
1900} // G4GMocrenFileSceneHandler::IsVisible()
1901
1902
1903//-----
1905{
1906 // This is typically called after an update and before drawing hits
1907 // of the next event. To simulate the clearing of "transients"
1908 // (hits, etc.) the detector is redrawn...
1909 if (fpViewer) {
1910 fpViewer -> SetView ();
1911 fpViewer -> ClearView ();
1912 fpViewer -> DrawView ();
1913 }
1914}
1915
1916//-----
1917void G4GMocrenFileSceneHandler::AddDetector(const G4VSolid & solid) {
1918
1919 Detector detector;
1920
1921 // detector name
1922 detector.name = solid.GetName();
1923 if(GFDEBUG_DET > 1)
1924 G4cout << "0 Detector name : " << detector.name << G4endl;
1925
1926 const G4VModel* pv_model = GetModel();
1927 if (!pv_model) { return ; }
1928 G4PhysicalVolumeModel* pPVModel =
1929 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
1930 if (!pPVModel) { return ; }
1931
1932 // edge points of the detector
1933 std::vector<G4float *> dedges;
1934 G4Polyhedron * poly = solid.CreatePolyhedron();
1935 detector.polyhedron = poly;
1936 detector.transform3D = fObjectTransformation;
1937
1938 // retrieve color
1939 unsigned char uccolor[3] = {30, 30, 30};
1940 if(pPVModel->GetCurrentLV()->GetVisAttributes()) {
1941 G4Color color = pPVModel->GetCurrentLV()->GetVisAttributes()->GetColor();
1942 uccolor[0] = (unsigned char)(color.GetRed()*255);
1943 uccolor[1] = (unsigned char)(color.GetGreen()*255);
1944 uccolor[2] = (unsigned char)(color.GetBlue()*255);
1945 //if(uccolor[0] < 2 && uccolor[1] < 2 && uccolor[2] < 2)
1946 //uccolor[0] = uccolor[1] = uccolor[2] = 30; // dark grey
1947 }
1948 for(G4int i = 0; i < 3; i++) detector.color[i] = uccolor[i];
1949 //
1950 kDetectors.push_back(detector);
1951
1952 if(GFDEBUG_DET > 1) {
1953 G4cout << "0 color: (" << (G4int)uccolor[0] << ", "
1954 << (G4int)uccolor[1] << ", " << (G4int)uccolor[2] << ")"
1955 << G4endl;
1956 }
1957
1958}
1959
1960//-----
1961void G4GMocrenFileSceneHandler::ExtractDetector() {
1962
1963 std::vector<Detector>::iterator itr = kDetectors.begin();
1964
1965 for(; itr != kDetectors.end(); itr++) {
1966
1967 // detector name
1968 G4String detname = itr->name;
1969 if(GFDEBUG_DET > 1)
1970 G4cout << "Detector name : " << detname << G4endl;
1971
1972 // edge points of the detector
1973 std::vector<G4float *> dedges;
1974 G4Polyhedron * poly = itr->polyhedron;
1975 poly->Transform(itr->transform3D);
1976 G4Transform3D invVolTrans = kVolumeTrans3D.inverse();
1977 poly->Transform(invVolTrans);
1978
1979 G4Point3D v1, v2;
1980 G4bool bnext = true;
1981 G4int next;
1982 G4int nedges = 0;
1983 //
1984 while(bnext) {
1985 if(!(poly->GetNextEdge(v1, v2, next))) bnext =false;
1986 G4float * edge = new G4float[6];
1987 edge[0] = v1.x()/mm;
1988 edge[1] = v1.y()/mm;
1989 edge[2] = v1.z()/mm;
1990 edge[3] = v2.x()/mm;
1991 edge[4] = v2.y()/mm;
1992 edge[5] = v2.z()/mm;
1993 dedges.push_back(edge);
1994 nedges++;
1995 }
1996 //delete poly;
1997 // detector color
1998 unsigned char uccolor[3] = {itr->color[0],
1999 itr->color[1],
2000 itr->color[2]};
2001 //
2002 kgMocrenIO->addDetector(detname, dedges, uccolor);
2003 for(G4int i = 0; i < nedges; i++) { // # of edges is 12.
2004 delete [] dedges[i];
2005 }
2006 dedges.clear();
2007
2008 if(GFDEBUG_DET > 1) {
2009 G4cout << " color: (" << (G4int)uccolor[0] << ", "
2010 << (G4int)uccolor[1] << ", " << (G4int)uccolor[2] << ")"
2011 << G4endl;
2012 }
2013 }
2014}
2015
2016void G4GMocrenFileSceneHandler::GetNestedVolumeIndex(G4int _idx, G4int _idx3d[3]) {
2017 if(kNestedVolumeDimension[0] == 0 ||
2018 kNestedVolumeDimension[1] == 0 ||
2019 kNestedVolumeDimension[2] == 0) {
2020 for(G4int i = 0; i < 3; i++) _idx3d[i] = 0;
2021 return;
2022 }
2023
2024
2025 if(kFlagParameterization == 0) {
2026
2027 G4int plane = kNestedVolumeDimension[2]*kNestedVolumeDimension[1];
2028 G4int line = kNestedVolumeDimension[2];
2029
2030 /*
2031 G4int idx3d[3];
2032 idx3d[0] = _idx/plane;
2033 idx3d[1] = (_idx%plane)/line;
2034 idx3d[2] = (_idx%plane)%line;
2035 _idx3d[0] = idx3d[kNestedVolumeDirAxis[0]];
2036 _idx3d[1] = idx3d[kNestedVolumeDirAxis[1]];
2037 _idx3d[2] = idx3d[kNestedVolumeDirAxis[2]];
2038 */
2039
2040 _idx3d[kNestedVolumeDirAxis[0]] = _idx/plane;
2041 _idx3d[kNestedVolumeDirAxis[1]] = (_idx%plane)/line;
2042 _idx3d[kNestedVolumeDirAxis[2]] = (_idx%plane)%line;
2043
2044
2045
2046 /*
2047
2048 G4cout << "G4GMocrenFileSceneHandler::GetNestedVolumeIndex : " << G4endl;
2049 G4cout << "(depi, depj, depk) : "
2050 << kNestedVolumeDirAxis[0] << ", "
2051 << kNestedVolumeDirAxis[1] << ", "
2052 << kNestedVolumeDirAxis[2] << G4endl;
2053 G4cout << "(ni, nj, nk) :"
2054 << kNestedVolumeDimension[0] << ", "
2055 << kNestedVolumeDimension[1] << ", "
2056 << kNestedVolumeDimension[2] << " - " << G4endl;
2057
2058 G4cout << " _idx = " << _idx << " : plane = "
2059 << plane << " , line = " << line << G4endl;
2060 G4cout << "(idx,idy,idz) + " << _idx3d[0] << ", "
2061 << _idx3d[1] << ", " << _idx3d[2] << " + " << G4endl;
2062
2063 */
2064
2065
2066
2067 } else {
2068
2069 G4int plane = kNestedVolumeDimension[0]*kNestedVolumeDimension[1];
2070 G4int line = kNestedVolumeDimension[0];
2071 _idx3d[kNestedVolumeDirAxis[2]] = _idx/plane;
2072 _idx3d[kNestedVolumeDirAxis[1]] = (_idx%plane)/line;
2073 _idx3d[kNestedVolumeDirAxis[0]] = (_idx%plane)%line;
2074
2075 }
2076
2077}
2078
2079
2080//-- --//
2081G4GMocrenFileSceneHandler::Detector::Detector()
2082 : polyhedron(0) {
2083 color[0] = color[1] = color[2] = 255;
2084}
2085G4GMocrenFileSceneHandler::Detector::~Detector() {
2086 if(!polyhedron) delete polyhedron;
2087}
2088void G4GMocrenFileSceneHandler::Detector::clear() {
2089 name.clear();
2090 if(!polyhedron) delete polyhedron;
2091 color[0] = color[1] = color[2] = 255;
2092 transform3D = G4Transform3D::Identity;
2093}
2094
2095//-- --//
2096G4GMocrenFileSceneHandler::Index3D::Index3D()
2097 : x(0), y(0), z(0) {
2098 ;
2099}
2100
2101G4GMocrenFileSceneHandler::Index3D::Index3D(const Index3D & _index3D)
2102 : x(_index3D.x), y(_index3D.y), z(_index3D.z) {
2103 //: x(_index3D.X()),
2104 //y(_index3D.Y()),
2105 //z(_index3D.Z()) {
2106 // : x(static_cast<Index3D>(_index3D).x),
2107 // y(static_cast<Index3D>(_index3D).y),
2108 // z(static_cast<Index3D>(_index3D).z) {
2109 ;
2110}
2111
2112G4GMocrenFileSceneHandler::Index3D::Index3D(G4int _x, G4int _y, G4int _z)
2113 : x(_x), y(_y), z(_z) {
2114 ;
2115}
2116G4bool G4GMocrenFileSceneHandler::Index3D::operator < (const Index3D & _right) const {
2117 if(z < static_cast<Index3D>(_right).z) {
2118 return true;
2119 } else if(z == _right.z) {
2120 if(y < static_cast<Index3D>(_right).y) return true;
2121 else if(y == _right.y)
2122 if(x < static_cast<Index3D>(_right).x) return true;
2123 }
2124 return false;
2125}
2126G4bool G4GMocrenFileSceneHandler::Index3D::operator == (const Index3D & _right) const {
2127 if(z == _right.z && y == _right.y && x == _right.x) return true;
2128 return false;
2129}
const int FR_MAX_FILE_NUM
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
const G4bool GFDEBUG
const G4bool GFDEBUG_TRK
const G4bool GFDEBUG_DIGI
const G4bool GFDEBUG_HIT
const G4int FR_MAX_FILE_NUM
const G4int GFDEBUG_DET
const char GDD_FILE_HEADER[]
const char DEFAULT_GDD_FILE_NAME[]
const G4int MAX_NUM_TRAJECTORIES
float G4float
Definition: G4Types.hh:84
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
Hep3Vector & transform(const HepRotation &)
Definition: ThreeVectorR.cc:20
void set(double x, double y, double z)
HepRotation inverse() const
bool isIdentity() const
Definition: Rotation.cc:167
Definition: G4Box.hh:56
G4double GetYHalfLength() const
G4double GetZHalfLength() const
G4double GetXHalfLength() const
G4Polyhedron * CreatePolyhedron() const
Definition: G4Box.cc:544
G4double GetBlue() const
Definition: G4Colour.hh:154
G4double GetRed() const
Definition: G4Colour.hh:152
G4double GetGreen() const
Definition: G4Colour.hh:153
Definition: G4Cons.hh:78
G4double GetDensity(G4int &_ct) const
G4GMocrenFileSceneHandler(G4GMocrenFile &system, G4GMocrenMessenger &messenger, const G4String &name="")
void AddCompound(const G4VTrajectory &traj)
virtual void BeginPrimitives(const G4Transform3D &objectTransformation)
void AddPrimitive(const G4Polyline &line)
void setModalityImageSize(int _size[3])
void setDoseDistUnit(std::string &_unit, int _num=0)
short convertDensityToHU(float &_dens)
void setDoseDistMinMax(short _minmax[2], int _num=0)
void setDoseDistScale(double &_scale, int _num=0)
void setModalityImageDensityMap(std::vector< float > &_map)
void translateDetector(std::vector< float > &_translate)
void clearDetector()
Definition: G4GMocrenIO.hh:451
void clearROIAll()
void setVoxelSpacing(float _spacing[3])
void setDoseDistName(std::string _name, int _num=0)
bool storeData(char *_filename)
Definition: G4GMocrenIO.cc:457
void setModalityImage(short *_image)
void setDoseDistSize(int _size[3], int _num=0)
void setModalityImageMinMax(short _minmax[2])
void newDoseDist()
void clearDoseDistAll()
void clearTracks()
Definition: G4GMocrenIO.hh:439
void addDetector(std::string &_name, std::vector< float * > &_det, unsigned char _color[3])
void setDoseDist(double *_image, int _num=0)
void translateTracks(std::vector< float > &_translateo)
void addTrack(float *_tracks)
virtual std::vector< G4String > getHitNames()
virtual std::vector< G4String > getHitScorerNames()
virtual G4bool getDrawVolumeGrid()
virtual G4String getVolumeName()
G4VSolid * GetSolid() const
const G4VisAttributes * GetVisAttributes() const
std::size_t GetNoDaughters() const
G4VPhysicalVolume * GetDaughter(const std::size_t i) const
G4double GetDensity() const
Definition: G4Material.hh:175
const G4String & GetName() const
Definition: G4Material.hh:172
Definition: G4Para.hh:79
virtual G4Material * ComputeMaterial(const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
G4double GetVoxelHalfZ() const
std::size_t GetNoVoxelsX() const
G4double GetVoxelHalfY() const
G4double GetVoxelHalfX() const
std::size_t GetNoVoxelsZ() const
std::size_t GetNoVoxelsY() const
const std::vector< G4PhysicalVolumeNodeID > & GetDrawnPVPath() const
G4VPhysicalVolume * GetCurrentPV() const
G4LogicalVolume * GetCurrentLV() const
G4VPhysicalVolume * GetTopPhysicalVolume() const
G4Material * GetCurrentMaterial() const
const std::vector< Model > & GetEndOfEventModelList() const
static G4ScoringManager * GetScoringManager()
G4VScoringMesh * FindMesh(G4VHitsCollection *map)
Definition: G4Text.hh:72
const G4VTrajectory * GetCurrentTrajectory() const
Definition: G4Trd.hh:63
Definition: G4Tubs.hh:75
G4Polyhedron * CreatePolyhedron() const
Definition: G4Tubs.cc:1759
Definition: G4VHit.hh:48
virtual std::vector< G4AttValue > * CreateAttValues() const
Definition: G4VHit.hh:64
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
virtual void ComputeTransformation(const G4int no, G4VPhysicalVolume *currentPV) const =0
virtual G4Material * ComputeMaterial(G4VPhysicalVolume *currentVol, const G4int repNo, const G4VTouchable *parentTouch=nullptr)=0
virtual G4bool IsNested() const
virtual G4bool IsReplicated() const =0
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetMultiplicity() const
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0
virtual G4int GetCopyNo() const =0
const G4String & GetName() const
virtual G4VPVParameterisation * GetParameterisation() const =0
G4ThreeVector GetObjectTranslation() const
virtual G4bool IsParameterised() const =0
virtual void BeginModeling()
G4VModel * GetModel() const
G4Scene * GetScene() const
G4Transform3D fObjectTransformation
virtual void EndPrimitives()
virtual void EndModeling()
G4VViewer * fpViewer
const G4String & GetName() const
const G4VisAttributes * fpVisAttribs
virtual void BeginPrimitives(const G4Transform3D &objectTransformation=G4Transform3D())
virtual void AddSolid(const G4Box &)
virtual void AddCompound(const G4VTrajectory &)
G4ThreeVector GetTranslation() const
G4ThreeVector GetSize() const
void GetNumberOfSegments(G4int nSegment[3])
G4RotationMatrix GetRotationMatrix() const
G4String GetName() const
virtual G4Polyhedron * CreatePolyhedron() const
Definition: G4VSolid.cc:700
virtual G4GeometryType GetEntityType() const =0
Map_t * GetMap() const
Definition: G4THitsMap.hh:155
virtual const G4ThreeVector GetPosition() const =0
virtual G4VTrajectoryPoint * GetPoint(G4int i) const =0
virtual G4int GetPointEntries() const =0
virtual G4String GetParticleName() const =0
virtual G4int GetTrackID() const =0
virtual G4ThreeVector GetInitialMomentum() const =0
virtual G4double GetCharge() const =0
virtual void DrawTrajectory() const
const G4VisAttributes * GetApplicableVisAttributes(const G4VisAttributes *) const
const G4Color & GetColor() const
G4bool IsVisible() const
static Verbosity GetVerbosity()
const G4VisAttributes * GetVisAttributes() const
static DLL_API const Transform3D Identity
Definition: Transform3D.h:196
CLHEP::HepRotation getRotation() const
CLHEP::Hep3Vector getTranslation() const
Transform3D inverse() const
Definition: Transform3D.cc:141
G4Point3D GetVertex(G4int index) const
HepPolyhedron & Transform(const G4Transform3D &t)
G4int GetNoFacets() const
G4bool GetNextVertexIndex(G4int &index, G4int &edgeFlag) const
G4int GetNoVertices() const
G4bool GetNextEdge(G4Point3D &p1, G4Point3D &p2, G4int &edgeFlag) const
EAxis
Definition: geomdefs.hh:54
@ kYAxis
Definition: geomdefs.hh:56
@ kXAxis
Definition: geomdefs.hh:55
@ kZAxis
Definition: geomdefs.hh:57
const char * name(G4int ptype)
#define DBL_MAX
Definition: templates.hh:62