Newer
Older
Markus Frank
committed
//==========================================================================
Markus Frank
committed
//--------------------------------------------------------------------------
// Copyright (C) Organisation europeenne pour la Recherche nucleaire (CERN)
Markus Frank
committed
// All rights reserved.
Markus Frank
committed
// For the licensing terms see $DD4hepINSTALL/LICENSE.
// For the list of contributors see $DD4hepINSTALL/doc/CREDITS.
Markus Frank
committed
// Author : M.Frank
//
//==========================================================================
Markus Frank
committed
// Framework include files
#include "DD4hep/Plugins.h"
#include "DD4hep/Shapes.h"
Markus Frank
committed
#include "DD4hep/Printout.h"
#include "DD4hep/DD4hepUnits.h"
#include "DD4hep/PropertyTable.h"
#include "DD4hep/detail/ShapesInterna.h"
#include "DD4hep/detail/ObjectsInterna.h"
#include "DD4hep/detail/DetectorInterna.h"
#include "DDG4/Geant4Field.h"
#include "DDG4/Geant4Converter.h"
#include "Geant4ShapeConverter.h"
Markus FRANK
committed
#include "TMath.h"
Markus FRANK
committed
//#include "TColor.h"
//#include "TGeoManager.h"
#include "G4Version.hh"
#include "G4VisAttributes.hh"
#include "G4ProductionCuts.hh"
#include "G4VUserRegionInformation.hh"
#include "G4Ellipsoid.hh"
#include "G4UnionSolid.hh"
#include "G4ReflectedSolid.hh"
#include "G4SubtractionSolid.hh"
#include "G4IntersectionSolid.hh"
#include "G4VSensitiveDetector.hh"
Markus FRANK
committed
#include "G4Element.hh"
Markus FRANK
committed
#include "G4Material.hh"
#include "G4UserLimits.hh"
#include "G4LogicalVolume.hh"
#include "G4ReflectionFactory.hh"
#include "G4OpticalSurface.hh"
#include "G4LogicalSkinSurface.hh"
Markus FRANK
committed
#include "G4ElectroMagneticField.hh"
#include "G4LogicalBorderSurface.hh"
#include "G4MaterialPropertiesTable.hh"
#include "G4MaterialPropertiesIndex.hh"
#if G4VERSION_NUMBER >= 1030
#include "G4ScaledSolid.hh"
#endif
#include "CLHEP/Units/SystemOfUnits.h"
Markus Frank
committed
// C/C++ include files
Markus Frank
committed
#include <iostream>
#include <iomanip>
#include <sstream>
namespace units = dd4hep;
using namespace dd4hep::detail;
using namespace dd4hep::sim;
using namespace dd4hep;
#include "DDG4/Geant4AssemblyVolume.h"
#include "DD4hep/DetectorTools.h"
static const double CM_2_MM = (CLHEP::centimeter/dd4hep::centimeter);
namespace {
static string indent = "";
static Double_t s_identity_rot[] = { 1., 0., 0., 0., 1., 0., 0., 0., 1. };
struct MyTransform3D : public G4Transform3D {
MyTransform3D(double XX, double XY, double XZ, double DX, double YX, double YY, double YZ, double DY, double ZX, double ZY,
Markus Frank
committed
double ZZ, double DZ)
: G4Transform3D(XX, XY, XZ, DX, YX, YY, YZ, DY, ZX, ZY, ZZ, DZ) {
Markus Frank
committed
MyTransform3D(const double* t, const double* r = s_identity_rot)
: G4Transform3D(r[0],r[1],r[2],t[0]*CM_2_MM,r[3],r[4],r[5],t[1]*CM_2_MM,r[6],r[7],r[8],t[2]*CM_2_MM) {
}
MyTransform3D(Transform3D&& copy) : Transform3D(copy) {}
bool is_left_handed(const TGeoMatrix* m) {
const Double_t* r = m->GetRotationMatrix();
if ( r ) {
Double_t det =
r[0]*r[4]*r[8] + r[3]*r[7]*r[2] + r[6]*r[1]*r[5] -
r[2]*r[4]*r[6] - r[5]*r[7]*r[0] - r[8]*r[1]*r[3];
return det < 0e0;
}
return false;
}
class G4UserRegionInformation : public G4VUserRegionInformation {
public:
Region region;
double threshold;
bool storeSecondaries;
G4UserRegionInformation()
Markus Frank
committed
: threshold(0.0), storeSecondaries(false) {
}
virtual ~G4UserRegionInformation() {
}
if (region.isValid())
printout(DEBUG, "Region", "Name:%s", region.name());
pair<double,double> g4PropertyConversion(int index) {
switch(index) {
case kRINDEX: return make_pair(CLHEP::keV/units::keV, 1.0);
case kREFLECTIVITY: return make_pair(CLHEP::keV/units::keV, 1.0);
case kREALRINDEX: return make_pair(CLHEP::keV/units::keV, 1.0);
case kIMAGINARYRINDEX: return make_pair(CLHEP::keV/units::keV, 1.0);
case kEFFICIENCY: return make_pair(CLHEP::keV/units::keV, 1.0);
case kTRANSMITTANCE: return make_pair(CLHEP::keV/units::keV, 1.0);
case kSPECULARLOBECONSTANT: return make_pair(CLHEP::keV/units::keV, 1.0);
case kSPECULARSPIKECONSTANT: return make_pair(CLHEP::keV/units::keV, 1.0);
case kBACKSCATTERCONSTANT: return make_pair(CLHEP::keV/units::keV, 1.0);
case kGROUPVEL: return make_pair(CLHEP::keV/units::keV, (CLHEP::m/CLHEP::s)/(units::m/units::s)); // meter/second
case kMIEHG: return make_pair(CLHEP::keV/units::keV, CLHEP::m/units::m);
case kRAYLEIGH: return make_pair(CLHEP::keV/units::keV, CLHEP::m/units::m); // ??? says its a length
case kWLSCOMPONENT: return make_pair(CLHEP::keV/units::keV, 1.0);
case kWLSABSLENGTH: return make_pair(CLHEP::keV/units::keV, CLHEP::m/units::m);
case kABSLENGTH: return make_pair(CLHEP::keV/units::keV, CLHEP::m/units::m);
#if G4VERSION_NUMBER >= 1100
case kWLSCOMPONENT2: return make_pair(CLHEP::keV/units::keV, 1.0);
case kWLSABSLENGTH2: return make_pair(CLHEP::keV/units::keV, CLHEP::m/units::m);
case kSCINTILLATIONCOMPONENT1: return make_pair(CLHEP::keV/units::keV, units::keV/CLHEP::keV);
case kSCINTILLATIONCOMPONENT2: return make_pair(CLHEP::keV/units::keV, units::keV/CLHEP::keV);
case kSCINTILLATIONCOMPONENT3: return make_pair(CLHEP::keV/units::keV, units::keV/CLHEP::keV);
#else
case kFASTCOMPONENT: return make_pair(CLHEP::keV/units::keV, 1.0);
case kSLOWCOMPONENT: return make_pair(CLHEP::keV/units::keV, 1.0);
case kPROTONSCINTILLATIONYIELD: return make_pair(CLHEP::keV/units::keV, units::keV/CLHEP::keV); // Yields: 1/energy
case kDEUTERONSCINTILLATIONYIELD: return make_pair(CLHEP::keV/units::keV, units::keV/CLHEP::keV);
case kTRITONSCINTILLATIONYIELD: return make_pair(CLHEP::keV/units::keV, units::keV/CLHEP::keV);
case kALPHASCINTILLATIONYIELD: return make_pair(CLHEP::keV/units::keV, units::keV/CLHEP::keV);
case kIONSCINTILLATIONYIELD: return make_pair(CLHEP::keV/units::keV, units::keV/CLHEP::keV);
case kELECTRONSCINTILLATIONYIELD: return make_pair(CLHEP::keV/units::keV, units::keV/CLHEP::keV);
default:
break;
}
printout(FATAL,"Geant4Converter", "+++ Cannot convert material property with index: %d", index);
#else
printout(FATAL,"Geant4Converter", "+++ Cannot convert material property with index: %d [Need Geant4 > 10.03]", index);
#endif
return make_pair(0e0,0e0);
}
double g4ConstPropertyConversion(int index) {
case kSURFACEROUGHNESS: return CLHEP::m/units::m; // Length
case kISOTHERMAL_COMPRESSIBILITY: return (CLHEP::m3/CLHEP::keV)/(units::m3/CLHEP::keV); // Volume/Energy
case kRS_SCALE_FACTOR: return 1.0; // ??
case kWLSMEANNUMBERPHOTONS: return 1.0; // ??
case kWLSTIMECONSTANT: return CLHEP::second/units::second; // Time
case kMIEHG_FORWARD: return 1.0;
case kMIEHG_BACKWARD: return 1.0;
case kMIEHG_FORWARD_RATIO: return 1.0;
case kSCINTILLATIONYIELD: return units::keV/CLHEP::keV; // Energy
case kRESOLUTIONSCALE: return 1.0;
case kFERMIPOT: return CLHEP::keV/units::keV; // Energy
case kDIFFUSION: return 1.0;
case kSPINFLIP: return 1.0;
case kLOSS: return 1.0; // ??
case kLOSSCS: return CLHEP::barn/units::barn; // ??
case kABSCS: return CLHEP::barn/units::barn; // ??
case kSCATCS: return CLHEP::barn/units::barn; // ??
case kMR_NBTHETA: return 1.0;
case kMR_NBE: return 1.0;
case kMR_RRMS: return 1.0; // ??
case kMR_CORRLEN: return CLHEP::m/units::m; // Length
case kMR_THETAMIN: return 1.0;
case kMR_THETAMAX: return 1.0;
case kMR_EMIN: return CLHEP::keV/units::keV; // Energy
case kMR_EMAX: return CLHEP::keV/units::keV; // Energy
case kMR_ANGNOTHETA: return 1.0;
case kMR_ANGNOPHI: return 1.0;
case kMR_ANGCUT: return 1.0;
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
#if G4VERSION_NUMBER >= 1100
case kSCINTILLATIONTIMECONSTANT1: return CLHEP::second/units::second; // Time
case kSCINTILLATIONTIMECONSTANT2: return CLHEP::second/units::second; // Time
case kSCINTILLATIONTIMECONSTANT3: return CLHEP::second/units::second; // Time
case kSCINTILLATIONRISETIME1: return CLHEP::second/units::second; // Time
case kSCINTILLATIONRISETIME2: return CLHEP::second/units::second; // Time
case kSCINTILLATIONRISETIME3: return CLHEP::second/units::second; // Time
case kSCINTILLATIONYIELD1: return 1.0;
case kSCINTILLATIONYIELD2: return 1.0;
case kSCINTILLATIONYIELD3: return 1.0;
case kPROTONSCINTILLATIONYIELD1: return 1.0;
case kPROTONSCINTILLATIONYIELD2: return 1.0;
case kPROTONSCINTILLATIONYIELD3: return 1.0;
case kDEUTERONSCINTILLATIONYIELD1: return 1.0;
case kDEUTERONSCINTILLATIONYIELD2: return 1.0;
case kDEUTERONSCINTILLATIONYIELD3: return 1.0;
case kALPHASCINTILLATIONYIELD1: return 1.0;
case kALPHASCINTILLATIONYIELD2: return 1.0;
case kALPHASCINTILLATIONYIELD3: return 1.0;
case kIONSCINTILLATIONYIELD1: return 1.0;
case kIONSCINTILLATIONYIELD2: return 1.0;
case kIONSCINTILLATIONYIELD3: return 1.0;
case kELECTRONSCINTILLATIONYIELD1: return 1.0;
case kELECTRONSCINTILLATIONYIELD2: return 1.0;
case kELECTRONSCINTILLATIONYIELD3: return 1.0;
#else
case kFASTTIMECONSTANT: return CLHEP::second/units::second; // Time
case kFASTSCINTILLATIONRISETIME: return CLHEP::second/units::second; // Time
case kSLOWTIMECONSTANT: return CLHEP::second/units::second; // Time
case kSLOWSCINTILLATIONRISETIME: return CLHEP::second/units::second; // Time
case kYIELDRATIO: return 1.0;
#endif
default:
break;
}
printout(FATAL,"Geant4Converter", "+++ Cannot convert CONST material property with index: %d", index);
#else
printout(FATAL,"Geant4Converter", "+++ Cannot convert material property with index: %d [Need Geant4 > 10.03]", index);
#endif
Geant4Converter::Geant4Converter(const Detector& description_ref)
: Geant4Mapping(description_ref), checkOverlaps(true) {
this->Geant4Mapping::init();
Markus Frank
committed
m_propagateRegions = true;
Markus Frank
committed
outputLevel = PrintLevel(printLevel() - 1);
Geant4Converter::Geant4Converter(const Detector& description_ref, PrintLevel level)
: Geant4Mapping(description_ref), checkOverlaps(true) {
Markus Frank
committed
m_propagateRegions = true;
Markus Frank
committed
outputLevel = level;
Markus Frank
committed
/// Standard destructor
Geant4Converter::~Geant4Converter() {
Markus Frank
committed
/// Handle the conversion of isotopes
void* Geant4Converter::handleIsotope(const string& /* name */, const TGeoIsotope* iso) const {
G4Isotope* g4i = data().g4Isotopes[iso];
double a_conv = (CLHEP::g / CLHEP::mole);
g4i = new G4Isotope(iso->GetName(), iso->GetZ(), iso->GetN(), iso->GetA()*a_conv);
printout(debugElements ? ALWAYS : outputLevel,
"Geant4Converter", "++ Created G4 Isotope %s from data: Z=%d N=%d A=%.3f [g/mole]",
iso->GetName(), iso->GetZ(), iso->GetN(), iso->GetA());
data().g4Isotopes[iso] = g4i;
}
return g4i;
}
/// Handle the conversion of elements
void* Geant4Converter::handleElement(const string& name, const Atom element) const {
G4Element* g4e = data().g4Elements[element];
PrintLevel lvl = debugElements ? ALWAYS : outputLevel;
if (element->GetNisotopes() > 0) {
g4e = new G4Element(name, element->GetTitle(), element->GetNisotopes());
for (int i = 0, n = element->GetNisotopes(); i < n; ++i) {
TGeoIsotope* iso = element->GetIsotope(i);
G4Isotope* g4iso = (G4Isotope*)handleIsotope(iso->GetName(), iso);
g4e->AddIsotope(g4iso, element->GetRelativeAbundance(i));
else {
// This adds in Geant4 the natural isotopes, which we normally do not want. We want to steer it outselves.
double a_conv = (CLHEP::g / CLHEP::mole);
g4e = new G4Element(element->GetTitle(), name, element->Z(), element->A()*a_conv);
printout(lvl, "Geant4Converter", "++ Created G4 Isotope %s from data: Z=%d N=%d A=%.3f [g/mole]",
element->GetName(), element->Z(), element->N(), element->A());
}
stringstream str;
str << (*g4e) << endl;
printout(lvl, "Geant4Converter", "++ Created G4 element %s", str.str().c_str());
data().g4Elements[element] = g4e;
}
return g4e;
}
/// Dump material in GDML format to output stream
void* Geant4Converter::handleMaterial(const string& name, Material medium) const {
Geant4GeometryInfo& info = data();
G4Material* mat = info.g4Materials[medium];
if ( !mat ) {
Markus Frank
committed
PrintLevel lvl = debugMaterials ? ALWAYS : outputLevel;
TGeoMaterial* material = medium->GetMaterial();
G4State state = kStateUndefined;
double density = material->GetDensity() * (CLHEP::gram / CLHEP::cm3);
if ( density < 1e-25 )
switch ( material->GetState() ) {
case TGeoMaterial::kMatStateSolid:
state = kStateSolid;
break;
case TGeoMaterial::kMatStateLiquid:
state = kStateLiquid;
break;
case TGeoMaterial::kMatStateGas:
state = kStateGas;
break;
default:
case TGeoMaterial::kMatStateUndefined:
state = kStateUndefined;
break;
}
if ( material->IsMixture() ) {
double A_total = 0.0;
double W_total = 0.0;
TGeoMixture* mix = (TGeoMixture*) material;
int nElements = mix->GetNelements();
mat = new G4Material(name, density, nElements, state,
material->GetTemperature(), material->GetPressure());
for (int i = 0; i < nElements; ++i) {
A_total += (mix->GetAmixt())[i];
W_total += (mix->GetWmixt())[i];
for (int i = 0; i < nElements; ++i) {
TGeoElement* e = mix->GetElement(i);
G4Element* g4e = (G4Element*) handleElement(e->GetName(), Atom(e));
if (!g4e) {
printout(ERROR, "Material",
"Missing component %s for material %s. A=%f W=%f",
e->GetName(), mix->GetName(), A_total, W_total);
Markus Frank
committed
}
//mat->AddElement(g4e, (mix->GetAmixt())[i] / A_total);
mat->AddElement(g4e, (mix->GetWmixt())[i] / W_total);
}
else {
double z = material->GetZ(), a = material->GetA();
if ( z < 1.0000001 ) z = 1.0;
if ( a < 0.5000001 ) a = 1.0;
mat = new G4Material(name, z, a, density, state,
material->GetTemperature(), material->GetPressure());
}
#if ROOT_VERSION_CODE >= ROOT_VERSION(6,17,0)
/// Attach the material properties if any
G4MaterialPropertiesTable* tab = 0;
TListIter propIt(&material->GetProperties());
for(TObject* obj=propIt.Next(); obj; obj = propIt.Next()) {
Markus Frank
committed
string exc_str;
TNamed* named = (TNamed*)obj;
TGDMLMatrix* matrix = info.manager->GetGDMLMatrix(named->GetTitle());
Geant4GeometryInfo::PropertyVector* v =
(Geant4GeometryInfo::PropertyVector*)handleMaterialProperties(matrix);
if ( 0 == v ) {
except("Geant4Converter", "++ FAILED to create G4 material %s [Cannot convert property:%s]",
material->GetName(), named->GetName());
if ( 0 == tab ) {
tab = new G4MaterialPropertiesTable();
mat->SetMaterialPropertiesTable(tab);
Markus Frank
committed
string err;
int idx = -1;
try {
idx = tab->GetPropertyIndex(named->GetName());
}
catch(const std::exception& e) {
exc_str = e.what();
idx = -1;
}
catch(...) {
idx = -1;
}
Markus Frank
committed
printout(ERROR, "Geant4Converter",
"++ UNKNOWN Geant4 CONST Property: %-20s %s [IGNORED]",
exc_str.c_str(), named->GetName());
continue;
}
// We need to convert the property from TGeo units to Geant4 units
auto conv = g4PropertyConversion(idx);
vector<double> bins(v->bins), vals(v->values);
for(size_t i=0, count=bins.size(); i<count; ++i)
bins[i] *= conv.first, vals[i] *= conv.second;
G4MaterialPropertyVector* vec = new G4MaterialPropertyVector(&bins[0], &vals[0], bins.size());
tab->AddProperty(named->GetName(), vec);
printout(lvl, "Geant4Converter", "++ Property: %-20s [%ld x %ld] -> %s ",
named->GetName(), matrix->GetRows(), matrix->GetCols(), named->GetTitle());
for(size_t i=0, count=v->bins.size(); i<count; ++i)
printout(lvl,named->GetName()," Geant4: %8.3g [MeV] TGeo: %8.3g [GeV] Conversion: %8.3g",
bins[i], v->bins[i], conv.first);
/// Attach the material properties if any
TListIter cpropIt(&material->GetConstProperties());
for(TObject* obj=cpropIt.Next(); obj; obj = cpropIt.Next()) {
Markus Frank
committed
string exc_str;
Bool_t err = kFALSE;
TNamed* named = (TNamed*)obj;
double v = info.manager->GetProperty(named->GetTitle(),&err);
if ( err != kFALSE ) {
except("Geant4Converter",
"++ FAILED to create G4 material %s "
"[Cannot convert const property: %s]",
material->GetName(), named->GetName());
}
if ( 0 == tab ) {
tab = new G4MaterialPropertiesTable();
mat->SetMaterialPropertiesTable(tab);
}
Markus Frank
committed
int idx = -1;
try {
idx = tab->GetPropertyIndex(named->GetName());
}
catch(const std::exception& e) {
exc_str = e.what();
idx = -1;
}
catch(...) {
idx = -1;
}
Markus Frank
committed
printout(ERROR, "Geant4Converter",
"++ UNKNOWN Geant4 CONST Property: %-20s %s [IGNORED]",
exc_str.c_str(), named->GetName());
continue;
}
// We need to convert the property from TGeo units to Geant4 units
double conv = g4ConstPropertyConversion(idx);
printout(lvl, "Geant4Converter", "++ CONST Property: %-20s %g ", named->GetName(), v);
tab->AddConstProperty(named->GetName(), v * conv);
}
#endif
auto* ionization = mat->GetIonisation();
stringstream str;
str << (*mat) << endl;
if ( ionization )
str << " log(MEE): " << ionization->GetLogMeanExcEnergy();
else
str << " log(MEE): UNKNOWN";
printout(lvl, "Geant4Converter", "++ Created G4 material %s", str.str().c_str());
info.g4Materials[medium] = mat;
void* Geant4Converter::handleSolid(const string& name, const TGeoShape* shape) const {
Markus Frank
committed
G4VSolid* solid = 0;
if ( shape ) {
if ( 0 != (solid = data().g4Solids[shape]) ) {
Markus Frank
committed
return solid;
}
TClass* isa = shape->IsA();
PrintLevel lvl = debugShapes ? ALWAYS : outputLevel;
if (isa == TGeoShapeAssembly::Class()) {
// Assemblies have no corresponding 'shape' in Geant4. Ignore the shape translation.
// It does not harm, since this 'shape' is never accessed afterwards.
data().g4Solids[shape] = solid = convertShape<TGeoShapeAssembly>(shape);
return solid;
Markus Frank
committed
}
else if (isa == TGeoBBox::Class())
solid = convertShape<TGeoBBox>(shape);
else if (isa == TGeoTube::Class())
solid = convertShape<TGeoTube>(shape);
else if (isa == TGeoTubeSeg::Class())
solid = convertShape<TGeoTubeSeg>(shape);
else if (isa == TGeoCtub::Class())
solid = convertShape<TGeoCtub>(shape);
else if (isa == TGeoEltu::Class())
solid = convertShape<TGeoEltu>(shape);
else if (isa == TwistedTubeObject::Class())
solid = convertShape<TwistedTubeObject>(shape);
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
else if (isa == TGeoTrd1::Class())
solid = convertShape<TGeoTrd1>(shape);
else if (isa == TGeoTrd2::Class())
solid = convertShape<TGeoTrd2>(shape);
else if (isa == TGeoHype::Class())
solid = convertShape<TGeoHype>(shape);
else if (isa == TGeoXtru::Class())
solid = convertShape<TGeoXtru>(shape);
else if (isa == TGeoPgon::Class())
solid = convertShape<TGeoPgon>(shape);
else if (isa == TGeoPcon::Class())
solid = convertShape<TGeoPcon>(shape);
else if (isa == TGeoCone::Class())
solid = convertShape<TGeoCone>(shape);
else if (isa == TGeoConeSeg::Class())
solid = convertShape<TGeoConeSeg>(shape);
else if (isa == TGeoParaboloid::Class())
solid = convertShape<TGeoParaboloid>(shape);
else if (isa == TGeoSphere::Class())
solid = convertShape<TGeoSphere>(shape);
else if (isa == TGeoTorus::Class())
solid = convertShape<TGeoTorus>(shape);
else if (isa == TGeoTrap::Class())
solid = convertShape<TGeoTrap>(shape);
else if (isa == TGeoArb8::Class())
solid = convertShape<TGeoArb8>(shape);
#if ROOT_VERSION_CODE > ROOT_VERSION(6,21,0)
else if (isa == TGeoTessellated::Class())
solid = convertShape<TGeoTessellated>(shape);
else if (isa == TGeoScaledShape::Class()) {
#if G4VERSION_NUMBER >= 1030
TGeoScaledShape* sh = (TGeoScaledShape*) shape;
const double* vals = sh->GetScale()->GetScale();
G4Scale3D scal(vals[0], vals[1], vals[2]);
G4VSolid* g4solid = (G4VSolid*)handleSolid(sol->GetName(), sol);
if ( scal.xx()>0e0 && scal.yy()>0e0 && scal.zz()>0e0 )
solid = new G4ScaledSolid(sh->GetName(), g4solid, scal);
else
solid = new G4ReflectedSolid(g4solid->GetName()+"_refl", g4solid, scal);
#else
except("Geant4Converter","++ TGeoScaledShape are only supported by Geant4 for versions >= 10.3");
#endif
}
else if ( isa == TGeoCompositeShape::Class() ) {
const TGeoCompositeShape* sh = (const TGeoCompositeShape*) shape;
const TGeoBoolNode* boolean = sh->GetBoolNode();
TGeoBoolNode::EGeoBoolType oper = boolean->GetBooleanOperator();
TGeoMatrix* matrix = boolean->GetRightMatrix();
G4VSolid* left = (G4VSolid*) handleSolid(name + "_left", boolean->GetLeftShape());
G4VSolid* right = (G4VSolid*) handleSolid(name + "_right", boolean->GetRightShape());
except("Geant4Converter","++ No left Geant4 Solid present for composite shape: %s",name.c_str());
except("Geant4Converter","++ No right Geant4 Solid present for composite shape: %s",name.c_str());
//specific case!
//Ellipsoid tag preparing
//if left == TGeoScaledShape AND right == TGeoBBox
// AND if TGeoScaledShape->GetShape == TGeoSphere
TGeoShape* ls = boolean->GetLeftShape();
TGeoShape* rs = boolean->GetRightShape();
if (strcmp(ls->ClassName(), "TGeoScaledShape") == 0 &&
Markus Frank
committed
strcmp(rs->ClassName(), "TGeoBBox") == 0) {
if (strcmp(((TGeoScaledShape *)ls)->GetShape()->ClassName(), "TGeoSphere") == 0) {
if (oper == TGeoBoolNode::kGeoIntersection) {
TGeoScaledShape* lls = (TGeoScaledShape *)ls;
TGeoBBox* rrs = (TGeoBBox*)rs;
Markus Frank
committed
double sx = lls->GetScale()->GetScale()[0];
double sy = lls->GetScale()->GetScale()[1];
Markus Frank
committed
double radius = ((TGeoSphere *)lls->GetShape())->GetRmax();
Markus Frank
committed
double dz = rrs->GetDZ();
double zorig = rrs->GetOrigin()[2];
double zcut2 = dz + zorig;
double zcut1 = 2 * zorig - zcut2;
Markus Frank
committed
solid = new G4Ellipsoid(name,
sx * radius * CM_2_MM,
sy * radius * CM_2_MM,
radius * CM_2_MM,
zcut1 * CM_2_MM,
zcut2 * CM_2_MM);
data().g4Solids[shape] = solid;
return solid;
}
}
if ( matrix->IsRotation() ) {
MyTransform3D transform(matrix->GetTranslation(),matrix->GetRotationMatrix());
if (oper == TGeoBoolNode::kGeoSubtraction)
solid = new G4SubtractionSolid(name, left, right, transform);
else if (oper == TGeoBoolNode::kGeoUnion)
solid = new G4UnionSolid(name, left, right, transform);
else if (oper == TGeoBoolNode::kGeoIntersection)
solid = new G4IntersectionSolid(name, left, right, transform);
const Double_t *t = matrix->GetTranslation();
G4ThreeVector transform(t[0] * CM_2_MM, t[1] * CM_2_MM, t[2] * CM_2_MM);
if (oper == TGeoBoolNode::kGeoSubtraction)
solid = new G4SubtractionSolid(name, left, right, 0, transform);
else if (oper == TGeoBoolNode::kGeoUnion)
solid = new G4UnionSolid(name, left, right, 0, transform);
else if (oper == TGeoBoolNode::kGeoIntersection)
solid = new G4IntersectionSolid(name, left, right, 0, transform);
except("Geant4Converter","++ Failed to handle unknown solid shape: %s of type %s",
name.c_str(), isa->GetName());
printout(lvl,"Geant4Converter","++ Successessfully converted shape [%p] of type:%s to %s.",
solid,isa->GetName(),typeName(typeid(*solid)).c_str());
data().g4Solids[shape] = solid;
}
return solid;
}
/// Dump logical volume in GDML format to output stream
void* Geant4Converter::handleVolume(const string& name, const TGeoVolume* volume) const {
Geant4GeometryInfo& info = data();
PrintLevel lvl = debugVolumes ? ALWAYS : outputLevel;
Markus Frank
committed
Geant4GeometryMaps::VolumeMap::const_iterator volIt = info.g4Volumes.find(volume);
Volume _v(volume);
if ( _v.testFlagBit(Volume::VETO_SIMU) ) {
printout(lvl, "Geant4Converter", "++ Volume %s not converted [Veto'ed for simulation]",volume->GetName());
return 0;
}
else if (volIt == info.g4Volumes.end() ) {
string n = volume->GetName();
TGeoMedium* med = volume->GetMedium();
TGeoShape* sh = volume->GetShape();
G4VSolid* solid = (G4VSolid*) handleSolid(sh->GetName(), sh);
bool is_assembly = sh->IsA() == TGeoShapeAssembly::Class() || volume->IsA() == TGeoVolumeAssembly::Class();
printout(lvl, "Geant4Converter", "++ Convert Volume %-32s: %p %s/%s assembly:%s",
n.c_str(), volume, sh->IsA()->GetName(), volume->IsA()->GetName(), yes_no(is_assembly));
if ( is_assembly ) {
return nullptr;
}
Region reg = _v.region();
LimitSet lim = _v.limitSet();
VisAttr vis = _v.visAttributes();
G4Region* region = reg.isValid() ? info.g4Regions[reg] : nullptr;
G4UserLimits* limits = lim.isValid() ? info.g4Limits[lim] : nullptr;
G4Material* medium = (G4Material*) handleMaterial(med->GetName(), Material(med));
/// Check all pre-conditions
if ( !solid ) {
except("G4Converter","++ No Geant4 Solid present for volume:" + n);
Markus Frank
committed
}
else if ( !medium ) {
except("G4Converter","++ No Geant4 material present for volume:" + n);
Markus Frank
committed
}
else if ( reg.isValid() && !region ) {
except("G4Cnv::volume["+name+"]"," ++ Failed to access Geant4 region %s.",reg.name());
Markus Frank
committed
}
else if ( lim.isValid() && !limits ) {
except("G4Cnv::volume["+name+"]","++ FATAL Failed to access Geant4 user limits %s.",lim.name());
}
else if ( limits ) {
printout(lvl, "Geant4Converter", "++ Volume + Apply LIMITS settings:%-24s to volume %s.",
Markus Frank
committed
lim.name(), _v.name());
Markus Frank
committed
}
G4VisAttributes* vattr = vis.isValid() ? (G4VisAttributes*)handleVis(vis.name(), vis) : nullptr;
G4LogicalVolume* g4vol = new G4LogicalVolume(solid, medium, n, 0, 0, limits);
if ( region ) {
printout(lvl, "Geant4Converter", "++ Volume + Apply REGION settings: %s to volume %s.",
Markus Frank
committed
reg.name(), _v.name());
g4vol->SetRegion(region);
region->AddRootLogicalVolume(g4vol);
if ( vattr ) {
g4vol->SetVisAttributes(vattr);
Markus Frank
committed
}
info.g4Volumes[volume] = g4vol;
printout(lvl, "Geant4Converter", "++ Volume + %s converted: %p ---> G4: %p", n.c_str(), volume, g4vol);
return nullptr;
/// Dump logical volume in GDML format to output stream
void* Geant4Converter::collectVolume(const string& /* name */, const TGeoVolume* volume) const {
Geant4GeometryInfo& info = data();
Volume _v(volume);
Region reg = _v.region();
LimitSet lim = _v.limitSet();
SensitiveDetector det = _v.sensitiveDetector();
if ( lim.isValid() )
info.limits[lim].insert(volume);
if ( reg.isValid() )
info.regions[reg].insert(volume);
if ( det.isValid() )
info.sensitives[det].insert(volume);
return (void*)volume;
/// Dump volume placement in GDML format to output stream
void* Geant4Converter::handleAssembly(const string& name, const TGeoNode* node) const {
TGeoVolume* mot_vol = node->GetVolume();
PrintLevel lvl = debugVolumes ? ALWAYS : outputLevel;
if ( mot_vol->IsA() != TGeoVolumeAssembly::Class() ) {
Volume _v(mot_vol);
if ( _v.testFlagBit(Volume::VETO_SIMU) ) {
printout(lvl, "Geant4Converter", "++ AssemblyNode %s not converted [Veto'ed for simulation]",node->GetName());
return 0;
}
Geant4GeometryInfo& info = data();
Geant4AssemblyVolume* g4 = info.g4AssemblyVolumes[node];
if ( g4 ) {
printout(ALWAYS, "Geant4Converter", "+++ Assembly: **** : Re-using existing assembly: %s",node->GetName());
}
if ( !g4 ) {
g4 = new Geant4AssemblyVolume();
for(Int_t i=0; i < mot_vol->GetNdaughters(); ++i) {
TGeoNode* dau = mot_vol->GetNode(i);
TGeoVolume* dau_vol = dau->GetVolume();
TGeoMatrix* tr = dau->GetMatrix();
MyTransform3D transform(tr->GetTranslation(),tr->IsRotation() ? tr->GetRotationMatrix() : s_identity_rot);
if ( is_left_handed(tr) ) {
G4Scale3D scale;
G4Rotate3D rot;
G4Translate3D trans;
transform.getDecomposition(scale, rot, trans);
printout(debugReflections ? ALWAYS : lvl, "Geant4Converter",
"Tr:x=%8.1f y=%8.1f z=%8.1f Scale:x=%4.2f y=%4.2f z=%4.2f",
dau_vol->GetName(), mot_vol->GetName(),
transform.dx(), transform.dy(), transform.dz(),
scale.xx(), scale.yy(), scale.zz());
if ( dau_vol->IsA() == TGeoVolumeAssembly::Class() ) {
Geant4GeometryMaps::AssemblyMap::iterator ia = info.g4AssemblyVolumes.find(dau);
if ( ia == info.g4AssemblyVolumes.end() ) {
Markus Frank
committed
printout(FATAL, "Geant4Converter", "+++ Invalid child assembly at %s : %d parent: %s child:%s",
__FILE__, __LINE__, name.c_str(), dau->GetName());
Markus Frank
committed
}
g4->placeAssembly(dau, (*ia).second, transform);
printout(lvl, "Geant4Converter", "+++ Assembly: AddPlacedAssembly : dau:%s "
Markus Frank
committed
"to mother %s Tr:x=%8.3f y=%8.3f z=%8.3f",
dau_vol->GetName(), mot_vol->GetName(),
transform.dx(), transform.dy(), transform.dz());
Geant4GeometryMaps::VolumeMap::iterator iv = info.g4Volumes.find(dau_vol);
if ( iv == info.g4Volumes.end() ) {
printout(FATAL,"Geant4Converter", "+++ Invalid child volume at %s : %d parent: %s child:%s",
__FILE__, __LINE__, name.c_str(), dau->GetName());
except("Geant4Converter", "+++ Invalid child volume at %s : %d parent: %s child:%s",
__FILE__, __LINE__, name.c_str(), dau->GetName());
Markus Frank
committed
}
g4->placeVolume(dau,(*iv).second, transform);
printout(lvl, "Geant4Converter", "+++ Assembly: AddPlacedVolume : dau:%s "
Markus Frank
committed
"to mother %s Tr:x=%8.3f y=%8.3f z=%8.3f",
dau_vol->GetName(), mot_vol->GetName(),
transform.dx(), transform.dy(), transform.dz());
info.g4AssemblyVolumes[node] = g4;
/// Dump volume placement in GDML format to output stream
void* Geant4Converter::handlePlacement(const string& name, const TGeoNode* node) const {
Geant4GeometryInfo& info = data();
PrintLevel lvl = debugPlacements ? ALWAYS : outputLevel;
Markus Frank
committed
Geant4GeometryMaps::PlacementMap::const_iterator g4it = info.g4Placements.find(node);
G4VPhysicalVolume* g4 = (g4it == info.g4Placements.end()) ? 0 : (*g4it).second;
TGeoVolume* vol = node->GetVolume();
Volume _v(vol);
if ( _v.testFlagBit(Volume::VETO_SIMU) ) {
printout(lvl, "Geant4Converter", "++ Placement %s not converted [Veto'ed for simulation]",node->GetName());
return 0;
}
TGeoVolume* mot_vol = node->GetMotherVolume();
TGeoMatrix* tr = node->GetMatrix();
except("Geant4Converter",
"++ Attempt to handle placement without transformation:%p %s of type %s vol:%p",
node, node->GetName(), node->IsA()->GetName(), vol);
except("Geant4Converter", "++ Unknown G4 volume:%p %s of type %s ptr:%p",
node, node->GetName(), node->IsA()->GetName(), vol);
int copy = node->GetNumber();
bool node_is_reflected = is_left_handed(tr);
bool node_is_assembly = vol->IsA() == TGeoVolumeAssembly::Class();
bool mother_is_assembly = mot_vol ? mot_vol->IsA() == TGeoVolumeAssembly::Class() : false;
MyTransform3D transform(tr->GetTranslation(),tr->IsRotation() ? tr->GetRotationMatrix() : s_identity_rot);
Markus Frank
committed
Geant4GeometryMaps::VolumeMap::const_iterator volIt = info.g4Volumes.find(mot_vol);
if ( mother_is_assembly ) {
Markus Frank
committed
//
// Mother is an assembly:
// Nothing to do here, because:
// -- placed volumes were already added before in "handleAssembly"
// -- imprint cannot be made, because this requires a logical volume as a mother
//
printout(lvl, "Geant4Converter", "+++ Assembly: **** : dau:%s "
Markus Frank
committed
"to mother %s Tr:x=%8.3f y=%8.3f z=%8.3f",
vol->GetName(), mot_vol->GetName(),
transform.dx(), transform.dy(), transform.dz());
Markus Frank
committed
}
G4Scale3D scale;
G4Rotate3D rot;
G4Translate3D trans;
transform.getDecomposition(scale, rot, trans);
if ( node_is_assembly ) {
Markus Frank
committed
//
// Node is an assembly:
// Imprint the assembly. The mother MUST already be transformed.
//
printout(lvl, "Geant4Converter", "++ Assembly: makeImprint: dau:%-12s %s in mother %-12s "
"Tr:x=%8.1f y=%8.1f z=%8.1f Scale:x=%4.2f y=%4.2f z=%4.2f",
node->GetName(), node_is_reflected ? "(REFLECTED)" : "",
mot_vol ? mot_vol->GetName() : "<unknown>",
transform.dx(), transform.dy(), transform.dz(),
scale.xx(), scale.yy(), scale.zz());
Markus Frank
committed
Geant4AssemblyVolume* ass = (Geant4AssemblyVolume*)info.g4AssemblyVolumes[node];
Geant4AssemblyVolume::Chain chain;
chain.emplace_back(node);
ass->imprint(*this, node, chain, ass, (*volIt).second, transform, copy, checkOverlaps);
else if ( node != info.manager->GetTopNode() && volIt == info.g4Volumes.end() ) {
Markus Frank
committed
throw logic_error("Geant4Converter: Invalid mother volume found!");
G4LogicalVolume* g4vol = info.g4Volumes[vol];
G4LogicalVolume* g4mot = info.g4Volumes[mot_vol];
G4PhysicalVolumesPair pvPlaced =
G4ReflectionFactory::Instance()->Place(transform, // no rotation
name, // its name
g4vol, // its logical volume
g4mot, // its mother (logical) volume
false, // no boolean operations
copy, // its copy number
checkOverlaps);
printout(debugReflections||debugPlacements ? ALWAYS : lvl, "Geant4Converter",
"++ Place %svolume %-12s in mother %-12s "
"Tr:x=%8.1f y=%8.1f z=%8.1f Scale:x=%4.2f y=%4.2f z=%4.2f",
node_is_reflected ? "REFLECTED " : "", _v.name(),
mot_vol ? mot_vol->GetName() : "<unknown>",
transform.dx(), transform.dy(), transform.dz(),
scale.xx(), scale.yy(), scale.zz());
// First 2 cases can be combined.
// Leave them separated for debugging G4ReflectionFactory for now...
if ( node_is_reflected && !pvPlaced.second )
return info.g4Placements[node] = pvPlaced.first;
else if ( !node_is_reflected && !pvPlaced.second )
return info.g4Placements[node] = pvPlaced.first;
// Now deal with valid pvPlaced.second ...
if ( node_is_reflected )
return info.g4Placements[node] = pvPlaced.first;
else if ( !node_is_reflected )
return info.g4Placements[node] = pvPlaced.first;
g4 = pvPlaced.second ? pvPlaced.second : pvPlaced.first;
info.g4Placements[node] = g4;
printout(ERROR, "Geant4Converter", "++ DEAD code. Should not end up here!");
/// Convert the geometry type region into the corresponding Geant4 object(s).
void* Geant4Converter::handleRegion(Region region, const set<const TGeoVolume*>& /* volumes */) const {
G4Region* g4 = data().g4Regions[region];
Markus Frank
committed
PrintLevel lvl = debugRegions ? ALWAYS : outputLevel;
Markus Frank
committed
g4 = new G4Region(r.name());
// create region info with storeSecondaries flag
Andre Sailer
committed
if( not r.wasThresholdSet() and r.storeSecondaries() ) {
throw runtime_error("G4Region: StoreSecondaries is True, but no explicit threshold set:");
}
printout(lvl, "Geant4Converter", "++ Setting up region: %s", r.name());
G4UserRegionInformation* info = new G4UserRegionInformation();
info->region = r;
info->threshold = r.threshold()*CLHEP::MeV/units::MeV;
info->storeSecondaries = r.storeSecondaries();
g4->SetUserInformation(info);
printout(lvl, "Geant4Converter", "++ Converted region settings of:%s.", r.name());
vector < string > &limits = r.limits();
G4ProductionCuts* cuts = 0;
// set production cut
if( not r.useDefaultCut() ) {
cuts = new G4ProductionCuts();
cuts->SetProductionCut(r.cut()*CLHEP::mm/units::mm);
printout(lvl, "Geant4Converter", "++ %s: Using default cut: %f [mm]",
r.name(), r.cut()*CLHEP::mm/units::mm);
}
for( const auto& nam : limits ) {
const LimitSet::Set& cts = ls.cuts();
for (const auto& c : cts ) {
if ( c.particles == "*" ) pid = -1;
else if ( c.particles == "e-" ) pid = idxG4ElectronCut;
else if ( c.particles == "e+" ) pid = idxG4PositronCut;
else if ( c.particles == "e[+-]" ) pid = -idxG4PositronCut-idxG4ElectronCut;
else if ( c.particles == "e[-+]" ) pid = -idxG4PositronCut-idxG4ElectronCut;
else if ( c.particles == "gamma" ) pid = idxG4GammaCut;
else if ( c.particles == "proton" ) pid = idxG4ProtonCut;
else throw runtime_error("G4Region: Invalid production cut particle-type:" + c.particles);
if ( !cuts ) cuts = new G4ProductionCuts();
if ( pid == -(idxG4PositronCut+idxG4ElectronCut) ) {
cuts->SetProductionCut(c.value*CLHEP::mm/units::mm, idxG4PositronCut);
cuts->SetProductionCut(c.value*CLHEP::mm/units::mm, idxG4ElectronCut);
}
else {
cuts->SetProductionCut(c.value*CLHEP::mm/units::mm, pid);
}
printout(lvl, "Geant4Converter", "++ %s: Set cut [%s/%d] = %f [mm]",
r.name(), c.particles.c_str(), pid, c.value*CLHEP::mm/units::mm);
}
const auto& lm = data().g4Limits;
for (const auto& j : lm ) {
if (nam == j.first->GetName()) {
g4->SetUserLimits(j.second);
printout(lvl, "Geant4Converter", "++ %s: Set limits %s to region type %s",
r.name(), nam.c_str(), j.second->GetType().c_str());
found = true;
break;
}
}
if ( found ) {
throw runtime_error("G4Region: Failed to resolve user limitset:" + nam);
/// Assign cuts to region if they were created
if ( cuts ) g4->SetProductionCuts(cuts);
data().g4Regions[region] = g4;
}
return g4;
}
/// Convert the geometry type LimitSet into the corresponding Geant4 object(s).
void* Geant4Converter::handleLimitSet(LimitSet limitset, const set<const TGeoVolume*>& /* volumes */) const {
G4UserLimits* g4 = data().g4Limits[limitset];
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
struct LimitPrint {
const LimitSet& ls;
LimitPrint(const LimitSet& lset) : ls(lset) {}
const LimitPrint& operator()(const std::string& pref, const Geant4UserLimits::Handler& h) const {
if ( !h.particleLimits.empty() ) {
printout(ALWAYS,"Geant4Converter",
"+++ LimitSet: Limit %s.%s applied for particles:",ls.name(), pref.c_str());
for(const auto& p : h.particleLimits)
printout(ALWAYS,"Geant4Converter","+++ LimitSet: Particle type: %-18s PDG: %-6d : %f",
p.first->GetParticleName().c_str(), p.first->GetPDGEncoding(), p.second);
}
else {
printout(ALWAYS,"Geant4Converter",
"+++ LimitSet: Limit %s.%s NOT APPLIED.",ls.name(), pref.c_str());
}
return *this;
}
};
Geant4UserLimits* limits = new Geant4UserLimits(limitset);
g4 = limits;
if ( debugRegions ) {
LimitPrint print(limitset);
print("maxTime", limits->maxTime)
("minEKine", limits->minEKine)
("minRange", limits->minRange)
("maxStepLength", limits->maxStepLength)
("maxTrackLength",limits->maxTrackLength);
}
data().g4Limits[limitset] = g4;
}
return g4;
}
/// Convert the geometry visualisation attributes to the corresponding Geant4 object(s).
void* Geant4Converter::handleVis(const string& /* name */, VisAttr attr) const {
Geant4GeometryInfo& info = data();
G4VisAttributes* g4 = info.g4Vis[attr];
if ( !g4 ) {