CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/SimG4Core/Geometry/src/DDG4SolidConverter.cc

Go to the documentation of this file.
00001 #include "SimG4Core/Geometry/interface/DDG4SolidConverter.h"
00002 #include "G4VSolid.hh"
00003 
00004 #include "DetectorDescription/Core/interface/DDSolid.h"
00005 
00006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00007 #include "FWCore/Utilities/interface/Exception.h"
00008 
00009 using namespace std;
00010 const vector<double> * DDG4SolidConverter::par_ = 0; 
00011 
00012 DDG4SolidConverter::DDG4SolidConverter()
00013 {
00014    // could also be done 'dynamically' from outside 
00015    // would then need to have a 'register' method ...
00016   par_=0;
00017   convDispatch_[ddbox]            = DDG4SolidConverter::box; 
00018   convDispatch_[ddtubs]           = DDG4SolidConverter::tubs;
00019   convDispatch_[ddtrap]           = DDG4SolidConverter::trap;
00020   convDispatch_[ddcons]           = DDG4SolidConverter::cons;
00021   convDispatch_[ddpolycone_rrz]   = DDG4SolidConverter::polycone_rrz;
00022   convDispatch_[ddpolycone_rz]    = DDG4SolidConverter::polycone_rz;
00023   convDispatch_[ddpolyhedra_rrz]  = DDG4SolidConverter::polyhedra_rrz;
00024   convDispatch_[ddpolyhedra_rz]   = DDG4SolidConverter::polyhedra_rz;   
00025   convDispatch_[ddtorus]          = DDG4SolidConverter::torus;   
00026   convDispatch_[ddreflected]      = DDG4SolidConverter::reflected;
00027   convDispatch_[ddunion]          = DDG4SolidConverter::unionsolid;
00028   convDispatch_[ddintersection]   = DDG4SolidConverter::intersection;
00029   convDispatch_[ddsubtraction]    = DDG4SolidConverter::subtraction;
00030   convDispatch_[ddpseudotrap]     = DDG4SolidConverter::pseudotrap;
00031   convDispatch_[ddtrunctubs]      = DDG4SolidConverter::trunctubs;
00032   convDispatch_[ddsphere]         = DDG4SolidConverter::sphere;   
00033   convDispatch_[ddorb]            = DDG4SolidConverter::orb;   
00034   convDispatch_[ddellipticaltube] = DDG4SolidConverter::ellipticaltube;   
00035   convDispatch_[ddellipsoid]      = DDG4SolidConverter::ellipsoid;   
00036   convDispatch_[ddparallelepiped] = DDG4SolidConverter::para;   
00037 }
00038 
00039 
00040 DDG4SolidConverter::~DDG4SolidConverter() { }
00041 
00042 G4VSolid * DDG4SolidConverter::convert(const DDSolid & s)
00043 {
00044   if ( !s ) {
00045     edm::LogError("SimG4CoreGeometry") <<" DDG4SolidConverter::convert(..) found an undefined DDSolid " << s.toString();
00046     throw cms::Exception("SimG4CoreGeometry", "DDG4SolidConverter::convert(..) found an undefined DDSolid " + s.toString());
00047   }
00048    G4VSolid * result = 0;
00049    par_ = &(s.parameters());
00050    map<DDSolidShape,FNPTR>::iterator it = convDispatch_.find(s.shape());
00051    if (it != convDispatch_.end()) {
00052      result = it->second(s);
00053    } 
00054    else {
00055      throw cms::Exception("DetectorDescriptionFault") 
00056        <<  "DDG4SolidConverter::convert: conversion failed for s=" << s
00057        << "\n solid.shape()=" << s.shape()
00058        << std::endl;
00059    }
00060    return result;
00061 }
00062 
00063 
00064 #include "G4Box.hh"
00065 G4VSolid * DDG4SolidConverter::box(const DDSolid & s)
00066 {
00067    LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: box = " << s ;
00068    return new G4Box(s.name().name(), (*par_)[0],(*par_)[1],(*par_)[2]);
00069 }
00070 
00071 
00072 #include "G4Tubs.hh"
00073 G4VSolid * DDG4SolidConverter::tubs(const DDSolid & s)
00074 {
00075    LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: tubs = " << s ;  
00076    return new G4Tubs(s.name().name(), (*par_)[1], // rmin
00077                                (*par_)[2], // rmax
00078                                (*par_)[0], // dzHalf
00079                                (*par_)[3], // phiStart
00080                                (*par_)[4]);// deltaPhi
00081 }
00082 
00083 
00084 #include "G4Trap.hh"
00085 G4VSolid * DDG4SolidConverter::trap(const DDSolid & s)
00086 {
00087    LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: trap = " << s ;
00088    return new G4Trap(s.name().name(), (*par_)[0],  // pDz
00089                                (*par_)[1],  // theta
00090                                (*par_)[2],  // phi
00091                                (*par_)[3],  // y1
00092                                (*par_)[4],  // x1
00093                                (*par_)[5],  // x2
00094                                (*par_)[6],  // alpha1
00095                                (*par_)[7],  // y2
00096                                (*par_)[8],  // x3
00097                                (*par_)[9],  // x4
00098                                (*par_)[10]);// alpha2
00099 }
00100 
00101 
00102 #include "G4Cons.hh"
00103 G4VSolid * DDG4SolidConverter::cons(const DDSolid & s) 
00104 {
00105    LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: cons = " << s ;  
00106    return new G4Cons(s.name().name(), (*par_)[1],  // rmin -z
00107                                (*par_)[2],  // rmax -z
00108                                (*par_)[3],  // rmin +z
00109                                (*par_)[4],  // rmax +z
00110                                (*par_)[0],  // zHalf
00111                                (*par_)[5],  // phistart
00112                                (*par_)[6]); // deltaphi
00113 }
00114 
00115             
00116 #include "G4Polycone.hh"
00117 G4VSolid * DDG4SolidConverter::polycone_rz(const DDSolid & s)
00118 {
00119   LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: pcon_rz = " << s ;  
00120   vector<double> r;
00121   vector<double> z;
00122   vector<double>::const_iterator i = (*par_).begin()+2;
00123   int count=0;
00124   for(; i!=(*par_).end(); ++i) {
00125     LogDebug("SimG4CoreGeometry") << "z=" << *i ;
00126     z.push_back(*i); ++i;
00127     LogDebug("SimG4CoreGeometry") << " r=" << *i ;
00128     r.push_back(*i);
00129     count++;
00130    }
00131    LogDebug("SimG4CoreGeometry") << "sp=" << (*par_)[0]/deg << " ep=" << (*par_)[1]/deg ;
00132    return new G4Polycone(s.name().name(), (*par_)[0], (*par_)[1], // start,delta-phi
00133                          count, // numRZ
00134                          &(r[0]),
00135                          &(z[0]));
00136 }
00137 
00138 
00139 G4VSolid * DDG4SolidConverter::polycone_rrz(const DDSolid & s)
00140 {
00141     LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: pcon_rrz = " << s ;  
00142     vector<double> z_p;
00143     vector<double> rmin_p;
00144     vector<double> rmax_p;
00145     vector<double>::const_iterator i = par_->begin()+2;
00146     int count = 0;
00147     for (; i!=par_->end(); ++i) {
00148         LogDebug("SimG4CoreGeometry") << "z=" << *i ;
00149         z_p.push_back(*i); ++i;
00150         LogDebug("SimG4CoreGeometry") << "rmin=" << *i ;
00151         rmin_p.push_back(*i); ++i;
00152         LogDebug("SimG4CoreGeometry") << "rmax=" << *i ;
00153         rmax_p.push_back(*i); 
00154       count++;
00155     }
00156     LogDebug("SimG4CoreGeometry") << "sp=" << (*par_)[0]/deg << " ep=" << (*par_)[1]/deg ;
00157     return new G4Polycone(s.name().name(), (*par_)[0], (*par_)[1], // start,delta-phi
00158                         count, // sections
00159                         &(z_p[0]),
00160                         &(rmin_p[0]),
00161                         &(rmax_p[0]));
00162                         
00163 }
00164 
00165 
00166 #include "G4Polyhedra.hh"
00167 G4VSolid * DDG4SolidConverter::polyhedra_rz(const DDSolid & s)
00168 {
00169     LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: phed_rz = " << s ;  
00170     vector<double> r;
00171     vector<double> z;
00172     vector<double>::const_iterator i = par_->begin()+3;
00173     int count=0;
00174     
00175     for(; i!=par_->end(); ++i) {
00176       z.push_back(*i); ++i;
00177       r.push_back(*i);
00178       count++;
00179     }
00180       
00181     return new G4Polyhedra(s.name().name(), (*par_)[1], (*par_)[2], int((*par_)[0]),// start,delta-phi;sides
00182                            count, // numRZ
00183                            &(r[0]),
00184                            &(z[0]));
00185 }
00186 
00187 
00188 G4VSolid * DDG4SolidConverter::polyhedra_rrz(const DDSolid & s)
00189 {
00190     LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: phed_rrz = " << s ;  
00191     vector<double> z_p;
00192     vector<double> rmin_p;
00193     vector<double> rmax_p;
00194     vector<double>::const_iterator i = par_->begin()+3;
00195     int count = 0;
00196     for (; i!=par_->end(); ++i) {
00197         LogDebug("SimG4CoreGeometry") << "z=" << *i ;
00198         z_p.push_back(*i); ++i;
00199         LogDebug("SimG4CoreGeometry") << "rmin=" << *i ;
00200         rmin_p.push_back(*i); ++i;
00201         LogDebug("SimG4CoreGeometry") << "rmax=" << *i ;
00202         rmax_p.push_back(*i); 
00203       count++;
00204     }
00205     LogDebug("SimG4CoreGeometry") << "sp=" << (*par_)[0]/deg << " ep=" << (*par_)[1]/deg ;
00206     return new G4Polyhedra(s.name().name(), (*par_)[1], (*par_)[2], int((*par_)[0]), // start,delta-phi,sides
00207                            count, // sections
00208                            &(z_p[0]),
00209                            &(rmin_p[0]),
00210                            &(rmax_p[0]));  
00211 }
00212 
00213 #include "G4Torus.hh"
00214 G4VSolid * DDG4SolidConverter::torus(const DDSolid & s)
00215 {
00216    LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: torus = " << s ;  
00217    return new G4Torus(s.name().name(), (*par_)[0], // rmin
00218                                (*par_)[1], // rmax
00219                                (*par_)[2], // Rtor
00220                                (*par_)[3], // phiStart
00221                                (*par_)[4]);// deltaPhi
00222 }
00223 
00224 
00225 #include "G4ReflectedSolid.hh"
00226 G4VSolid * DDG4SolidConverter::reflected(const DDSolid & s)
00227 {
00228   LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: reflected = " << s ;  
00229   G4ReflectedSolid * rs = 0;
00230   DDReflectionSolid rfs(s); 
00231   if (rfs) {    
00232     static /* G4Transform3D */ HepGeom::ReflectZ3D z_reflection; // = HepGeom::ReflectZ3D;      
00233     rs = new G4ReflectedSolid(s.name().name(), 
00234                               DDG4SolidConverter().convert(rfs.unreflected()), 
00235                               z_reflection);
00236     
00237   } // else ?
00238   return rs;
00239 }
00240 
00241 
00242 #include "G4UnionSolid.hh"
00243 G4VSolid * DDG4SolidConverter::unionsolid(const DDSolid & s)
00244 {
00245   LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: unionsolid = " << s.name() ;
00246   G4UnionSolid * us = 0;
00247   DDBooleanSolid bs(s);
00248   if (bs) {
00249     LogDebug("SimG4CoreGeometry") << "SolidA=" << bs.solidA();
00250     G4VSolid * sa = DDG4SolidConverter().convert(bs.solidA());
00251     LogDebug("SimG4CoreGeometry") << "SolidB=" << bs.solidB();
00252     G4VSolid * sb = DDG4SolidConverter().convert(bs.solidB());
00253     LogDebug("SimG4CoreGeometry") << " name:" << s.name() << " t=" << bs.translation() << flush;
00254     LogDebug("SimG4CoreGeometry") << " " << bs.rotation().rotation()->Inverse() << flush;
00255     std::vector<double> tdbl(9);
00256     bs.rotation().rotation()->Inverse().GetComponents(tdbl.begin(), tdbl.end());
00257     CLHEP::HepRep3x3 temprep(tdbl[0], tdbl[1], tdbl[2], tdbl[3], tdbl[4], tdbl[5], tdbl[6], tdbl[7], tdbl[8]);
00258     CLHEP::Hep3Vector temphvec(bs.translation().X(), bs.translation().Y(), bs.translation().Z()); 
00259     us = new G4UnionSolid(s.name().name(),
00260                           sa,
00261                           sb,
00262                           new CLHEP::HepRotation(temprep),
00263                           temphvec);
00264     
00265   } // else?
00266   return us;       
00267 }
00268 
00269 
00270 #include "G4SubtractionSolid.hh"
00271 #include <sstream>
00272 G4VSolid * DDG4SolidConverter::subtraction(const DDSolid & s)
00273 {
00274   LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: subtraction = " << s ;
00275   G4SubtractionSolid * us = 0;
00276   DDBooleanSolid bs(s);
00277   if (bs) {
00278     G4VSolid * sa = DDG4SolidConverter().convert(bs.solidA());
00279     G4VSolid * sb = DDG4SolidConverter().convert(bs.solidB());
00280     LogDebug("SimG4CoreGeometry") << " name:" << s.name() << " t=" << bs.translation() << flush;
00281     //       stringstream sst;
00282     //       bs.rotation().rotation()->inverse().print(sst);
00283     //       LogDebug("SimG4CoreGeometry") << " " << sst.str() << flush;
00284     LogDebug("SimG4CoreGeometry") << " " << bs.rotation().rotation()->Inverse() << flush;
00285     std::vector<double> tdbl(9);
00286     bs.rotation().rotation()->Inverse().GetComponents(tdbl.begin(), tdbl.end());
00287     CLHEP::HepRep3x3 temprep(tdbl[0], tdbl[1], tdbl[2], tdbl[3], tdbl[4], tdbl[5], tdbl[6], tdbl[7], tdbl[8]);
00288     CLHEP::Hep3Vector temphvec(bs.translation().X(), bs.translation().Y(), bs.translation().Z()); 
00289     us = new G4SubtractionSolid(s.name().name(),
00290                                 sa,
00291                                 sb,
00292                                 new CLHEP::HepRotation(temprep),
00293                                 temphvec);
00294   }
00295   return us;       
00296 }
00297 
00298 
00299 #include "G4IntersectionSolid.hh"
00300 G4VSolid * DDG4SolidConverter::intersection(const DDSolid & s)
00301 {
00302   LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: intersection = " << s ;
00303   G4IntersectionSolid * us = 0;
00304   DDBooleanSolid bs(s);
00305   if (bs) {
00306     G4VSolid * sa = DDG4SolidConverter().convert(bs.solidA());
00307     G4VSolid * sb = DDG4SolidConverter().convert(bs.solidB());
00308     LogDebug("SimG4CoreGeometry") << " name:" << s.name() << " t=" << bs.translation() << flush;
00309     LogDebug("SimG4CoreGeometry") << " " << bs.rotation().rotation()->Inverse() << flush;
00310     std::vector<double> tdbl(9);
00311     bs.rotation().rotation()->Inverse().GetComponents(tdbl.begin(), tdbl.end());
00312     CLHEP::HepRep3x3 temprep(tdbl[0], tdbl[1], tdbl[2], tdbl[3], tdbl[4], tdbl[5], tdbl[6], tdbl[7], tdbl[8]);
00313     CLHEP::Hep3Vector temphvec(bs.translation().X(), bs.translation().Y(), bs.translation().Z()); 
00314     us = new G4IntersectionSolid(s.name().name(),
00315                                  sa,
00316                                  sb,
00317                                  new CLHEP::HepRotation(temprep),
00318                                  temphvec);
00319   }
00320   return us;       
00321 }
00322 
00323 
00324 #include "G4Trd.hh"
00325 G4VSolid * DDG4SolidConverter::pseudotrap(const DDSolid & s)
00326 {
00327   static G4RotationMatrix * rot = 0;
00328   static bool firstTime=true;
00329   if (firstTime) {
00330     firstTime=false;
00331     rot = new G4RotationMatrix;
00332     rot->rotateX(90.*deg);
00333     
00334   }
00335   LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: pseudoTrap = " << s ;
00336   G4Trd * trap = 0;
00337   G4Tubs * tubs = 0;
00338   G4VSolid * result = 0;
00339   DDPseudoTrap pt(s); // pt...PseudoTrap
00340   double r = pt.radius();
00341   bool atMinusZ = pt.atMinusZ();
00342   double x = 0;
00343   double h = 0;
00344   bool intersec = false; // union or intersection solid
00345   if (pt.atMinusZ()) {
00346     x = pt.x1(); // tubs radius
00347   } 
00348   else {
00349     x = pt.x2(); // tubs radius
00350   }
00351   double openingAngle = 2.*asin(x/abs(r));
00352   //trap = new G4Trd(s.name().name(), 
00353   double displacement=0;
00354   double startPhi=0;
00355   /* calculate the displacement of the tubs w.r.t. to the trap,
00356      determine the opening angle of the tubs */
00357   double delta = sqrt(r*r-x*x);
00358   if (r < 0 && abs(r) >= x) {
00359     intersec = true; // intersection solid
00360     h = pt.y1() < pt.y2() ? pt.y2() : pt.y1(); // tubs half height
00361     h += h/20.; // enlarge a bit - for subtraction solid
00362     if (atMinusZ) {
00363       displacement = - pt.halfZ() - delta; 
00364       startPhi = 270.*deg - openingAngle/2.;
00365     }
00366     else {
00367       displacement =   pt.halfZ() + delta;
00368       startPhi = 90.*deg - openingAngle/2.;
00369     }
00370   }
00371   else if ( r > 0 && abs(r) >= x )
00372     {
00373       if (atMinusZ) {
00374         displacement = - pt.halfZ() + delta;
00375         startPhi = 90.*deg - openingAngle/2.;
00376         h = pt.y1();
00377       }
00378       else {
00379         displacement =   pt.halfZ() - delta; 
00380         startPhi = 270.*deg - openingAngle/2.;
00381         h = pt.y2();
00382       }    
00383     }
00384   else {
00385     throw cms::Exception("DetectorDescriptionFault", "Check parameters of the PseudoTrap! name=" + pt.name().name());
00386    }
00387   G4ThreeVector displ(0.,0.,displacement); // displacement of the tubs w.r.t. trap
00388   LogDebug("SimG4CoreGeometry") << "DDSolidConverter::pseudotrap(): displacement=" << displacement 
00389                                 << " openingAngle=" << openingAngle/deg << " x=" << x << " h=" << h;
00390     
00391   // Now create two solids (trd & tubs), and a boolean solid out of them 
00392   string name=pt.name().name();
00393   trap = new G4Trd(name, pt.x1(), pt.x2(), pt.y1(), pt.y2(), pt.halfZ());
00394   tubs = new G4Tubs(name, 
00395                     0., // rMin
00396                     abs(r), // rMax
00397                     h, // half height
00398                     startPhi, // start angle
00399                     openingAngle);
00400   if (intersec) {
00401     result = new G4SubtractionSolid(name, trap, tubs, rot, displ);
00402   }
00403   else {
00405     G4VSolid * tubicCap = new G4SubtractionSolid(name, 
00406                                                  tubs, 
00407                                                  new G4Box(name, 1.1*x, sqrt(r*r-x*x), 1.1*h),  
00408                                                  0, 
00409                                                  G4ThreeVector());
00410     result = new G4UnionSolid(name, trap, tubicCap, rot, displ);
00411             
00412     // approximative implementation - also fails to visualize due to G4/Iguana limitations
00413     /*
00414       delete tubs;
00415       tubs = new G4Tubs(name, 
00416       sqrt(r*r-x*x), // rMin-approximation!
00417       abs(r), // rMax
00418       h, // half height
00419       startPhi, // start angle
00420       openingAngle);
00421       result = new G4UnionSolid(name, trap, tubs, rot, displ);
00422     */
00423   }                                
00424   return result;
00425 }
00426 
00427 
00428 G4VSolid * DDG4SolidConverter::trunctubs(const DDSolid & s)
00429 {
00430   // truncated tube-section: a boolean subtraction solid:
00431   //                         from a tube-section a box is subtracted according to the  
00432   //                         given parameters
00433   LogDebug("SimG4CoreGeometry") << "MantisConverter: solidshape=" << DDSolidShapesName::name(s.shape()) << " " << s;
00434   LogDebug("SimG4CoreGeometry") << "before";
00435   DDTruncTubs tt(s);
00436   LogDebug("SimG4CoreGeometry") << "after";
00437   double rIn(tt.rIn()), rOut(tt.rOut()), zHalf(tt.zHalf()),
00438     startPhi(tt.startPhi()), deltaPhi(tt.deltaPhi()), 
00439     cutAtStart(tt.cutAtStart()), cutAtDelta(tt.cutAtDelta());
00440   bool cutInside(bool(tt.cutInside()));
00441   string name=tt.name().name();
00442 
00443   // check the parameters
00444   if (rIn <= 0 || rOut <=0 || cutAtStart <=0 || cutAtDelta <= 0) {
00445     throw cms::Exception("DetectorDescriptionFault", "TruncTubs " + string(tt.name().fullname()) + ": 0 <= rIn,cutAtStart,rOut,cutAtDelta,rOut violated!");
00446   }
00447   if (rIn >= rOut) {
00448     throw cms::Exception("DetectorDescriptionFault", "TruncTubs " + string(tt.name().fullname()) + ": rIn<rOut violated!");
00449   }
00450   if (startPhi != 0.) {
00451     throw cms::Exception("DetectorDescriptionFault", "TruncTubs " + string(tt.name().fullname()) + ": startPhi != 0 not supported!");
00452   }
00453   //     if (cutInside != false) {
00454   //       throw cms::Exception("DetectorDescriptionFault", "TruncTubs " + string(tt.name()) + " cutInside == true not supported!");
00455   //     }
00456 
00457   startPhi=0.;
00458   double r(cutAtStart), R(cutAtDelta);
00459   G4VSolid * result(0);
00460   G4VSolid * tubs = new G4Tubs(name,rIn,rOut,zHalf,startPhi,deltaPhi);
00461   LogDebug("SimG4CoreGeometry") << "G4Tubs: " << rIn << ' ' << rOut << ' ' << zHalf << ' ' << startPhi/deg << ' ' << deltaPhi/deg;
00462   LogDebug("SimG4CoreGeometry") << s;
00463   // length & hight of the box 
00464   double boxX(30.*rOut), boxY(20.*rOut); // exaggerate dimensions - does not matter, it's subtracted!
00465    
00466   // width of the box > width of the tubs
00467   double boxZ(1.1*zHalf);
00468    
00469   // angle of the box w.r.t. tubs
00470   double cath = r-R*cos(deltaPhi);
00471   double hypo = sqrt(r*r+R*R-2.*r*R*cos(deltaPhi));
00472   double cos_alpha = cath/hypo;
00473 
00474   double alpha = -acos(cos_alpha);
00475   LogDebug("SimG4CoreGeometry") << "cath=" << cath/m;
00476   LogDebug("SimG4CoreGeometry") << "hypo=" << hypo/m;
00477   LogDebug("SimG4CoreGeometry") << "al=" << acos(cath/hypo)/deg;
00478   LogDebug("SimG4CoreGeometry") << "deltaPhi=" << deltaPhi/deg << "\n"
00479                                 << "r=" << r/m << "\n"
00480                                 <<  "R=" << R/m;
00481 
00482   LogDebug("SimG4CoreGeometry") << "alpha=" << alpha/deg;
00483     
00484   // rotationmatrix of box w.r.t. tubs
00485   G4RotationMatrix * rot = new G4RotationMatrix;
00486   rot->rotateZ(-alpha);
00487   LogDebug("SimG4CoreGeometry") << (*rot);
00488 
00489   // center point of the box
00490   double xBox;
00491   if (!cutInside) {
00492     xBox = r+boxY/sin(abs(alpha));
00493   } else {
00494     xBox = -(boxY/sin(abs(alpha))-r);
00495   }
00496 
00497   G4ThreeVector trans(xBox,0.,0.);
00498   LogDebug("SimG4CoreGeometry") << "trans=" << trans;
00499 
00500   G4VSolid * box = new G4Box(name,boxX,boxY,boxZ);
00501   result = new G4SubtractionSolid(name,tubs,box,rot,trans);
00502       
00503   return result;
00504 
00505 }
00506 
00507 #include "G4Sphere.hh"
00508 G4VSolid * DDG4SolidConverter::sphere(const DDSolid & s) 
00509 {
00510    LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: sphere = " << s ;  
00511    DDSphere sp(s);
00512    return new G4Sphere(s.name().name(), sp.innerRadius(),
00513                        sp.outerRadius(),
00514                        sp.startPhi(),
00515                        sp.deltaPhi(),
00516                        sp.startTheta(),
00517                        sp.deltaTheta());
00518 }
00519 
00520 #include "G4Orb.hh"
00521 G4VSolid * DDG4SolidConverter::orb(const DDSolid & s) 
00522 {
00523    LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: orb = " << s ;  
00524    DDOrb sp(s);
00525    return new G4Orb(s.name().name(), sp.radius());
00526 }
00527 
00528 #include "G4EllipticalTube.hh"
00529 G4VSolid * DDG4SolidConverter::ellipticaltube(const DDSolid & s) 
00530 {
00531    LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: ellipticaltube = " << s ;  
00532    DDEllipticalTube sp(s);
00533    return new G4EllipticalTube(s.name().name(),
00534                                sp.xSemiAxis(),
00535                                sp.ySemiAxis(),
00536                                sp.zHeight());
00537 }
00538 
00539 #include "G4Ellipsoid.hh"
00540 G4VSolid * DDG4SolidConverter::ellipsoid(const DDSolid & s) 
00541 {
00542    LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: ellipsoid = " << s ;  
00543    DDEllipsoid sp(s);
00544    return new G4Ellipsoid(s.name().name(),
00545                           sp.xSemiAxis(),
00546                           sp.ySemiAxis(),
00547                           sp.zSemiAxis(),
00548                           sp.zBottomCut(),
00549                           sp.zTopCut());
00550 }
00551 
00552 #include "G4Para.hh"
00553 G4VSolid * DDG4SolidConverter::para(const DDSolid & s) 
00554 {
00555    LogDebug("SimG4CoreGeometry") << "DDG4SolidConverter: parallelepiped = " << s ;  
00556    DDParallelepiped sp(s);
00557    return new G4Para(s.name().name(),
00558                      sp.xHalf(),
00559                      sp.yHalf(),
00560                      sp.zHalf(),
00561                      sp.alpha(),
00562                      sp.theta(),
00563                      sp.phi());
00564 }
00565