source: trunk/source/visualization/gMocren/src/G4GMocrenFileSceneHandler.cc @ 1343

Last change on this file since 1343 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

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