00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038 #include <fstream>
00039 #include <cstdlib>
00040 #include <cstring>
00041 #include <sstream>
00042
00043 #include "globals.hh"
00044 #include "G4PhysicalConstants.hh"
00045 #include "G4SystemOfUnits.hh"
00046 #include "G4VisManager.hh"
00047
00048 #include "G4GMocrenFile.hh"
00049 #include "G4GMocrenFileSceneHandler.hh"
00050 #include "G4GMocrenFileViewer.hh"
00051 #include "G4Point3D.hh"
00052 #include "G4VisAttributes.hh"
00053 #include "G4Scene.hh"
00054 #include "G4Transform3D.hh"
00055 #include "G4Polyhedron.hh"
00056 #include "G4Box.hh"
00057 #include "G4Cons.hh"
00058 #include "G4Polyline.hh"
00059 #include "G4Trd.hh"
00060 #include "G4Tubs.hh"
00061 #include "G4Trap.hh"
00062 #include "G4Torus.hh"
00063 #include "G4Sphere.hh"
00064 #include "G4Para.hh"
00065 #include "G4Text.hh"
00066 #include "G4Circle.hh"
00067 #include "G4Square.hh"
00068 #include "G4VPhysicalVolume.hh"
00069 #include "G4PhysicalVolumeModel.hh"
00070 #include "G4LogicalVolume.hh"
00071 #include "G4Material.hh"
00072
00073 #include "G4VPVParameterisation.hh"
00074 #include "G4VVolumeMaterialScanner.hh"
00075 #include "G4VisTrajContext.hh"
00076 #include "G4TrajectoriesModel.hh"
00077 #include "G4VTrajectoryModel.hh"
00078 #include "G4TrajectoryDrawByCharge.hh"
00079 #include "G4HitsModel.hh"
00080 #include "G4GMocrenMessenger.hh"
00081
00082 #include "G4GMocrenIO.hh"
00083 #include "G4VNestedParameterisation.hh"
00084 #include "G4GMocrenTouchable.hh"
00085 #include "G4GMocrenFileCTtoDensityMap.hh"
00086 #include "G4PhantomParameterisation.hh"
00087
00088 #include "G4ScoringManager.hh"
00089 #include "G4ScoringBox.hh"
00090
00091
00092 const char GDD_FILE_HEADER [] = "g4_";
00093 const char DEFAULT_GDD_FILE_NAME[] = "g4_00.gdd";
00094
00095 const G4int FR_MAX_FILE_NUM = 100 ;
00096 const G4int MAX_NUM_TRAJECTORIES = 100000;
00097
00098
00099 const G4bool GFDEBUG = false;
00100 const G4bool GFDEBUG_TRK = false;
00101 const G4bool GFDEBUG_HIT = false;
00102 const G4bool GFDEBUG_DIGI = false;
00103 const G4int GFDEBUG_DET = 0;
00104
00106
00108
00109
00110 G4int G4GMocrenFileSceneHandler::kSceneIdCount = 0;
00111
00113
00115
00116
00117
00118 G4GMocrenFileSceneHandler::G4GMocrenFileSceneHandler(G4GMocrenFile& system,
00119 G4GMocrenMessenger & messenger,
00120 const G4String& name)
00121 : G4VSceneHandler(system, kSceneIdCount++, name),
00122 kSystem(system),
00123 kMessenger(messenger),
00124 kgMocrenIO(new G4GMocrenIO()),
00125 kbSetModalityVoxelSize(false),
00126 kbModelingTrajectory(false),
00127
00128 kFlagInModeling(false),
00129 kFlagSaving_g4_gdd(false),
00130 kFlagParameterization(0),
00131 kFLagProcessedInteractiveScorer(false) {
00132
00133
00134 if(getenv("G4GMocrenFile_DEST_DIR") == NULL) {
00135 kGddDestDir[0] = '\0';
00136
00137
00138 std::strncpy(kGddFileName, DEFAULT_GDD_FILE_NAME, std::strlen(DEFAULT_GDD_FILE_NAME));
00139 } else {
00140 const char * env = getenv("G4GMocrenFile_DEST_DIR");
00141 std::strncpy(kGddDestDir, env, std::strlen(env));
00142 std::strncpy(kGddFileName, DEFAULT_GDD_FILE_NAME, std::strlen(DEFAULT_GDD_FILE_NAME));
00143 }
00144
00145
00146 kMaxFileNum = FR_MAX_FILE_NUM ;
00147 if ( std::getenv( "G4GMocrenFile_MAX_FILE_NUM" ) != NULL ) {
00148
00149 std::sscanf( getenv("G4GMocrenFile_MAX_FILE_NUM"), "%d", &kMaxFileNum ) ;
00150
00151 } else {
00152 kMaxFileNum = FR_MAX_FILE_NUM ;
00153 }
00154 if( kMaxFileNum < 1 ) { kMaxFileNum = 1 ; }
00155
00156 InitializeParameters();
00157
00158 }
00159
00160
00161
00162 G4GMocrenFileSceneHandler::~G4GMocrenFileSceneHandler ()
00163 {
00164 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00165 G4cout << "***** ~G4GMocrenFileSceneHandler" << G4endl;
00166
00167 if(kGddDest) {
00168
00169
00170 GFEndModeling();
00171 }
00172 }
00173
00174
00175 void G4GMocrenFileSceneHandler::InitializeParameters() {
00176
00177 kbSetModalityVoxelSize = false;
00178
00179 for(G4int i = 0; i < 3; i++) {
00180 kModalitySize[i] = 0;
00181 kNestedVolumeDimension[i] = 0;
00182 kNestedVolumeDirAxis[i] = -1;
00183 }
00184
00185
00186
00187 }
00188
00189
00190 void G4GMocrenFileSceneHandler::SetGddFileName()
00191 {
00192
00193 const G4int MAX_FILE_INDEX = kMaxFileNum - 1 ;
00194
00195
00196 std::strncpy(kGddFileName, kGddDestDir, std::strlen(kGddDestDir));
00197
00198
00199 std::strncat ( kGddFileName, DEFAULT_GDD_FILE_NAME, std::strlen(DEFAULT_GDD_FILE_NAME));
00200
00201
00202 static G4int currentNumber = 0;
00203 for( G4int i = currentNumber ; i < kMaxFileNum ; i++) {
00204
00205
00206 if( i == MAX_FILE_INDEX )
00207 {
00208 if (G4VisManager::GetVerbosity() >= G4VisManager::warnings) {
00209 G4cout << "===========================================" << G4endl;
00210 G4cout << "WARNING MESSAGE from GMocrenFile driver: " << G4endl;
00211 G4cout << " This file name is the final one in the " << G4endl;
00212 G4cout << " automatic updation of the output file name." << G4endl;
00213 G4cout << " You may overwrite existing files, i.e. " << G4endl;
00214 G4cout << " g4_XX.gdd." << G4endl;
00215 G4cout << "===========================================" << G4endl;
00216 }
00217 }
00218
00219
00220 if( i >= 0 && i <= 9 ) {
00221 std::sprintf( kGddFileName, "%s%s%s%d.gdd" , kGddDestDir, GDD_FILE_HEADER, "0", i );
00222 } else {
00223 std::sprintf( kGddFileName, "%s%s%d.gdd" , kGddDestDir, GDD_FILE_HEADER, i );
00224 }
00225
00226
00227 std::ifstream fin(kGddFileName);
00228 if(GFDEBUG)
00229 G4cout << "FILEOPEN: " << i << " : " << kGddFileName << fin.fail() << G4endl;
00230 if(!fin) {
00231
00232 fin.close();
00233 currentNumber = i+1;
00234 break;
00235 } else {
00236
00237 fin.close();
00238 }
00239
00240 }
00241
00242 G4cout << "======================================================================" << G4endl;
00243 G4cout << "Output file: " << kGddFileName << G4endl;
00244 G4cout << "Destination directory (current dir if NULL): " << kGddDestDir << G4endl;
00245 G4cout << "Maximum number of files in the destination directory: " << kMaxFileNum << G4endl;
00246 G4cout << "Note:" << G4endl;
00247 G4cout << " * The maximum number is customizable as: " << G4endl;
00248 G4cout << " % setenv G4GMocrenFile_MAX_FILE_NUM number " << G4endl;
00249 G4cout << " * The destination directory is customizable as:" << G4endl;
00250 G4cout << " % setenv G4GMocrenFile_DEST_DIR dir_name/ " << G4endl;
00251 G4cout << " ** Do not forget \"/\" at the end of the dir_name, e.g. \"./tmp/\"." << G4endl;
00252
00253 G4cout << G4endl;
00254 G4cout << "Maximum number of trajectories is set to " << MAX_NUM_TRAJECTORIES << "."<< G4endl;
00255 G4cout << "======================================================================" << G4endl;
00256
00257 }
00258
00259
00260
00261 void G4GMocrenFileSceneHandler::BeginSavingGdd( void )
00262 {
00263 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00264 G4cout << "***** BeginSavingGdd (called)" << G4endl;
00265
00266 if( !IsSavingGdd() ) {
00267
00268 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations) {
00269 G4cout << "***** (started) " ;
00270 G4cout << "(open g4.gdd, ##)" << G4endl;
00271 }
00272
00273 SetGddFileName() ;
00274 kFlagSaving_g4_gdd = true;
00275
00276
00277 G4GMocrenFileCTtoDensityMap ctdens;
00278 short minmax[2];
00279 minmax[0] = ctdens.GetMinCT();
00280 minmax[1] = ctdens.GetMaxCT();
00281 kgMocrenIO->setModalityImageMinMax(minmax);
00282 std::vector<G4float> map;
00283 G4float dens;
00284 for(G4int i = minmax[0]; i <= minmax[1]; i++) {
00285 dens = ctdens.GetDensity(i);
00286 map.push_back(dens);
00287 }
00288 kgMocrenIO->setModalityImageDensityMap(map);
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315
00316 kgMocrenIO->clearDoseDistAll();
00317 kgMocrenIO->clearROIAll();
00318 kgMocrenIO->clearTracks();
00319 kgMocrenIO->clearDetector();
00320 std::vector<Detector>::iterator itr = kDetectors.begin();
00321 for(; itr != kDetectors.end(); itr++) {
00322 itr->clear();
00323 }
00324 kDetectors.clear();
00325
00326 kNestedHitsList.clear();
00327 kNestedVolumeNames.clear();
00328
00329 }
00330 }
00331
00332 void G4GMocrenFileSceneHandler::EndSavingGdd ( void )
00333 {
00334 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00335 G4cout << "***** EndSavingGdd (called)" << G4endl;
00336
00337 if(IsSavingGdd()) {
00338 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00339 G4cout << "***** (started) (close "
00340 << kGddFileName << ")" << G4endl;
00341
00342 if(kGddDest) kGddDest.close();
00343 kFlagSaving_g4_gdd = false;
00344
00345 std::map<Index3D, G4float>::iterator itr = kNestedModality.begin();
00346 G4int xmax=0, ymax=0, zmax=0;
00347 for(; itr != kNestedModality.end(); itr++) {
00348 if(itr->first.x > xmax) xmax = itr->first.x;
00349 if(itr->first.y > ymax) ymax = itr->first.y;
00350 if(itr->first.z > zmax) zmax = itr->first.z;
00351 }
00352
00353 kModalitySize[0] = xmax+1;
00354 kModalitySize[1] = ymax+1;
00355 kModalitySize[2] = zmax+1;
00356 kgMocrenIO->setModalityImageSize(kModalitySize);
00357 if(GFDEBUG) G4cout << "gMocren-file driver : modality size : "
00358 << kModalitySize[0] << " x "
00359 << kModalitySize[1] << " x "
00360 << kModalitySize[2] << G4endl;
00361
00362 G4int nxy = kModalitySize[0]*kModalitySize[1];
00363
00364 for(G4int z = 0; z < kModalitySize[2]; z++) {
00365 short * modality = new short[nxy];
00366 for(G4int y = 0; y < kModalitySize[1]; y++) {
00367 for(G4int x = 0; x < kModalitySize[0]; x++) {
00368
00369
00370
00371 G4int ixy = x + y*kModalitySize[0];
00372 Index3D idx(x,y,z);
00373 itr = kNestedModality.find(idx);
00374 if(itr != kNestedModality.end()) {
00375
00376 modality[ixy] = kgMocrenIO->convertDensityToHU(itr->second);
00377 } else {
00378 modality[ixy] = -1024;
00379 }
00380
00381 }
00382 }
00383 kgMocrenIO->setModalityImage(modality);
00384 }
00385
00386
00387 size_t nhits = kNestedHitsList.size();
00388 if(GFDEBUG) G4cout << "gMocren-file driver : # hits = " << nhits << G4endl;
00389
00390 std::map<Index3D, G4double>::iterator hitsItr;
00391 std::map<G4String, std::map<Index3D, G4double> >::iterator hitsListItr = kNestedHitsList.begin();
00392
00393 for(G4int n = 0; hitsListItr != kNestedHitsList.end(); hitsListItr++, n++) {
00394
00395 kgMocrenIO->newDoseDist();
00396 kgMocrenIO->setDoseDistName(hitsListItr->first, n);
00397 kgMocrenIO->setDoseDistSize(kModalitySize, n);
00398
00399 G4double minmax[2] = {DBL_MAX, -DBL_MAX};
00400 for(G4int z = 0 ; z < kModalitySize[2]; z++) {
00401 G4double * values = new G4double[nxy];
00402 for(G4int y = 0; y < kModalitySize[1]; y++) {
00403 for(G4int x = 0; x < kModalitySize[0]; x++) {
00404
00405 G4int ixy = x + y*kModalitySize[0];
00406 Index3D idx(x,y,z);
00407 hitsItr = hitsListItr->second.find(idx);
00408 if(hitsItr != hitsListItr->second.end()) {
00409
00410 values[ixy] = hitsItr->second;
00411 } else {
00412 values[ixy] = 0.;
00413 }
00414 if(values[ixy] < minmax[0]) minmax[0] = values[ixy];
00415 if(values[ixy] > minmax[1]) minmax[1] = values[ixy];
00416 }
00417 }
00418 kgMocrenIO->setDoseDist(values, n);
00419 }
00420 kgMocrenIO->setDoseDistMinMax(minmax, n);
00421 G4double lower = 0.;
00422 if(minmax[0] < 0) lower = minmax[0];
00423 G4double scale = (minmax[1]-lower)/25000.;
00424 kgMocrenIO->setDoseDistScale(scale, n);
00425 G4String sunit("unit?");
00426 kgMocrenIO->setDoseDistUnit(sunit, n);
00427 }
00428
00429
00430
00431 if(false) {
00432 G4ThreeVector trans;
00433 G4RotationMatrix rot;
00434 trans = kVolumeTrans3D.getTranslation();
00435 rot = kVolumeTrans3D.getRotation().inverse();
00436
00437 std::vector<G4float *> tracks;
00438 unsigned char colors[3];
00439 G4float * trk = new G4float[6];
00440 tracks.push_back(trk);
00441
00442 G4ThreeVector orig(0.,0.,0), xa(2000.,0.,0.), ya(0.,2000.,0.), za(0.,0.,2000.);
00443 orig -= trans;
00444 orig.transform(rot);
00445 xa -= trans;
00446 xa.transform(rot);
00447 ya -= trans;
00448 ya.transform(rot);
00449 za -= trans;
00450 za.transform(rot);
00451 for(G4int i = 0; i < 3; i++) trk[i] = orig[i];
00452 for(G4int i = 0; i < 3; i++) trk[i+3] = xa[i];
00453 colors[0] = 255; colors[1] = 0; colors[2] = 0;
00454 kgMocrenIO->addTrack(tracks, colors);
00455
00456 for(G4int i = 0; i < 3; i++) trk[i+3] = ya[i];
00457 colors[0] = 0; colors[1] = 255; colors[2] = 0;
00458 kgMocrenIO->addTrack(tracks, colors);
00459
00460 for(G4int i = 0; i < 3; i++) trk[i+3] = za[i];
00461 colors[0] = 0; colors[1] = 0; colors[2] = 255;
00462 kgMocrenIO->addTrack(tracks, colors);
00463 }
00464
00465
00466 ExtractDetector();
00467
00468
00469 if(GFDEBUG_DET) G4cout << ">>>>>>>>>>>>>>>>>>>>>> (";
00470 std::vector<G4float> transformObjects;
00471 for(G4int i = 0; i < 3; i++) {
00472
00473 transformObjects.push_back((kVolumeSize[i]/2. - kVoxelDimension[i]/2.));
00474 if(GFDEBUG_DET) G4cout << transformObjects[i] << ", ";
00475 }
00476 if(GFDEBUG_DET) G4cout << ")" << G4endl;
00477
00478
00479 kgMocrenIO->translateTracks(transformObjects);
00480 kgMocrenIO->translateDetector(transformObjects);
00481
00482
00483 kgMocrenIO->storeData(kGddFileName);
00484 }
00485
00486 }
00487
00488
00489
00490 void G4GMocrenFileSceneHandler::GFBeginModeling( void )
00491 {
00492 G4VSceneHandler::BeginModeling();
00493
00494 if( !GFIsInModeling() ) {
00495
00496
00497 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00498 G4cout << "***** G4GMocrenFileSceneHandler::GFBeginModeling (called & started)" << G4endl;
00499
00500
00501 BeginSavingGdd();
00502
00503 kFlagInModeling = true ;
00504
00505
00506
00507
00508
00509 if(GFDEBUG_HIT) {
00510 G4Scene * scene = GetScene();
00511 std::vector<G4Scene::Model> vmodel = scene->GetEndOfEventModelList();
00512 std::vector<G4Scene::Model>::iterator itr = vmodel.begin();
00513 for(; itr != vmodel.end(); itr++) {
00514 G4cout << " IIIIII model name: " << itr->fpModel->GetGlobalTag() << G4endl;
00515 }
00516 }
00517
00518 }
00519 }
00520
00521
00522
00523
00524
00525 void G4GMocrenFileSceneHandler::AddPrimitive (const G4Polyline& polyline)
00526 {
00527 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00528 G4cout << "***** AddPrimitive" << G4endl;
00529
00530 if (fProcessing2D) {
00531 static G4bool warned = false;
00532 if (!warned) {
00533 warned = true;
00534 G4Exception
00535 ("G4GMocrenFileSceneHandler::AddPrimitive (const G4Polyline&)",
00536 "gMocren1001", JustWarning,
00537 "2D polylines not implemented. Ignored.");
00538 }
00539 return;
00540 }
00541
00542
00543 GFBeginModeling();
00544
00545 static G4int numTrajectories = 0;
00546 if(numTrajectories >= MAX_NUM_TRAJECTORIES) return;
00547
00548
00549 if(kbModelingTrajectory) {
00550
00551 G4TrajectoriesModel * pTrModel = dynamic_cast<G4TrajectoriesModel*>(fpModel);
00552 if (!pTrModel) {
00553 G4Exception
00554 ("G4VSceneHandler::AddCompound(const G4Polyline&)",
00555 "gMocren0002", FatalException, "Not a G4TrajectoriesModel.");
00556 }
00557
00558 G4ThreeVector trans;
00559 G4RotationMatrix rot;
00560 trans = kVolumeTrans3D.getTranslation();
00561 rot = kVolumeTrans3D.getRotation().inverse();
00562
00563 if(GFDEBUG_TRK) G4cout << " trajectory points : " << G4endl;
00564 std::vector<G4float *> trajectory;
00565 if(polyline.size() < 2) return;
00566 G4Polyline::const_iterator preitr = polyline.begin();
00567 G4Polyline::const_iterator postitr = preitr; postitr++;
00568 for(; postitr != polyline.end(); preitr++, postitr++) {
00569 G4ThreeVector prePts(preitr->x(), preitr->y(), preitr->z());
00570 prePts -= trans;
00571 prePts.transform(rot);
00572 G4ThreeVector postPts(postitr->x(), postitr->y(), postitr->z());
00573 postPts -= trans;
00574 postPts.transform(rot);
00575 G4float * stepPts = new G4float[6];
00576 stepPts[0] = prePts.x();
00577 stepPts[1] = prePts.y();
00578 stepPts[2] = prePts.z();
00579 stepPts[3] = postPts.x();
00580 stepPts[4] = postPts.y();
00581 stepPts[5] = postPts.z();
00582 trajectory.push_back(stepPts);
00583
00584 if(GFDEBUG_TRK) {
00585 G4cout << " ("
00586 << stepPts[0] << ", "
00587 << stepPts[1] << ", "
00588 << stepPts[2] << ") - ("
00589 << stepPts[3] << ", "
00590 << stepPts[4] << ", "
00591 << stepPts[5] << ")" << G4endl;
00592 }
00593 }
00594
00595 const G4VisAttributes * att = polyline.GetVisAttributes();
00596 G4Color color = att->GetColor();
00597 unsigned char trkcolor[3];
00598 trkcolor[0] = (unsigned char)(color.GetRed()*255);
00599 trkcolor[1] = (unsigned char)(color.GetGreen()*255);
00600 trkcolor[2] = (unsigned char)(color.GetBlue()*255);
00601 if(GFDEBUG_TRK) {
00602 G4cout << " color : ["
00603 << color.GetRed() << ", "
00604 << color.GetGreen() << ", "
00605 << color.GetBlue() << "]" << G4endl;
00606 }
00607
00608 kgMocrenIO->addTrack(trajectory, trkcolor);
00609
00610 numTrajectories++;
00611 }
00612
00613 }
00614
00615
00616
00617 void G4GMocrenFileSceneHandler::AddPrimitive (const G4NURBS&)
00618 {
00619
00620 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00621 G4cout << "***** AddPrimitive( G4NURBS )" << G4endl;
00622
00623
00624 GFBeginModeling();
00625
00626 }
00627
00628
00629
00630
00631 void G4GMocrenFileSceneHandler::AddPrimitive ( const G4Text& text )
00632 {
00633 if (fProcessing2D) {
00634 static G4bool warned = false;
00635 if (!warned) {
00636 warned = true;
00637 G4Exception
00638 ("G4GMocrenFileSceneHandler::AddPrimitive (const G4Text&)",
00639 "gMocren1002", JustWarning,
00640 "2D text not implemented. Ignored.");
00641 }
00642 return;
00643 }
00644
00645
00646 G4Text dummytext = text;
00647
00648
00649 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00650 G4cout << "***** AddPrimitive( G4Text )" << G4endl;
00651
00652
00653 GFBeginModeling();
00654
00655 }
00656
00657
00658
00659 void G4GMocrenFileSceneHandler::AddPrimitive ( const G4Circle& mark_circle )
00660 {
00661
00662 G4Circle dummycircle = mark_circle;
00663
00664 if (fProcessing2D) {
00665 static G4bool warned = false;
00666 if (!warned) {
00667 warned = true;
00668 G4Exception
00669 ("G4GMocrenFileSceneHandler::AddPrimitive (const G4Circle&)",
00670 "gMocren1003", JustWarning,
00671 "2D circles not implemented. Ignored.");
00672 }
00673 return;
00674 }
00675
00676
00677 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00678 G4cout << "***** AddPrimitive( G4Circle )" << G4endl;
00679
00680
00681 GFBeginModeling();
00682
00683
00684 }
00685
00686
00687
00688 void G4GMocrenFileSceneHandler::AddPrimitive (const G4Square& mark_square )
00689 {
00690
00691 G4Square dummysquare = mark_square;
00692
00693 if (fProcessing2D) {
00694 static G4bool warned = false;
00695 if (!warned) {
00696 warned = true;
00697 G4Exception
00698 ("G4GMocrenFileSceneHandler::AddPrimitive (const G4Square&)",
00699 "gMocren1004", JustWarning,
00700 "2D squares not implemented. Ignored.");
00701 }
00702 return;
00703 }
00704
00705
00706 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00707 G4cout << "***** AddPrimitive( G4Square )" << G4endl;
00708
00709
00710 GFBeginModeling();
00711
00712 }
00713
00714
00715
00716 void G4GMocrenFileSceneHandler::AddPrimitive ( const G4Polyhedron& polyhedron )
00717 {
00718
00719 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00720 G4cout << "***** AddPrimitive( G4Polyhedron )" << G4endl;
00721
00722
00723 if (polyhedron.GetNoFacets() == 0) return;
00724
00725 if (fProcessing2D) {
00726 static G4bool warned = false;
00727 if (!warned) {
00728 warned = true;
00729 G4Exception
00730 ("G4GMocrenFileSceneHandler::AddPrimitive (const G4Polyhedron&)",
00731 "gMocren1005", JustWarning,
00732 "2D polyhedra not implemented. Ignored.");
00733 }
00734 return;
00735 }
00736
00737
00738 GFBeginModeling();
00739
00740
00741 for (G4int f = polyhedron.GetNoFacets(); f; f--){
00742 G4bool notLastEdge = true;
00743 G4int index = -1;
00744 G4int edgeFlag = 1;
00745
00746
00747 G4int i = 0;
00748 do {
00749
00750 notLastEdge = polyhedron.GetNextVertexIndex(index, edgeFlag);
00751
00752 i++;
00753 }while (notLastEdge);
00754 switch (i){
00755 case 3:
00756
00757 break;
00758 case 4:
00759
00760 break;
00761 default:
00762 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00763 G4cout <<
00764 "ERROR G4GMocrenFileSceneHandler::AddPrimitive(G4Polyhedron)" << G4endl;
00765 G4PhysicalVolumeModel* pPVModel =
00766 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
00767 if (pPVModel)
00768 if(G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00769 G4cout << "Volume " << pPVModel->GetCurrentPV()->GetName() <<
00770 ", Solid " << pPVModel->GetCurrentLV()->GetSolid()->GetName() <<
00771 " (" << pPVModel->GetCurrentLV()->GetSolid()->GetEntityType();
00772
00773 if(G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00774 G4cout <<
00775 "\nG4Polyhedron facet with " << i << " edges" << G4endl;
00776 }
00777 }
00778
00779 }
00780
00781
00782
00783 void G4GMocrenFileSceneHandler::GFEndModeling ()
00784 {
00785 G4VSceneHandler::EndModeling();
00786
00787
00788 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00789 G4cout << "***** GFEndModeling (called)" << G4endl;
00790
00791 if( GFIsInModeling() ) {
00792
00793 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations) {
00794 G4cout << "***** GFEndModeling (started) " ;
00795 G4cout << "(/EndModeling, /DrawAll, /CloseDevice)" << G4endl;
00796 }
00797
00798
00799 EndSavingGdd() ;
00800
00801
00802 kFlagInModeling = false ;
00803
00804 }
00805
00806 }
00807
00808
00809
00810 void G4GMocrenFileSceneHandler::BeginPrimitives (const G4Transform3D& objectTransformation)
00811 {
00812 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00813 G4cout << "***** BeginPrimitives " << G4endl;
00814
00815 GFBeginModeling();
00816
00817 G4VSceneHandler::BeginPrimitives (objectTransformation);
00818
00819 }
00820
00821
00822
00823 void G4GMocrenFileSceneHandler::EndPrimitives ()
00824 {
00825 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00826 G4cout << "***** EndPrimitives " << G4endl;
00827
00828 G4VSceneHandler::EndPrimitives ();
00829 }
00830
00831
00832
00833
00834
00835 void G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )
00836 {
00837 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
00838 G4cout << "***** AddSolid ( box )" << G4endl;
00839
00840 if(GFDEBUG_DET > 0)
00841 G4cout << "G4GMocrenFileSceneHandler::AddSolid(const G4Box&) : "
00842 << box.GetName() << G4endl;
00843
00844
00845 if( !IsVisible() ) { return ; }
00846
00847
00848 GFBeginModeling();
00849
00850
00851
00852 if(GFDEBUG_DET > 1) {
00853 G4cout << "-------" << G4endl;
00854 G4cout << " " << box.GetName() << G4endl;
00855 G4Polyhedron * poly = box.CreatePolyhedron();
00856 poly->Transform(fObjectTransformation);
00857
00858 G4Point3D v1, v2;
00859 G4int next;
00860
00861 for(G4int i = 0; i < 12; i++) {
00862 poly->GetNextEdge(v1, v2, next);
00863 if(next == 0) break;
00864 G4cout << " (" << v1.x() << ", "
00865 << v1.y() << ", "
00866 << v1.z() << ") - ("
00867 << v2.x() << ", "
00868 << v2.y() << ", "
00869 << v2.z() << ") [" << next << "]"
00870 << G4endl;
00871 }
00872 delete poly;
00873 }
00874
00875
00876
00877 G4String volName = kMessenger.getVolumeName();
00878
00879 if(kFlagParameterization != 2) {
00880 G4ScoringManager * pScrMan = G4ScoringManager::GetScoringManager();
00881 if(pScrMan) {
00882 G4ScoringBox * pScBox = dynamic_cast<G4ScoringBox*>(pScrMan->FindMesh(volName));
00883 G4bool bMesh = false;
00884 if(pScBox != NULL) bMesh = true;
00885 if(bMesh) kFlagParameterization = 2;
00886 if(GFDEBUG_DET > 0) G4cout << " G4ScoringManager::FindMesh() : "
00887 << volName << " - " << bMesh << G4endl;
00888 }
00889 }
00890
00891 const G4VModel* pv_model = GetModel();
00892 if (!pv_model) { return ; }
00893 G4PhysicalVolumeModel* pPVModel =
00894 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
00895 if (!pPVModel) { return ; }
00896
00897
00898
00899 if(GFDEBUG_DET > 0) {
00900 G4Material * mat = pPVModel->GetCurrentMaterial();
00901 G4String name = mat->GetName();
00902 G4double dens = mat->GetDensity()/(g/cm3);
00903 G4int copyNo = pPVModel->GetCurrentPV()->GetCopyNo();
00904 G4int depth = pPVModel->GetCurrentDepth();
00905 G4cout << " copy no.: " << copyNo << G4endl;
00906 G4cout << " depth : " << depth << G4endl;
00907 G4cout << " density : " << dens << " [g/cm3]" << G4endl;
00908 G4cout << " location: " << pPVModel->GetCurrentPV()->GetObjectTranslation() << G4endl;
00909 G4cout << " Multiplicity : " << pPVModel->GetCurrentPV()->GetMultiplicity() << G4endl;
00910 G4cout << " Is replicated? : " << pPVModel->GetCurrentPV()->IsReplicated() << G4endl;
00911 G4cout << " Is parameterised? : " << pPVModel->GetCurrentPV()->IsParameterised() << G4endl;
00912 G4cout << " top phys. vol. name : " << pPVModel->GetTopPhysicalVolume()->GetName() << G4endl;
00913 }
00914
00915
00916 if(box.GetName() == volName) {
00917
00918 kVolumeTrans3D = fObjectTransformation;
00919
00920 G4ThreeVector raxis(1., 0., 0.), dummy(0.,0.,0.);
00921 G4RotationMatrix rot(raxis, pi*rad);
00922 G4Transform3D trot(rot, dummy);
00923 if(GFDEBUG_DET) {
00924 G4ThreeVector trans1 = kVolumeTrans3D.getTranslation();
00925 G4RotationMatrix rot1 = kVolumeTrans3D.getRotation().inverse();
00926 G4cout << "kVolumeTrans3D: " << trans1 << G4endl << rot1 << G4endl;
00927 }
00928 kVolumeTrans3D = kVolumeTrans3D*trot;
00929 if(GFDEBUG_DET) G4cout << " Parameterised volume : " << box.GetName() << G4endl;
00930
00931
00932
00933
00934 G4VPhysicalVolume * pv[3] = {0,0,0};
00935 pv[0] = pPVModel->GetCurrentPV()->GetLogicalVolume()->GetDaughter(0);
00936 if(!pv[0]) {
00937 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
00938 "gMocren0003", FatalException, "Unexpected volume.");
00939 }
00940 G4int dirAxis[3] = {-1,-1,-1};
00941 G4int nDaughters[3] = {0,0,0};
00942
00943 EAxis axis; G4int nReplicas; G4double width; G4double offset; G4bool consuming;
00944 pv[0]->GetReplicationData(axis, nReplicas, width, offset, consuming);
00945 nDaughters[0] = nReplicas;
00946 switch(axis) {
00947 case kXAxis: dirAxis[0] = 0; break;
00948 case kYAxis: dirAxis[0] = 1; break;
00949 case kZAxis: dirAxis[0] = 2; break;
00950 default:
00951 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
00952 "gMocren0004", FatalException, "Error.");
00953 }
00954 kNestedVolumeNames.push_back(pv[0]->GetName());
00955 if(GFDEBUG_DET)
00956 G4cout << " daughter name : " << pv[0]->GetName()
00957 << " # : " << nDaughters[0] << G4endl;
00958
00959
00960 if(GFDEBUG_DET) {
00961 if(pv[0]->GetLogicalVolume()->GetNoDaughters()) {
00962 G4cout << "# of daughters : "
00963 << pv[0]->GetLogicalVolume()->GetNoDaughters() << G4endl;
00964 } else {
00965
00966
00967 }
00968 }
00969
00970
00971 if(GFDEBUG_DET) G4cout << "# of daughters : "
00972 << pv[0]->GetLogicalVolume()->GetNoDaughters() << G4endl;
00973 if(pv[0]->GetLogicalVolume()->GetNoDaughters() == 0) {
00974 kFlagParameterization = 1;
00975
00976
00977 }
00978
00979 if(kFlagParameterization == 0) {
00980
00981 pv[1] = pv[0]->GetLogicalVolume()->GetDaughter(0);
00982 if(pv[1]) {
00983 pv[1]->GetReplicationData(axis, nReplicas, width, offset, consuming);
00984 nDaughters[1] = nReplicas;
00985 switch(axis) {
00986 case kXAxis: dirAxis[1] = 0; break;
00987 case kYAxis: dirAxis[1] = 1; break;
00988 case kZAxis: dirAxis[1] = 2; break;
00989 default:
00990 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
00991 "gMocren0007", FatalException, "Error.");
00992 }
00993 kNestedVolumeNames.push_back(pv[1]->GetName());
00994 if(GFDEBUG_DET)
00995 G4cout << " sub-daughter name : " << pv[1]->GetName()
00996 << " # : " << nDaughters[1]<< G4endl;
00997
00998
00999 pv[2] = pv[1]->GetLogicalVolume()->GetDaughter(0);
01000 if(pv[2]) {
01001 nDaughters[2] = pv[2]->GetMultiplicity();
01002 kNestedVolumeNames.push_back(pv[2]->GetName());
01003 if(GFDEBUG_DET)
01004 G4cout << " sub-sub-daughter name : " << pv[2]->GetName()
01005 << " # : " << nDaughters[2] << G4endl;
01006
01007 if(nDaughters[2] > 1) {
01008 G4VNestedParameterisation * nestPara
01009 = dynamic_cast<G4VNestedParameterisation*>(pv[2]->GetParameterisation());
01010 if(nestPara == NULL)
01011 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
01012 "gMocren0008", FatalException, "Non-nested parameterisation");
01013
01014 nestPara->ComputeTransformation(0, pv[2]);
01015 G4ThreeVector trans0 = pv[2]->GetObjectTranslation();
01016 nestPara->ComputeTransformation(1, pv[2]);
01017 G4ThreeVector trans1 = pv[2]->GetObjectTranslation();
01018 G4ThreeVector diff(trans0 - trans1);
01019 if(GFDEBUG_DET)
01020 G4cout << trans0 << " - " << trans1 << " - " << diff << G4endl;
01021
01022 if(diff.x() != 0.) dirAxis[2] = 0;
01023 else if(diff.y() != 0.) dirAxis[2] = 1;
01024 else if(diff.z() != 0.) dirAxis[2] = 2;
01025 else
01026 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
01027 "gMocren0009", FatalException, "Unexpected nested parameterisation");
01028 }
01029 }
01030 }
01031
01032 for(G4int i = 0; i < 3; i++) {
01033 kNestedVolumeDimension[i] = nDaughters[i];
01034
01035 kNestedVolumeDirAxis[i] = dirAxis[i];
01036 }
01037
01038
01039
01040
01041 G4VNestedParameterisation * nestPara
01042 = dynamic_cast<G4VNestedParameterisation*>(pv[2]->GetParameterisation());
01043 if(nestPara != NULL) {
01044 G4double prexyz[3] = {0.,0.,0.}, xyz[3] = {0.,0.,0.};
01045 for(G4int n0 = 0; n0 < nDaughters[0]; n0++) {
01046 for(G4int n1 = 0; n1 < nDaughters[1]; n1++) {
01047 for(G4int n2 = 0; n2 < nDaughters[2]; n2++) {
01048
01049 G4GMocrenTouchable * touch = new G4GMocrenTouchable(n1, n0);
01050 if(GFDEBUG_DET)
01051 G4cout << " retrieve volume : copy # : " << n0
01052 << ", " << n1 << ", " << n2 << G4endl;
01053 G4Material * mat = nestPara->ComputeMaterial(pv[2], n2, touch);
01054 delete touch;
01055 G4double dens = mat->GetDensity()/(g/cm3);
01056
01057 if(GFDEBUG_DET)
01058 G4cout << " density :" << dens << " [g/cm3]" << G4endl;
01059
01060 G4Box tbox(box);
01061 nestPara->ComputeDimensions(tbox, n2, pv[2]);
01062 xyz[0] = tbox.GetXHalfLength()/mm;
01063 xyz[1] = tbox.GetYHalfLength()/mm;
01064 xyz[2] = tbox.GetZHalfLength()/mm;
01065 if(n0 != 0 || n1 != 0 || n2 != 0) {
01066 for(G4int i = 0; i < 3; i++) {
01067 if(xyz[i] != prexyz[i])
01068 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
01069 "gMocren0010", FatalException, "Unsupported parameterisation");
01070 }
01071 }
01072 if(GFDEBUG_DET)
01073 G4cout << " size : " << tbox.GetXHalfLength()/mm << " x "
01074 << tbox.GetYHalfLength()/mm << " x "
01075 << tbox.GetZHalfLength()/mm << " [mm3]" << G4endl;
01076
01077 G4int idx[3];
01078 idx[dirAxis[0]] = n0;
01079 idx[dirAxis[1]] = n1;
01080 idx[dirAxis[2]] = n2;
01081 Index3D i3d(idx[0],idx[1],idx[2]);
01082 kNestedModality[i3d] = dens;
01083 if(GFDEBUG_DET)
01084 G4cout << " index: " << idx[0] << ", " << idx[1] << ", " << idx[2]
01085 << " density: " << dens << G4endl;
01086
01087 for(G4int i = 0; i < 3; i++) prexyz[i] = xyz[i];
01088 }
01089 }
01090 }
01091
01092 kVolumeSize.set(box.GetXHalfLength()*2/mm,
01093 box.GetYHalfLength()*2/mm,
01094 box.GetZHalfLength()*2/mm);
01095
01096 if(!kbSetModalityVoxelSize) {
01097 G4float spacing[3] = {static_cast<G4float>(2*xyz[0]),
01098 static_cast<G4float>(2*xyz[1]),
01099 static_cast<G4float>(2*xyz[2])};
01100 kgMocrenIO->setVoxelSpacing(spacing);
01101 kVoxelDimension.set(spacing[0], spacing[1], spacing[2]);
01102 kbSetModalityVoxelSize = true;
01103 }
01104
01105 } else {
01106 if(GFDEBUG_DET)
01107 G4cout << pv[2]->GetName() << G4endl;
01108 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
01109 "gMocren0011", FatalException, "Non-nested parameterisation");
01110 }
01111
01112
01113
01114
01115 if(GFDEBUG_DET > 1) {
01116 if(pPVModel->GetCurrentPV()->IsParameterised()) {
01117 G4VPVParameterisation * para = pPVModel->GetCurrentPV()->GetParameterisation();
01118 G4cout << " Is nested parameterisation? : " << para->IsNested() << G4endl;
01119
01120
01121 G4int npvp = pPVModel->GetDrawnPVPath().size();
01122 G4cout << " physical volume node id : "
01123 << "size: " << npvp << ", PV name: ";
01124 for(G4int i = 0; i < npvp; i++) {
01125 G4cout << pPVModel->GetDrawnPVPath()[i].GetPhysicalVolume()->GetName()
01126 << " [param:"
01127 << pPVModel->GetDrawnPVPath()[i].GetPhysicalVolume()->IsParameterised()
01128 << ",rep:"
01129 << pPVModel->GetDrawnPVPath()[i].GetPhysicalVolume()->IsReplicated();
01130 if(pPVModel->GetDrawnPVPath()[i].GetPhysicalVolume()->GetParameterisation()) {
01131 G4cout << ",nest:"
01132 << pPVModel->GetDrawnPVPath()[i].GetPhysicalVolume()->GetParameterisation()->IsNested();
01133 }
01134 G4cout << ",copyno:"
01135 << pPVModel->GetDrawnPVPath()[i].GetPhysicalVolume()->GetCopyNo();
01136 G4cout << "] - ";
01137 }
01138 G4cout << G4endl;
01139
01140
01141 pPVModel->GetCurrentPV()->GetReplicationData(axis, nReplicas, width, offset, consuming);
01142 G4cout << " # replicas : " << nReplicas << G4endl;
01143 G4double pareDims[3] = {0.,0.,0.};
01144 G4Box * pbox = dynamic_cast<G4Box *>(pPVModel->GetDrawnPVPath()[npvp-2].GetPhysicalVolume()->GetLogicalVolume()->GetSolid());
01145 if(pbox) {
01146 pareDims[0] = 2.*pbox->GetXHalfLength()*mm;
01147 pareDims[1] = 2.*pbox->GetYHalfLength()*mm;
01148 pareDims[2] = 2.*pbox->GetZHalfLength()*mm;
01149 G4cout << " mother size ["
01150 << pPVModel->GetDrawnPVPath()[npvp-2].GetPhysicalVolume()->GetName()
01151 << "] : "
01152 << pareDims[0] << " x "
01153 << pareDims[1] << " x "
01154 << pareDims[2] << " [mm3]"
01155 << G4endl;
01156 }
01157 G4double paraDims[3];
01158 G4Box * boxP = dynamic_cast<G4Box *>(pPVModel->GetDrawnPVPath()[npvp-1].GetPhysicalVolume()->GetLogicalVolume()->GetSolid());
01159 if(boxP) {
01160 paraDims[0] = 2.*boxP->GetXHalfLength()*mm;
01161 paraDims[1] = 2.*boxP->GetYHalfLength()*mm;
01162 paraDims[2] = 2.*boxP->GetZHalfLength()*mm;
01163 G4cout << " parameterised volume? ["
01164 << pPVModel->GetDrawnPVPath()[npvp-1].GetPhysicalVolume()->GetName()
01165 << "] : "
01166 << paraDims[0] << " x "
01167 << paraDims[1] << " x "
01168 << paraDims[2] << " [mm3] : "
01169 << G4int(pareDims[0]/paraDims[0]) << " x "
01170 << G4int(pareDims[1]/paraDims[1]) << " x "
01171 << G4int(pareDims[2]/paraDims[2]) << G4endl;
01172 } else {
01173 G4cout << pPVModel->GetDrawnPVPath()[npvp-2].GetPhysicalVolume()->GetName()
01174 << " isn't a G4Box." << G4endl;
01175 }
01176 }
01177 }
01178
01179
01180 } else if(kFlagParameterization == 1) {
01181
01182
01183 G4PhantomParameterisation * phantomPara
01184 = dynamic_cast<G4PhantomParameterisation*>(pv[0]->GetParameterisation());
01185 if(phantomPara == NULL) {
01186 G4Exception("G4GMocrenFileSceneHandler::AddSolid( const G4Box& box )",
01187 "gMocren0012", FatalException, "no G4PhantomParameterisation");
01188 } else {
01189 ;
01190 }
01191
01192 kNestedVolumeDimension[0] = phantomPara->GetNoVoxelX();
01193 kNestedVolumeDimension[1] = phantomPara->GetNoVoxelY();
01194 kNestedVolumeDimension[2] = phantomPara->GetNoVoxelZ();
01195 kNestedVolumeDirAxis[0] = 0;
01196 kNestedVolumeDirAxis[1] = 1;
01197 kNestedVolumeDirAxis[2] = 2;
01198
01199
01200 G4int nX = kNestedVolumeDimension[0];
01201 G4int nXY = kNestedVolumeDimension[0]*kNestedVolumeDimension[1];
01202
01203 for(G4int n0 = 0; n0 < kNestedVolumeDimension[0]; n0++) {
01204 for(G4int n1 = 0; n1 < kNestedVolumeDimension[1]; n1++) {
01205 for(G4int n2 = 0; n2 < kNestedVolumeDimension[2]; n2++) {
01206
01207 G4int repNo = n0 + n1*nX + n2*nXY;
01208 G4Material * mat = phantomPara->ComputeMaterial(repNo, pv[0]);
01209 G4double dens = mat->GetDensity()/(g/cm3);
01210
01211
01212 G4int idx[3];
01213 idx[kNestedVolumeDirAxis[0]] = n0;
01214 idx[kNestedVolumeDirAxis[1]] = n1;
01215 idx[kNestedVolumeDirAxis[2]] = n2;
01216 Index3D i3d(idx[0],idx[1],idx[2]);
01217 kNestedModality[i3d] = dens;
01218
01219 if(GFDEBUG_DET)
01220 G4cout << " index: " << idx[0] << ", " << idx[1] << ", " << idx[2]
01221 << " density: " << dens << G4endl;
01222
01223 }
01224 }
01225 }
01226
01227 kVolumeSize.set(box.GetXHalfLength()*2/mm,
01228 box.GetYHalfLength()*2/mm,
01229 box.GetZHalfLength()*2/mm);
01230
01231
01232 if(!kbSetModalityVoxelSize) {
01233 G4float spacing[3] = {static_cast<G4float>(2*phantomPara->GetVoxelHalfX()),
01234 static_cast<G4float>(2*phantomPara->GetVoxelHalfY()),
01235 static_cast<G4float>(2*phantomPara->GetVoxelHalfZ())};
01236 kgMocrenIO->setVoxelSpacing(spacing);
01237 kVoxelDimension.set(spacing[0], spacing[1], spacing[2]);
01238 kbSetModalityVoxelSize = true;
01239 }
01240 }
01241
01242 }
01243
01244
01245
01246 if(!kFLagProcessedInteractiveScorer) {
01247
01248
01249
01250 G4ScoringManager * pScrMan = G4ScoringManager::GetScoringManager();
01251 if(!pScrMan) return;
01252 G4ScoringBox * scoringBox
01253 = dynamic_cast<G4ScoringBox*>(pScrMan->FindMesh(volName));
01254 if(scoringBox == NULL) return;
01255
01256
01257 G4int nVoxels[3];
01258 scoringBox->GetNumberOfSegments(nVoxels);
01259
01260 kNestedVolumeDimension[0] = nVoxels[2];
01261 kNestedVolumeDimension[1] = nVoxels[1];
01262 kNestedVolumeDimension[2] = nVoxels[0];
01263 kNestedVolumeDirAxis[0] = 2;
01264 kNestedVolumeDirAxis[1] = 1;
01265 kNestedVolumeDirAxis[2] = 0;
01266
01267
01268 for(G4int n0 = 0; n0 < kNestedVolumeDimension[0]; n0++) {
01269 for(G4int n1 = 0; n1 < kNestedVolumeDimension[1]; n1++) {
01270 for(G4int n2 = 0; n2 < kNestedVolumeDimension[2]; n2++) {
01271
01272 G4double dens = 0.*(g/cm3);
01273
01274 G4int idx[3];
01275 idx[kNestedVolumeDirAxis[0]] = n0;
01276 idx[kNestedVolumeDirAxis[1]] = n1;
01277 idx[kNestedVolumeDirAxis[2]] = n2;
01278 Index3D i3d(idx[0],idx[1],idx[2]);
01279 kNestedModality[i3d] = dens;
01280
01281 }
01282 }
01283 }
01284
01285 G4ThreeVector boxSize = scoringBox->GetSize();
01286 if(GFDEBUG_DET > 1) {
01287 G4cout << "Interactive Scorer : size - "
01288 << boxSize.x()/cm << " x "
01289 << boxSize.y()/cm << " x "
01290 << boxSize.z()/cm << " [cm3]" << G4endl;
01291 G4cout << "Interactive Scorer : # voxels - "
01292 << nVoxels[0] << " x "
01293 << nVoxels[1] << " x "
01294 << nVoxels[2] << G4endl;
01295 }
01296 kVolumeSize.set(boxSize.x()*2,
01297 boxSize.y()*2,
01298 boxSize.z()*2);
01299
01300
01301 if(!kbSetModalityVoxelSize) {
01302 G4float spacing[3] = {static_cast<G4float>(boxSize.x()*2/nVoxels[0]),
01303 static_cast<G4float>(boxSize.y()*2/nVoxels[1]),
01304 static_cast<G4float>(boxSize.z()*2/nVoxels[2])};
01305
01306 kgMocrenIO->setVoxelSpacing(spacing);
01307 kVoxelDimension.set(spacing[0], spacing[1], spacing[2]);
01308 kbSetModalityVoxelSize = true;
01309
01310 }
01311
01312
01313 kVolumeTrans3D = fObjectTransformation;
01314
01315
01316 G4ThreeVector sbth = scoringBox->GetTranslation();
01317 G4Translate3D sbtranslate(sbth);
01318 kVolumeTrans3D = kVolumeTrans3D*sbtranslate;
01319
01320
01321 G4RotationMatrix sbrm;
01322 sbrm = scoringBox->GetRotationMatrix();
01323 if(!sbrm.isIdentity()) {
01324 G4ThreeVector sbdummy(0.,0.,0.);
01325 G4Transform3D sbrotate(sbrm.inverse(), sbdummy);
01326 kVolumeTrans3D = kVolumeTrans3D*sbrotate;
01327 }
01328
01329
01330
01331 G4ThreeVector raxisY(0., 1., 0.), dummyY(0.,0.,0.);
01332 G4RotationMatrix rotY(raxisY, pi*rad);
01333 G4Transform3D trotY(rotY, dummyY);
01334 G4ThreeVector raxisZ(0., 0., 1.), dummyZ(0.,0.,0.);
01335 G4RotationMatrix rotZ(raxisZ, pi*rad);
01336 G4Transform3D trotZ(rotZ, dummyZ);
01337
01338 kVolumeTrans3D = kVolumeTrans3D*trotY*trotZ;
01339
01340
01341
01342 kFLagProcessedInteractiveScorer = true;
01343 }
01344
01345
01346
01347
01348 G4bool bAddDet = true;
01349 if(!kMessenger.getDrawVolumeGrid()) {
01350
01351 if(kFlagParameterization == 0) {
01352
01353 if(volName == box.GetName()) {
01354 bAddDet = false;
01355 }
01356
01357 std::vector<G4String>::iterator itr = kNestedVolumeNames.begin();
01358 for(; itr != kNestedVolumeNames.end(); itr++) {
01359 if(*itr == box.GetName()) {
01360 bAddDet = false;
01361 break;
01362 }
01363 }
01364 } else if(kFlagParameterization == 1) {
01365
01366 if(volName != box.GetName()) {
01367 bAddDet = false;
01368 }
01369
01370 } else if(kFlagParameterization == 2) {
01371 ;
01372 }
01373
01374 }
01375 if(bAddDet) AddDetector(box);
01376
01377
01378 }
01379
01380
01381
01382 void
01383 G4GMocrenFileSceneHandler::AddSolid( const G4Tubs& tubes )
01384 {
01385 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
01386 G4cout << "***** AddSolid ( tubes )" << G4endl;
01387
01388
01389 if( !IsVisible() ) { return ; }
01390
01391
01392 GFBeginModeling();
01393
01394
01395 AddDetector(tubes);
01396
01397
01398
01399 if(GFDEBUG_DET > 0) {
01400 G4cout << "-------" << G4endl;
01401 G4cout << " " << tubes.GetName() << G4endl;
01402 G4Polyhedron * poly = tubes.CreatePolyhedron();
01403 G4int nv = poly->GetNoVertices();
01404 for(G4int i = 0; i < nv; i++) {
01405 G4cout << " (" << poly->GetVertex(i).x() << ", "
01406 << poly->GetVertex(i).y() << ", "
01407 << poly->GetVertex(i).z() << ")" << G4endl;
01408 }
01409 delete poly;
01410 }
01411
01412 const G4VModel* pv_model = GetModel();
01413 if (!pv_model) { return ; }
01414 G4PhysicalVolumeModel* pPVModel =
01415 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
01416 if (!pPVModel) { return ; }
01417 G4Material * mat = pPVModel->GetCurrentMaterial();
01418 G4String name = mat->GetName();
01419
01420 }
01421
01422
01423
01424
01425 void
01426 G4GMocrenFileSceneHandler::AddSolid( const G4Cons& cons )
01427 {
01428 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
01429 G4cout << "***** AddSolid ( cons )" << G4endl;
01430
01431
01432 if( !IsVisible() ) { return ; }
01433
01434
01435 GFBeginModeling();
01436
01437
01438 AddDetector(cons);
01439
01440 }
01441
01442
01443
01444 void G4GMocrenFileSceneHandler::AddSolid ( const G4Trd& trd )
01445 {
01446 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
01447 G4cout << "***** AddSolid ( trd )" << G4endl;
01448
01449
01450
01451 if( !IsVisible() ) { return ; }
01452
01453
01454 GFBeginModeling();
01455
01456
01457 AddDetector(trd);
01458
01459 }
01460
01461
01462
01463 void G4GMocrenFileSceneHandler::AddSolid ( const G4Sphere& sphere )
01464 {
01465 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
01466 G4cout << "***** AddSolid ( sphere )" << G4endl;
01467
01468
01469 if( !IsVisible() ) { return ; }
01470
01471
01472 GFBeginModeling();
01473
01474
01475 AddDetector(sphere);
01476
01477 }
01478
01479
01480
01481 void G4GMocrenFileSceneHandler::AddSolid (const G4Para& para)
01482 {
01483 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
01484 G4cout << "***** AddSolid ( para )" << G4endl;
01485
01486
01487 if( !IsVisible() ) { return ; }
01488
01489
01490 GFBeginModeling();
01491
01492
01493 AddDetector(para);
01494
01495 }
01496
01497
01498
01499 void G4GMocrenFileSceneHandler::AddSolid (const G4Trap& trap)
01500 {
01501 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
01502 G4cout << "***** AddSolid ( trap )" << G4endl;
01503
01504
01505 if( !IsVisible() ) { return ; }
01506
01507
01508 GFBeginModeling();
01509
01510
01511 AddDetector(trap);
01512
01513 }
01514
01515
01516
01517 void
01518 G4GMocrenFileSceneHandler::AddSolid( const G4Torus& torus )
01519 {
01520 if(GFDEBUG || G4VisManager::GetVerbosity() >= G4VisManager::confirmations)
01521 G4cout << "***** AddSolid ( torus )" << G4endl;
01522
01523
01524 if( !IsVisible() ) { return ; }
01525
01526
01527 GFBeginModeling();
01528
01529
01530 AddDetector(torus);
01531
01532 }
01533
01534
01535
01536
01537 void G4GMocrenFileSceneHandler::AddSolid ( const G4VSolid& solid )
01538 {
01539
01540 if( !IsVisible() ) { return ; }
01541
01542
01543 GFBeginModeling();
01544
01545
01546 AddDetector(solid);
01547
01548
01549 G4VSceneHandler::AddSolid( solid ) ;
01550
01551 }
01552
01553 #include "G4TrajectoriesModel.hh"
01554 #include "G4VTrajectory.hh"
01555 #include "G4VTrajectoryPoint.hh"
01556
01557
01558 void G4GMocrenFileSceneHandler::AddCompound(const G4VTrajectory & traj) {
01559
01560 kbModelingTrajectory = true;
01561
01562 G4VSceneHandler::AddCompound(traj);
01563
01564 if(GFDEBUG_TRK) {
01565 G4cout << " ::AddCompound(const G4VTrajectory&) >>>>>>>>> " << G4endl;
01566 G4TrajectoriesModel * pTrModel = dynamic_cast<G4TrajectoriesModel*>(fpModel);
01567 if (!pTrModel) {
01568 G4Exception
01569 ("G4VSceneHandler::AddCompound(const G4VTrajectory&)",
01570 "gMocren0013", FatalException, "Not a G4TrajectoriesModel.");
01571 } else {
01572 traj.DrawTrajectory(pTrModel->GetDrawingMode());
01573
01574 const G4VTrajectory * trj = pTrModel->GetCurrentTrajectory();
01575 G4cout << "------ track" << G4endl;
01576 G4cout << " name: " << trj->GetParticleName() << G4endl;
01577 G4cout << " id: " << trj->GetTrackID() << G4endl;
01578 G4cout << " charge: " << trj->GetCharge() << G4endl;
01579 G4cout << " momentum: " << trj->GetInitialMomentum() << G4endl;
01580
01581 G4int nPnt = trj->GetPointEntries();
01582 G4cout << " point: ";
01583 for(G4int i = 0; i < nPnt; i++) {
01584 G4cout << trj->GetPoint(i)->GetPosition() << ", ";
01585 }
01586 G4cout << G4endl;
01587 }
01588 G4cout << G4endl;
01589 }
01590
01591 kbModelingTrajectory = false;
01592 }
01593
01594 #include <vector>
01595 #include "G4VHit.hh"
01596 #include "G4AttValue.hh"
01597
01598 void G4GMocrenFileSceneHandler::AddCompound( const G4VHit & hit) {
01599 if(GFDEBUG_HIT) G4cout << " ::AddCompound(const G4VHit&) >>>>>>>>> " << G4endl;
01600
01601 G4VSceneHandler::AddCompound(hit);
01602
01603
01604
01605
01606
01607
01608
01609
01610
01611
01612
01613 std::vector<G4String> hitNames = kMessenger.getHitNames();
01614 if(GFDEBUG_HIT) {
01615 std::vector<G4String>::iterator itr = hitNames.begin();
01616 for(; itr != hitNames.end(); itr++)
01617 G4cout << " hit name : " << *itr << G4endl;
01618 }
01619
01620 std::vector<G4AttValue> * attval = hit.CreateAttValues();
01621 if(!attval) {G4cout << "0 empty " << (unsigned long)attval << G4endl;}
01622 else {
01623
01624 G4bool bid[3] = {false, false, false};
01625 Index3D id;
01626
01627 std::vector<G4AttValue>::iterator itr;
01628
01629 for(itr = attval->begin(); itr != attval->end(); itr++) {
01630 std::string stmp = itr->GetValue();
01631 std::istringstream sval(stmp.c_str());
01632
01633 if(itr->GetName() == G4String("XID")) {
01634 sval >> id.x;
01635 bid[0] = true;
01636 continue;
01637 }
01638 if(itr->GetName() == G4String("YID")) {
01639 sval >> id.y;
01640 bid[1] = true;
01641 continue;
01642 }
01643 if(itr->GetName() == G4String("ZID")) {
01644 sval >> id.z;
01645 bid[2] = true;
01646 continue;
01647 }
01648 }
01649
01650 G4int nhitname = (G4int)hitNames.size();
01651
01652 if(bid[0] && bid[1] && bid[2]) {
01653
01654 if(GFDEBUG_HIT)
01655 G4cout << " Hit : index(" << id.x << ", " << id.y << ", "
01656 << id.z << ")" << G4endl;
01657
01658
01659 for(itr = attval->begin(); itr != attval->end(); itr++) {
01660 for(G4int i = 0; i < nhitname; i++) {
01661 if(itr->GetName() == hitNames[i]) {
01662
01663 std::string stmp = itr->GetValue();
01664 std::istringstream sval(stmp.c_str());
01665 G4double value;
01666 G4String unit;
01667 sval >> value >> unit;
01668
01669 std::map<G4String, std::map<Index3D, G4double> >::iterator kNestedHitsListItr;
01670 kNestedHitsListItr = kNestedHitsList.find(hitNames[i]);
01671 if(kNestedHitsListItr != kNestedHitsList.end()) {
01672
01673
01674 kNestedHitsListItr->second[id] = value;
01675 } else {
01676 std::map<Index3D, G4double> hits;
01677 hits.insert(std::map<Index3D, G4double>::value_type(id, value));
01678 kNestedHitsList[hitNames[i]] = hits;
01679 }
01680
01681
01682 if(GFDEBUG_HIT)
01683 G4cout << " : " << hitNames[i] << " -> " << value
01684 << " [" << unit << "]" << G4endl;
01685 }
01686 }
01687 }
01688 } else {
01689 G4Exception("G4GMocrenFileSceneHandler::AddCompound(const G4VHit &)",
01690 "gMocren0014", FatalException, "Error");
01691 }
01692
01693 delete attval;
01694 }
01695
01696 }
01697
01698 void G4GMocrenFileSceneHandler::AddCompound( const G4VDigi & digi) {
01699 if(GFDEBUG_DIGI) G4cout << " ::AddCompound(const G4VDigi&) >>>>>>>>> " << G4endl;
01700 G4VSceneHandler::AddCompound(digi);
01701 }
01702
01703 void G4GMocrenFileSceneHandler::AddCompound(const G4THitsMap<G4double> & hits) {
01704 if(GFDEBUG_HIT)
01705 G4cout << " ::AddCompound(const std::map<G4int, G4double*> &) >>>>>>>>> " << G4endl;
01706
01707
01708 std::vector<G4String> hitScorerNames = kMessenger.getHitScorerNames();
01709 G4int nhitname = (G4int)hitScorerNames.size();
01710 G4String scorername = static_cast<G4VHitsCollection>(hits).GetName();
01711
01712
01713
01714
01715
01716
01717
01718
01719
01720
01721
01722 {
01723
01724
01725
01726 G4int idx[3];
01727 std::map<G4int, G4double*> * map = hits.GetMap();
01728 std::map<G4int, G4double*>::const_iterator itr = map->begin();
01729 for(; itr != map->end(); itr++) {
01730 GetNestedVolumeIndex(itr->first, idx);
01731 Index3D id(idx[0], idx[1], idx[2]);
01732
01733 std::map<G4String, std::map<Index3D, G4double> >::iterator nestedHitsListItr;
01734 nestedHitsListItr = kNestedHitsList.find(scorername);
01735 if(nestedHitsListItr != kNestedHitsList.end()) {
01736 nestedHitsListItr->second[id] = *itr->second;
01737 } else {
01738 std::map<Index3D, G4double> hit;
01739 hit.insert(std::map<Index3D, G4double>::value_type(id, *itr->second));
01740 kNestedHitsList[scorername] = hit;
01741 }
01742 }
01743
01744
01745
01746
01747 }
01748
01749 if(GFDEBUG_HIT) {
01750 G4String meshname = static_cast<G4VHitsCollection>(hits).GetSDname();
01751 G4cout << " >>>>> " << meshname << " : " << scorername << G4endl;
01752
01753 for(G4int i = 0; i < nhitname; i++)
01754 if(scorername == hitScorerNames[i])
01755 G4cout << " !!!! Hit scorer !!!! " << scorername << G4endl;
01756
01757 G4cout << " dimension: "
01758 << kNestedVolumeDimension[0] << " x "
01759 << kNestedVolumeDimension[1] << " x "
01760 << kNestedVolumeDimension[2] << G4endl;
01761
01762 G4int id[3];
01763 std::map<G4int, G4double*> * map = hits.GetMap();
01764 std::map<G4int, G4double*>::const_iterator itr = map->begin();
01765 for(; itr != map->end(); itr++) {
01766 GetNestedVolumeIndex(itr->first, id);
01767 G4cout << "[" << itr->first << "] "
01768 << "("<< id[0] << "," << id[1] << "," << id[2] << ")"
01769 << *itr->second << ", ";
01770 }
01771 G4cout << G4endl;
01772 }
01773
01774 }
01775
01776
01777 G4bool G4GMocrenFileSceneHandler::IsVisible()
01778 {
01779
01780 G4bool visibility = true ;
01781
01782 const G4VisAttributes* pVisAttribs =
01783 fpViewer->GetApplicableVisAttributes( fpVisAttribs );
01784
01785 if(pVisAttribs) {
01786 visibility = pVisAttribs->IsVisible();
01787 }
01788
01789 return visibility ;
01790
01791 }
01792
01793
01794
01795 void G4GMocrenFileSceneHandler::ClearTransientStore()
01796 {
01797
01798
01799
01800 if (fpViewer) {
01801 fpViewer -> SetView ();
01802 fpViewer -> ClearView ();
01803 fpViewer -> DrawView ();
01804 }
01805 }
01806
01807
01808 void G4GMocrenFileSceneHandler::AddDetector(const G4VSolid & solid) {
01809
01810 Detector detector;
01811
01812
01813 detector.name = solid.GetName();
01814 if(GFDEBUG_DET > 1)
01815 G4cout << "0 Detector name : " << detector.name << G4endl;
01816
01817 const G4VModel* pv_model = GetModel();
01818 if (!pv_model) { return ; }
01819 G4PhysicalVolumeModel* pPVModel =
01820 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
01821 if (!pPVModel) { return ; }
01822
01823
01824 std::vector<G4float *> dedges;
01825 G4Polyhedron * poly = solid.CreatePolyhedron();
01826 detector.polyhedron = poly;
01827 detector.transform3D = fObjectTransformation;
01828
01829
01830 unsigned char uccolor[3] = {30, 30, 30};
01831 if(pPVModel->GetCurrentLV()->GetVisAttributes()) {
01832 G4Color color = pPVModel->GetCurrentLV()->GetVisAttributes()->GetColor();
01833 uccolor[0] = (unsigned char)(color.GetRed()*255);
01834 uccolor[1] = (unsigned char)(color.GetGreen()*255);
01835 uccolor[2] = (unsigned char)(color.GetBlue()*255);
01836
01837
01838 }
01839 for(G4int i = 0; i < 3; i++) detector.color[i] = uccolor[i];
01840
01841 kDetectors.push_back(detector);
01842
01843 if(GFDEBUG_DET > 1) {
01844 G4cout << "0 color: (" << (G4int)uccolor[0] << ", "
01845 << (G4int)uccolor[1] << ", " << (G4int)uccolor[2] << ")"
01846 << G4endl;
01847 }
01848
01849 }
01850
01851
01852 void G4GMocrenFileSceneHandler::ExtractDetector() {
01853
01854 std::vector<Detector>::iterator itr = kDetectors.begin();
01855
01856 for(; itr != kDetectors.end(); itr++) {
01857
01858
01859 G4String detname = itr->name;
01860 if(GFDEBUG_DET > 1)
01861 G4cout << "Detector name : " << detname << G4endl;
01862
01863
01864 std::vector<G4float *> dedges;
01865 G4Polyhedron * poly = itr->polyhedron;
01866 poly->Transform(itr->transform3D);
01867 G4Transform3D invVolTrans = kVolumeTrans3D.inverse();
01868 poly->Transform(invVolTrans);
01869
01870 G4Point3D v1, v2;
01871 G4bool bnext = true;
01872 G4int next;
01873 G4int nedges = 0;
01874
01875 while(bnext) {
01876 if(!(poly->GetNextEdge(v1, v2, next))) bnext =false;
01877 G4float * edge = new G4float[6];
01878 edge[0] = v1.x()/mm;
01879 edge[1] = v1.y()/mm;
01880 edge[2] = v1.z()/mm;
01881 edge[3] = v2.x()/mm;
01882 edge[4] = v2.y()/mm;
01883 edge[5] = v2.z()/mm;
01884 dedges.push_back(edge);
01885 nedges++;
01886 }
01887
01888
01889 unsigned char uccolor[3] = {itr->color[0],
01890 itr->color[1],
01891 itr->color[2]};
01892
01893 kgMocrenIO->addDetector(detname, dedges, uccolor);
01894 for(G4int i = 0; i < nedges; i++) {
01895 delete [] dedges[i];
01896 }
01897 dedges.clear();
01898
01899 if(GFDEBUG_DET > 1) {
01900 G4cout << " color: (" << (G4int)uccolor[0] << ", "
01901 << (G4int)uccolor[1] << ", " << (G4int)uccolor[2] << ")"
01902 << G4endl;
01903 }
01904 }
01905 }
01906
01907 void G4GMocrenFileSceneHandler::GetNestedVolumeIndex(G4int _idx, G4int _idx3d[3]) {
01908 if(kNestedVolumeDimension[0] == 0 ||
01909 kNestedVolumeDimension[1] == 0 ||
01910 kNestedVolumeDimension[2] == 0) {
01911 for(G4int i = 0; i < 3; i++) _idx3d[i] = 0;
01912 return;
01913 }
01914
01915
01916 if(kFlagParameterization == 0) {
01917
01918 G4int plane = kNestedVolumeDimension[2]*kNestedVolumeDimension[1];
01919 G4int line = kNestedVolumeDimension[2];
01920
01921
01922
01923
01924
01925
01926
01927
01928
01929
01930
01931 _idx3d[kNestedVolumeDirAxis[0]] = _idx/plane;
01932 _idx3d[kNestedVolumeDirAxis[1]] = (_idx%plane)/line;
01933 _idx3d[kNestedVolumeDirAxis[2]] = (_idx%plane)%line;
01934
01935
01936
01937
01938
01939
01940
01941
01942
01943
01944
01945
01946
01947
01948
01949
01950
01951
01952
01953
01954
01955
01956
01957
01958 } else {
01959
01960 G4int plane = kNestedVolumeDimension[0]*kNestedVolumeDimension[1];
01961 G4int line = kNestedVolumeDimension[0];
01962 _idx3d[kNestedVolumeDirAxis[2]] = _idx/plane;
01963 _idx3d[kNestedVolumeDirAxis[1]] = (_idx%plane)/line;
01964 _idx3d[kNestedVolumeDirAxis[0]] = (_idx%plane)%line;
01965
01966 }
01967
01968 }
01969
01970
01971
01972 G4GMocrenFileSceneHandler::Detector::Detector()
01973 : polyhedron(0) {
01974 color[0] = color[1] = color[2] = 255;
01975 }
01976 G4GMocrenFileSceneHandler::Detector::~Detector() {
01977 if(!polyhedron) delete polyhedron;
01978 }
01979 void G4GMocrenFileSceneHandler::Detector::clear() {
01980 name.clear();
01981 if(!polyhedron) delete polyhedron;
01982 color[0] = color[1] = color[2] = 255;
01983 transform3D = G4Transform3D::Identity;
01984 }
01985
01986
01987 G4GMocrenFileSceneHandler::Index3D::Index3D()
01988 : x(0), y(0), z(0) {
01989 ;
01990 }
01991
01992 G4GMocrenFileSceneHandler::Index3D::Index3D(const Index3D & _index3D)
01993 : x(_index3D.x), y(_index3D.y), z(_index3D.z) {
01994
01995
01996
01997
01998
01999
02000 ;
02001 }
02002
02003 G4GMocrenFileSceneHandler::Index3D::Index3D(G4int _x, G4int _y, G4int _z)
02004 : x(_x), y(_y), z(_z) {
02005 ;
02006 }
02007 G4bool G4GMocrenFileSceneHandler::Index3D::operator < (const Index3D & _right) const {
02008 if(z < static_cast<Index3D>(_right).z) {
02009 return true;
02010 } else if(z == _right.z) {
02011 if(y < static_cast<Index3D>(_right).y) return true;
02012 else if(y == _right.y)
02013 if(x < static_cast<Index3D>(_right).x) return true;
02014 }
02015 return false;
02016 }
02017 G4bool G4GMocrenFileSceneHandler::Index3D::operator == (const Index3D & _right) const {
02018 if(z == _right.z && y == _right.y && x == _right.x) return true;
02019 return false;
02020 }