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
Markus Frank
committed
#include <DD4hep/Shapes.h>
#include <DD4hep/Volumes.h>
Markus Frank
committed
#include <DD4hep/Plugins.h>
Markus Frank
committed
#include <DD4hep/Printout.h>
Markus Frank
committed
#include <DD4hep/Detector.h>
Markus Frank
committed
#include <DD4hep/DD4hepUnits.h>
#include <DD4hep/PropertyTable.h>
Markus Frank
committed
#include <DD4hep/DetectorTools.h>
Markus Frank
committed
#include <DD4hep/detail/ShapesInterna.h>
#include <DD4hep/detail/ObjectsInterna.h>
#include <DD4hep/detail/DetectorInterna.h>
Markus Frank
committed
#include <DDG4/Geant4Field.h>
#include <DDG4/Geant4Helpers.h>
Markus Frank
committed
#include <DDG4/Geant4Converter.h>
#include <DDG4/Geant4UserLimits.h>
Markus Frank
committed
#include <DDG4/Geant4AssemblyVolume.h>
#include <DDG4/Geant4PlacementParameterisation.h>
#include "Geant4ShapeConverter.h"
#include <TClass.h>
Markus Frank
committed
#include <TGeoBoolNode.h>
Markus Frank
committed
#include <G4Version.hh>
#include <G4VisAttributes.hh>
#include <G4PVParameterised.hh>
Markus Frank
committed
#include <G4ProductionCuts.hh>
#include <G4VUserRegionInformation.hh>
Markus Frank
committed
#include <G4Box.hh>
#include <G4Tubs.hh>
#include <G4Ellipsoid.hh>
#include <G4UnionSolid.hh>
#include <G4ReflectedSolid.hh>
#include <G4SubtractionSolid.hh>
#include <G4IntersectionSolid.hh>
#include <G4VSensitiveDetector.hh>
Markus Frank
committed
#include <G4Region.hh>
#include <G4Element.hh>
#include <G4Isotope.hh>
#include <G4Material.hh>
#include <G4UserLimits.hh>
#include <G4RegionStore.hh>
Markus Frank
committed
#include <G4FieldManager.hh>
#include <G4LogicalVolume.hh>
#include <G4OpticalSurface.hh>
#include <G4ReflectionFactory.hh>
Markus Frank
committed
#include <G4LogicalSkinSurface.hh>
#include <G4ElectroMagneticField.hh>
#include <G4LogicalBorderSurface.hh>
#include <G4MaterialPropertiesTable.hh>
Markus Frank
committed
#include <G4MaterialPropertiesIndex.hh>
Markus Frank
committed
#include <G4ScaledSolid.hh>
#include <CLHEP/Units/SystemOfUnits.h>
Markus Frank
committed
// C/C++ include files
Markus Frank
committed
#include <iostream>
#include <iomanip>
#include <sstream>
Markus Frank
committed
namespace units = dd4hep;
using namespace dd4hep::sim;
using namespace dd4hep;
Markus Frank
committed
namespace {
Markus Frank
committed
  static constexpr const double CM_2_MM = (CLHEP::centimeter/dd4hep::centimeter);
  static constexpr const char* GEANT4_TAG_IGNORE = "Geant4-ignore";
  static constexpr const char* GEANT4_TAG_PLUGIN = "Geant4-plugin";
  static constexpr const char* GEANT4_TAG_BIRKSCONSTANT    = "BirksConstant";
  static constexpr const char* GEANT4_TAG_MEE              = "MeanExcitationEnergy";
  static constexpr const char* GEANT4_TAG_ENE_PER_ION_PAIR = "MeanEnergyPerIonPair";
Markus Frank
committed
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
  template <typename O, typename C, typename F> void handleRefs(const O* o, const C& c, F pmf) {
    for (typename C::const_iterator i = c.begin(); i != c.end(); ++i) {
      //(o->*pmf)((*i)->GetName(), *i);
      (o->*pmf)("", *i);
    }
  }
  template <typename O, typename C, typename F> void handle(const O* o, const C& c, F pmf) {
    for (typename C::const_iterator i = c.begin(); i != c.end(); ++i) {
      (o->*pmf)((*i)->GetName(), *i);
    }
  }
  template <typename O, typename F> void handleArray(const O* o, const TObjArray* c, F pmf) {
    TObjArrayIter arr(c);
    for(TObject* i = arr.Next(); i; i=arr.Next())
      (o->*pmf)(i);
  }
  template <typename O, typename C, typename F> void handleMap(const O* o, const C& c, F pmf) {
    for (typename C::const_iterator i = c.begin(); i != c.end(); ++i)
      (o->*pmf)((*i).first, (*i).second);
  }
  template <typename O, typename C, typename F> void handleRMap(const O* o, const C& c, F pmf) {
    for (typename C::const_reverse_iterator i = c.rbegin(); i != c.rend(); ++i)  {
      //cout << "Handle RMAP [ " << (*i).first << " ]" << endl;
      handle(o, (*i).second, pmf);
    }
  }
  template <typename O, typename C, typename F> void handleRMap_(const O* o, const C& c, F pmf) {
    for (typename C::const_iterator i = c.begin(); i != c.end(); ++i)  {
      const auto& cc = (*i).second;
      for (const auto& j : cc)   {
        (o->*pmf)(j);
      }
    }
  }
  string make_NCName(const string& in)   {
    string res = detail::str_replace(in, "/", "_");
    res = detail::str_replace(res, "#", "_");
    return res;
  }
  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;
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
#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 )  {
    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;
    }
Markus Frank
committed
    printout(lvl,"Geant4Material","+++ Setting up material %s", name.c_str());
    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) {
Markus Frank
committed
          printout(ERROR, name, 
                   "Missing element 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());
    }
Markus Frank
committed
    string plugin_name { };
    double value = 0e0;
    double ionisation_mee = -2e100;
    double ionisation_birks_constant = -2e100;
    double ionisation_ene_per_ion_pair = -2e100;
#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());
      const char*  cptr   = ::strstr(matrix->GetName(), GEANT4_TAG_IGNORE);
Markus Frank
committed
        printout(INFO,name,"++ Ignore property %s [%s]. Not Suitable for Geant4.",
                 matrix->GetName(), matrix->GetTitle());
        continue;
      }
      cptr = ::strstr(matrix->GetTitle(), GEANT4_TAG_IGNORE);
Markus Frank
committed
        printout(INFO,name,"++ Ignore property %s [%s]. Not Suitable for Geant4.",
                 matrix->GetName(), matrix->GetTitle());
        continue;
      Geant4GeometryInfo::PropertyVector* v =
        (Geant4GeometryInfo::PropertyVector*)handleMaterialProperties(matrix);
        except("Geant4Converter", "++ FAILED to create G4 material %s [Cannot convert property:%s]",
               material->GetName(), named->GetName());
        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",
Markus Frank
committed
                 "++ UNKNOWN Geant4 Property: %-20s %s [IGNORED]",
Markus Frank
committed
                 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(std::size_t i=0, count=bins.size(); i<count; ++i)
        bins[i] *= conv.first, vals[i] *= conv.second;
      G4MaterialPropertyVector* vec =
Markus Frank
committed
        new G4MaterialPropertyVector(&bins[0], &vals[0], bins.size());
      tab->AddProperty(named->GetName(), vec);
Markus Frank
committed
      printout(lvl, name, "++      Property: %-20s [%ld x %ld] -> %s ",
               named->GetName(), matrix->GetRows(), matrix->GetCols(), named->GetTitle());
      for(std::size_t i=0, count=v->bins.size(); i<count; ++i)
Markus Frank
committed
        printout(lvl, name, "  Geant4: %s %8.3g [MeV]  TGeo: %8.3g [GeV] Conversion: %8.3g",
                 named->GetName(), 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;
      const char*  cptr = ::strstr(named->GetName(), GEANT4_TAG_IGNORE);
Markus Frank
committed
        printout(INFO, name, "++ Ignore CONST property %s [%s].",
                 named->GetName(), named->GetTitle());
        continue;
      }
      cptr = ::strstr(named->GetTitle(), GEANT4_TAG_IGNORE);
Markus Frank
committed
        printout(INFO, name,"++ Ignore CONST property %s [%s].",
                 named->GetName(), named->GetTitle());
        continue;
      cptr = ::strstr(named->GetName(), GEANT4_TAG_PLUGIN);
        printout(INFO, name, "++ Ignore CONST property %s [%s]  --> Plugin.",
                 named->GetName(), named->GetTitle());
Markus Frank
committed
        plugin_name = named->GetTitle();
        continue;
      }
      cptr = ::strstr(named->GetName(), GEANT4_TAG_BIRKSCONSTANT);
Markus Frank
committed
        err = kFALSE;
        value = material->GetConstProperty(GEANT4_TAG_BIRKSCONSTANT,&err);
        if ( err == kFALSE ) ionisation_birks_constant = value * (CLHEP::mm/CLHEP::MeV)/(units::mm/units::MeV);
        continue;
      }
      cptr = ::strstr(named->GetName(), GEANT4_TAG_MEE);
Markus Frank
committed
        err = kFALSE;
        value = material->GetConstProperty(GEANT4_TAG_MEE, &err);
        if ( err == kFALSE ) ionisation_mee = value * (CLHEP::MeV/units::MeV);
        continue;
      }
      cptr = ::strstr(named->GetName(), GEANT4_TAG_ENE_PER_ION_PAIR);
Markus Frank
committed
        err = kFALSE;
        value = material->GetConstProperty(GEANT4_TAG_ENE_PER_ION_PAIR,&err);
        if ( err == kFALSE ) ionisation_ene_per_ion_pair = value * (CLHEP::MeV/units::MeV);
        continue;
      }
      err = kFALSE;
Markus Frank
committed
      value = info.manager->GetProperty(named->GetTitle(), &err);
Markus Frank
committed
        except(name,
               "++ FAILED to create G4 material %s [Cannot convert const property: %s]",
               material->GetName(), named->GetName());
      }
        tab = new G4MaterialPropertiesTable();
        mat->SetMaterialPropertiesTable(tab);
      }
Markus Frank
committed
      int idx = -1;
      try   {
        idx = tab->GetConstPropertyIndex(named->GetName());
Markus Frank
committed
      }
      catch(const std::exception& e)   {
        exc_str = e.what();
        idx = -1;
      }
      catch(...)   {
        idx = -1;
      }
Markus Frank
committed
        printout(ERROR, name,
Markus Frank
committed
                 "++ 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, name, "++      CONST Property: %-20s %g ", named->GetName(), value);
      tab->AddConstProperty(named->GetName(), value * conv);
    // Set Birk's constant if it was supplied in the material table of the TGeoMaterial
    auto* ionisation = mat->GetIonisation();
    str << (*mat);
    if ( ionisation )   {
      if ( ionisation_birks_constant > 0e0 )   {
Markus Frank
committed
        ionisation->SetBirksConstant(ionisation_birks_constant);
      }
      if ( ionisation_mee > -1e100 )   {
Markus Frank
committed
        ionisation->SetMeanExcitationEnergy(ionisation_mee);
      }
      if ( ionisation_ene_per_ion_pair > 0e0 )   {
Markus Frank
committed
        ionisation->SetMeanEnergyPerIonPair(ionisation_ene_per_ion_pair);
      }
      str << "          log(MEE): " << std::setprecision(4) << ionisation->GetLogMeanExcEnergy();
      if ( ionisation_birks_constant > 0e0 )
Markus Frank
committed
        str << "  Birk's constant: " << std::setprecision(4) << ionisation->GetBirksConstant() << " [mm/MeV]";
      if ( ionisation_ene_per_ion_pair > 0e0 )
Markus Frank
committed
        str << "  Mean Energy Per Ion Pair: " << std::setprecision(4) << ionisation->GetMeanEnergyPerIonPair()/CLHEP::eV << " [eV]";
    }
    else  {
      str << "          No ionisation parameters availible.";
    }
Markus Frank
committed
    printout(lvl, name, "++ Created G4 material %s", str.str().c_str());
    if ( !plugin_name.empty() )    {
      // Call plugin to create extended material if requested
      Detector* det = const_cast<Detector*>(&m_detDesc);
      G4Material* extended_mat = PluginService::Create<G4Material*>(plugin_name, det, medium, mat);
      if ( !extended_mat )   {
Markus Frank
committed
        except("G4Cnv::material["+name+"]","++ FATAL Failed to call plugin to create material.");
      }
      mat = extended_mat;
    }
    info.g4Materials[medium] = mat;
void* Geant4Converter::handleSolid(const string& name, const TGeoShape* shape) const {
  G4VSolid* solid = nullptr;
  if ( shape ) {
    if ( nullptr != (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);
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
    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);
    else if (isa == TGeoPara::Class())
      solid = convertShape<TGeoPara>(shape);
#if ROOT_VERSION_CODE > ROOT_VERSION(6,21,0)
    else if (isa == TGeoTessellated::Class()) 
      solid = convertShape<TGeoTessellated>(shape);
    else if (isa == TGeoScaledShape::Class())  {
      TGeoScaledShape* sh   = (TGeoScaledShape*) shape;
      if ( sol->IsA() == TGeoShapeAssembly::Class() )  {
Markus Frank
committed
        return solid;
      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 )
Markus Frank
committed
        solid = new G4ScaledSolid(sh->GetName(), g4solid, scal);
Markus Frank
committed
        solid = new G4ReflectedSolid(g4solid->GetName()+"_refl", g4solid, scal);
    }
    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());
      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() ) {
        G4Transform3D transform;
Markus Frank
committed
        g4Transform(matrix, transform);
        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 {
  Volume _v(volume);
  Geant4GeometryInfo& info = data();
  PrintLevel lvl = debugVolumes ? ALWAYS : outputLevel;
Markus Frank
committed
  Geant4GeometryMaps::VolumeMap::const_iterator volIt = info.g4Volumes.find(volume);
  if ( _v.testFlagBit(Volume::VETO_SIMU) )  {
    printout(lvl, "Geant4Converter", "++ Volume %s not converted [Veto'ed for simulation]",volume->GetName());
    return nullptr;
  }
  else if (volIt == info.g4Volumes.end() ) {
    const char*  vnam = volume->GetName();
    TGeoMedium*  med  = volume->GetMedium();
    Solid        sh   = volume->GetShape();
    bool         is_assembly = sh->IsA() == TGeoShapeAssembly::Class() || volume->IsAssembly();
    printout(lvl, "Geant4Converter", "++ Convert Volume %-32s: %p %s/%s assembly:%s",
             vnam, volume, sh.type(), _v.type(), yes_no(is_assembly));
    if ( is_assembly ) {
      return nullptr;
    }
    Region        reg      = _v.region();
    LimitSet      lim      = _v.limitSet();
    VisAttr       vis      = _v.visAttributes();
    G4Region*     g4region = reg.isValid() ? info.g4Regions[reg] : nullptr;
    G4UserLimits* g4limits = lim.isValid() ? info.g4Limits[lim]  : nullptr;
    G4VSolid*     g4solid  = (G4VSolid*)   handleSolid(sh->GetName(), sh);
    G4Material*   g4medium = (G4Material*) handleMaterial(med->GetName(), Material(med));
    /// Check all pre-conditions
    if ( !g4solid )   {
      except("G4Converter","++ No Geant4 Solid present for volume: %s", vnam);
Markus Frank
committed
    }
    else if ( !g4medium )   {
      except("G4Converter","++ No Geant4 material present for volume: %s", vnam);
Markus Frank
committed
    }
    else if ( reg.isValid() && !g4region )  {
      except("G4Cnv::volume["+name+"]"," ++ Failed to access Geant4 region %s.", reg.name());
Markus Frank
committed
    }
    else if ( lim.isValid() && !g4limits )  {
      except("G4Cnv::volume["+name+"]","++ FATAL Failed to access Geant4 user limits %s.", lim.name());
    else if ( g4limits )   {
      printout(lvl, "Geant4Converter", "++ Volume     + Apply LIMITS settings: %-24s to volume %s.",
               lim.name(), vnam);
Markus Frank
committed
    }
    G4LogicalVolume* g4vol = nullptr;
    if ( _v.hasProperties() && !_v.getProperty(GEANT4_TAG_PLUGIN,"").empty() )   {
      Detector* det = const_cast<Detector*>(&m_detDesc); 
      string plugin = _v.getProperty(GEANT4_TAG_PLUGIN,"");
      g4vol = PluginService::Create<G4LogicalVolume*>(plugin, det, _v, g4solid, g4medium);
      if ( !g4vol )    {
Markus Frank
committed
        except("G4Cnv::volume["+name+"]","++ FATAL Failed to call plugin to create logical volume.");
      }
    }
    else  {
      g4vol = new G4LogicalVolume(g4solid, g4medium, vnam, nullptr, nullptr, nullptr);
    }
    if ( g4limits )   {
      g4vol->SetUserLimits(g4limits);
    }
    if ( g4region )   {
      PrintLevel plevel = (debugVolumes||debugRegions||debugLimits) ? ALWAYS : outputLevel;
      printout(plevel, "Geant4Converter", "++ Volume     + Apply REGION settings: %-24s to volume %s.",
Markus Frank
committed
               reg.name(), vnam);
      // Handle the region settings for the world volume seperately.
      // Geant4 does NOT WANT any regions assigned to the workd volume.
      // The world's region is created in the G4RunManagerKernel!
      if ( _v == m_detDesc.worldVolume() )   {
Markus Frank
committed
        const char* wrd_nam = "DefaultRegionForTheWorld";
        const char* src_nam = g4region->GetName().c_str();
        auto* world_region  = G4RegionStore::GetInstance()->GetRegion(wrd_nam, false);
        if ( auto* cuts = g4region->GetProductionCuts() )   {
          world_region->SetProductionCuts(cuts);
          printout(plevel, "Geant4Converter",
                   "++ Volume %s Region: %s. Apply production cuts from %s", 
                   vnam, wrd_nam, src_nam);
        }
        if ( auto* lims = g4region->GetUserLimits() )   {
          world_region->SetUserLimits(lims);
          printout(plevel, "Geant4Converter",
                   "++ Volume %s Region: %s. Apply user limits from %s", 
                   vnam, wrd_nam, src_nam);
        }
Markus Frank
committed
        g4vol->SetRegion(g4region);
        g4region->AddRootLogicalVolume(g4vol);
    G4VisAttributes* g4vattr = vis.isValid()
      ? (G4VisAttributes*)handleVis(vis.name(), vis) : nullptr;
    if ( g4vattr )   {
      g4vol->SetVisAttributes(g4vattr);
Markus Frank
committed
    }
    info.g4Volumes[volume] = g4vol;
    printout(lvl, "Geant4Converter",
Markus Frank
committed
             "++ Volume     + %s converted: %p ---> G4: %p", vnam, 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();
  bool              world = (volume == m_detDesc.worldVolume().ptr());
  if ( !world )   {
    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() )    {
    return nullptr;
  Volume _v(mot_vol);
  if ( _v.testFlagBit(Volume::VETO_SIMU) )  {
    printout(lvl, "Geant4Converter", "++ AssemblyNode %s not converted [Veto'ed for simulation]",node->GetName());
    return nullptr;
  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();
      G4Transform3D transform;
      g4Transform(tr, transform);
      if ( is_left_handed(tr) )   {
Markus Frank
committed
        G4Scale3D     scale;
        G4Rotate3D    rot;
        G4Translate3D trans;
        transform.getDecomposition(scale, rot, trans);
        printout(debugReflections ? ALWAYS : lvl, "Geant4Converter",
                 "++ Placing reflected ASSEMBLY. dau:%s to mother %s "
                 "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(),
Markus Frank
committed
                 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());
          return nullptr;
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 nullptr;
    TGeoVolume* mot_vol = node->GetMotherVolume();
    TGeoMatrix* tr = node->GetMatrix();
      except("Geant4Converter",
Markus Frank
committed
             "++ Attempt to handle placement without transformation:%p %s of type %s vol:%p",
             node, node->GetName(), node->IsA()->GetName(), vol);
    else if (nullptr == vol) {
      except("Geant4Converter", "++ Unknown G4 volume:%p %s of type %s ptr:%p",
Markus Frank
committed
             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;
      G4Transform3D transform;
Markus Frank
committed
      Geant4GeometryMaps::VolumeMap::const_iterator volIt = info.g4Volumes.find(mot_vol);
      g4Transform(tr, transform);
      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)" : "",