CMS 3D CMS Logo

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

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