File indexing completed on 2026-04-10 07:50:34
0001 #pragma once
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014 #include "G4String.hh"
0015 #include "G4OpticalSurface.hh"
0016 #include "G4MaterialPropertiesTable.hh"
0017 #include "G4LogicalBorderSurface.hh"
0018 #include "G4LogicalSkinSurface.hh"
0019 #include "G4Version.hh"
0020 #include "G4Material.hh"
0021
0022 #include "U4SurfaceType.h"
0023 #include "U4OpticalSurfaceFinish.h"
0024 #include "U4OpticalSurfaceModel.h"
0025
0026 #include "U4SurfacePerfect.h"
0027 #include "S4.h"
0028
0029 #include "sdomain.h"
0030 #include "sproplist.h"
0031
0032 #if G4VERSION_NUMBER >= 1070
0033 #include "SNameOrder.h"
0034 #endif
0035
0036
0037
0038 struct NPFold ;
0039
0040
0041 enum {
0042 U4Surface_UNSET,
0043 U4Surface_PerfectAbsorber,
0044 U4Surface_PerfectDetector
0045 };
0046
0047 struct U4Surface
0048 {
0049 static constexpr const char* PerfectAbsorber = "PerfectAbsorber" ;
0050 static constexpr const char* PerfectDetector = "PerfectDetector" ;
0051 static unsigned Type(const char* type_);
0052
0053 static G4OpticalSurface* MakeOpticalSurface( const char* name_ );
0054
0055 static G4LogicalBorderSurface* MakeBorderSurface(const char* name_, const char* type_, const char* pv1_, const char* pv2_, const G4VPhysicalVolume* start_pv );
0056 static G4LogicalBorderSurface* MakePerfectAbsorberBorderSurface(const char* name_, const char* pv1, const char* pv2, const G4VPhysicalVolume* start_pv );
0057 static G4LogicalBorderSurface* MakePerfectDetectorBorderSurface(const char* name_, const char* pv1, const char* pv2, const G4VPhysicalVolume* start_pv );
0058
0059 static G4LogicalBorderSurface* MakeBorderSurface(const char* name_, const char* type_, const G4VPhysicalVolume* pv1, const G4VPhysicalVolume* pv2 );
0060 static G4LogicalBorderSurface* MakePerfectAbsorberBorderSurface(const char* name_, const G4VPhysicalVolume* pv1, const G4VPhysicalVolume* pv2 );
0061 static G4LogicalBorderSurface* MakePerfectDetectorBorderSurface(const char* name_, const G4VPhysicalVolume* pv1, const G4VPhysicalVolume* pv2 );
0062
0063
0064 static const std::vector<G4LogicalBorderSurface*>* PrepareBorderSurfaceVector(const G4LogicalBorderSurfaceTable* tab );
0065 static const std::vector<G4LogicalSkinSurface*>* PrepareSkinSurfaceVector(const G4LogicalSkinSurfaceTable* tab );
0066
0067 static void Collect( std::vector<const G4LogicalSurface*>& surfaces );
0068 static void CollectRawNames( std::vector<std::string>& rawnames, const std::vector<const G4LogicalSurface*>& surfaces );
0069
0070 static NPFold* MakeFold(const std::vector<const G4LogicalSurface*>& surfaces, const std::vector<std::string>& keys );
0071 static NPFold* MakeFold();
0072 static G4LogicalSurface* Find( const G4VPhysicalVolume* thePrePV, const G4VPhysicalVolume* thePostPV ) ;
0073
0074 };
0075
0076
0077 #include "U4Material.hh"
0078 #include "U4MaterialPropertiesTable.h"
0079 #include "U4Volume.h"
0080 #include "NPFold.h"
0081
0082
0083 inline unsigned U4Surface::Type(const char* type_)
0084 {
0085 unsigned type = U4Surface_UNSET ;
0086 if(strcmp(type_, PerfectAbsorber) == 0) type = U4Surface_PerfectAbsorber ;
0087 if(strcmp(type_, PerfectDetector) == 0) type = U4Surface_PerfectDetector ;
0088 return type ;
0089 }
0090
0091
0092
0093 inline G4OpticalSurface* U4Surface::MakeOpticalSurface( const char* name_ )
0094 {
0095 G4String name = name_ ;
0096 G4OpticalSurfaceModel model = glisur ;
0097 G4OpticalSurfaceFinish finish = polished ;
0098 G4SurfaceType type = dielectric_dielectric ;
0099 G4double value = 1.0 ;
0100 G4OpticalSurface* os = new G4OpticalSurface(name, model, finish, type, value );
0101 return os ;
0102 }
0103
0104
0105
0106
0107
0108
0109
0110
0111
0112
0113
0114
0115
0116
0117
0118
0119 inline G4LogicalBorderSurface* U4Surface::MakeBorderSurface(const char* name_, const char* type_, const char* pv1_, const char* pv2_, const G4VPhysicalVolume* start_pv )
0120 {
0121 const G4VPhysicalVolume* pv1 = U4Volume::FindPV( start_pv, pv1_ );
0122 const G4VPhysicalVolume* pv2 = U4Volume::FindPV( start_pv, pv2_ );
0123 return ( pv1 && pv2 ) ? MakeBorderSurface(name_, type_, pv1, pv2 ) : nullptr ;
0124 }
0125
0126 inline G4LogicalBorderSurface* U4Surface::MakePerfectAbsorberBorderSurface(const char* name_, const char* pv1, const char* pv2, const G4VPhysicalVolume* start_pv)
0127 {
0128 return MakeBorderSurface(name_, PerfectAbsorber, pv1, pv2, start_pv );
0129 }
0130 inline G4LogicalBorderSurface* U4Surface::MakePerfectDetectorBorderSurface(const char* name_, const char* pv1, const char* pv2, const G4VPhysicalVolume* start_pv)
0131 {
0132 return MakeBorderSurface(name_, PerfectDetector, pv1, pv2, start_pv );
0133 }
0134
0135
0136
0137
0138 inline G4LogicalBorderSurface* U4Surface::MakeBorderSurface(const char* name_, const char* type_, const G4VPhysicalVolume* pv1, const G4VPhysicalVolume* pv2)
0139 {
0140 unsigned type = Type(type_);
0141
0142 G4OpticalSurface* os = MakeOpticalSurface( name_ );
0143 G4MaterialPropertiesTable* mpt = new G4MaterialPropertiesTable ;
0144 os->SetMaterialPropertiesTable(mpt);
0145
0146 G4MaterialPropertyVector* rindex = U4Material::MakeProperty(1.);
0147 mpt->AddProperty("RINDEX", rindex );
0148
0149 G4MaterialPropertyVector* reflectivity = U4Material::MakeProperty(0.);
0150 mpt->AddProperty("REFLECTIVITY",reflectivity );
0151
0152
0153 if( type == U4Surface_PerfectAbsorber )
0154 {
0155 }
0156 else if( type == U4Surface_PerfectDetector )
0157 {
0158 G4MaterialPropertyVector* efficiency = U4Material::MakeProperty(1.);
0159 mpt->AddProperty("EFFICIENCY",efficiency );
0160 }
0161
0162 G4String name = name_ ;
0163
0164 G4VPhysicalVolume* pv1_ = const_cast<G4VPhysicalVolume*>(pv1);
0165 G4VPhysicalVolume* pv2_ = const_cast<G4VPhysicalVolume*>(pv2);
0166 G4LogicalBorderSurface* bs = new G4LogicalBorderSurface(name, pv1_, pv2_, os );
0167 return bs ;
0168 }
0169
0170 inline G4LogicalBorderSurface* U4Surface::MakePerfectAbsorberBorderSurface(const char* name_, const G4VPhysicalVolume* pv1, const G4VPhysicalVolume* pv2)
0171 {
0172 return MakeBorderSurface(name_, PerfectAbsorber, pv1, pv2 );
0173 }
0174 inline G4LogicalBorderSurface* U4Surface::MakePerfectDetectorBorderSurface(const char* name_, const G4VPhysicalVolume* pv1, const G4VPhysicalVolume* pv2)
0175 {
0176 return MakeBorderSurface(name_, PerfectDetector, pv1, pv2 );
0177 }
0178
0179
0180
0181
0182
0183
0184
0185
0186
0187
0188
0189
0190
0191
0192
0193
0194
0195
0196
0197
0198
0199
0200
0201
0202
0203
0204
0205
0206
0207
0208
0209
0210 inline const std::vector<G4LogicalBorderSurface*>* U4Surface::PrepareBorderSurfaceVector(const G4LogicalBorderSurfaceTable* tab )
0211 {
0212 typedef std::vector<G4LogicalBorderSurface*> VBS ;
0213 #if G4VERSION_NUMBER >= 1070
0214 typedef std::pair<const G4VPhysicalVolume*, const G4VPhysicalVolume*> PPV ;
0215 typedef std::map<PPV, G4LogicalBorderSurface*>::const_iterator IT ;
0216
0217 VBS* vec = new VBS ;
0218 for(IT it=tab->begin() ; it != tab->end() ; it++ )
0219 {
0220 G4LogicalBorderSurface* bs = it->second ;
0221 vec->push_back(bs);
0222 const PPV ppv = it->first ;
0223 assert( ppv.first == bs->GetVolume1());
0224 assert( ppv.second == bs->GetVolume2());
0225 }
0226
0227 {
0228 bool reverse = false ;
0229 const char* tail = "0x" ;
0230 SNameOrder<G4LogicalBorderSurface>::Sort( *vec, reverse, tail );
0231 std::cout << "U4Surface::PrepareBorderSurfaceVector\n" << SNameOrder<G4LogicalBorderSurface>::Desc( *vec ) << std::endl ;
0232 }
0233
0234 #else
0235 const VBS* vec = tab ;
0236 #endif
0237 return vec ;
0238 }
0239
0240
0241
0242
0243
0244
0245
0246
0247
0248
0249
0250 inline const std::vector<G4LogicalSkinSurface*>* U4Surface::PrepareSkinSurfaceVector(const G4LogicalSkinSurfaceTable* tab )
0251 {
0252 typedef std::vector<G4LogicalSkinSurface*> VKS ;
0253 #if G4VERSION_NUMBER >= 1122
0254 typedef std::map<const G4LogicalVolume*,G4LogicalSkinSurface*>::const_iterator IT ;
0255 VKS* vec = new VKS ;
0256
0257 for(IT it=tab->begin() ; it != tab->end() ; it++ )
0258 {
0259 G4LogicalSkinSurface* ks = it->second ;
0260 vec->push_back(ks);
0261 }
0262
0263 {
0264 bool reverse = false ;
0265 const char* tail = "0x" ;
0266 SNameOrder<G4LogicalSkinSurface>::Sort( *vec, reverse, tail );
0267 std::cout << "U4Surface::PrepareSkinSurfaceVector\n" << SNameOrder<G4LogicalSkinSurface>::Desc( *vec ) << std::endl ;
0268 }
0269
0270 #else
0271 const VKS* vec = tab ;
0272 #endif
0273 return vec ;
0274 }
0275
0276
0277
0278
0279
0280
0281
0282
0283
0284
0285
0286
0287 inline void U4Surface::Collect( std::vector<const G4LogicalSurface*>& surfaces )
0288 {
0289 const G4LogicalBorderSurfaceTable* border_ = G4LogicalBorderSurface::GetSurfaceTable() ;
0290 const std::vector<G4LogicalBorderSurface*>* border = PrepareBorderSurfaceVector(border_);
0291
0292 for(unsigned i=0 ; i < border->size() ; i++)
0293 {
0294 G4LogicalBorderSurface* bs = (*border)[i] ;
0295 surfaces.push_back(bs) ;
0296 }
0297
0298 const G4LogicalSkinSurfaceTable* skin_ = G4LogicalSkinSurface::GetSurfaceTable() ;
0299 const std::vector<G4LogicalSkinSurface*>* skin = PrepareSkinSurfaceVector(skin_);
0300 for(unsigned i=0 ; i < skin->size() ; i++)
0301 {
0302 G4LogicalSkinSurface* ks = (*skin)[i] ;
0303 surfaces.push_back(ks) ;
0304 }
0305 }
0306
0307 inline void U4Surface::CollectRawNames( std::vector<std::string>& rawnames, const std::vector<const G4LogicalSurface*>& surfaces )
0308 {
0309 for(unsigned i=0 ; i < surfaces.size() ; i++)
0310 {
0311 const G4LogicalSurface* ls = surfaces[i] ;
0312 const G4String& name = ls->GetName() ;
0313 const char* raw = name.c_str() ;
0314 rawnames.push_back(raw);
0315 }
0316 }
0317
0318
0319
0320
0321
0322
0323
0324
0325
0326
0327
0328 inline NPFold* U4Surface::MakeFold(const std::vector<const G4LogicalSurface*>& surfaces, const std::vector<std::string>& keys )
0329 {
0330 assert( surfaces.size() == keys.size() );
0331
0332 NPFold* fold = new NPFold ;
0333 for(unsigned i=0 ; i < surfaces.size() ; i++)
0334 {
0335 const G4LogicalSurface* ls = surfaces[i] ;
0336 [[maybe_unused]] const char* rawname = ls->GetName().c_str() ;
0337 const char* key = keys[i].c_str() ;
0338
0339 G4OpticalSurface* os = dynamic_cast<G4OpticalSurface*>(ls->GetSurfaceProperty());
0340
0341 G4SurfaceType theType = os->GetType();
0342 G4OpticalSurfaceModel theModel = os->GetModel();
0343 G4OpticalSurfaceFinish theFinish = os->GetFinish();
0344
0345
0346 G4double ModelValue = theModel == glisur ? os->GetPolish() : os->GetSigmaAlpha() ;
0347 assert( ModelValue >= 0. && ModelValue <= 1. );
0348
0349 const char* osn = os->GetName().c_str() ;
0350 G4MaterialPropertiesTable* mpt = os->GetMaterialPropertiesTable() ;
0351
0352 assert(mpt);
0353 NPFold* sub = U4MaterialPropertiesTable::MakeFold(mpt) ;
0354
0355 sub->set_meta<std::string>("rawname", rawname) ;
0356 sub->set_meta<std::string>("OpticalSurfaceName", osn) ;
0357 sub->set_meta<std::string>("TypeName", U4SurfaceType::Name(theType)) ;
0358 sub->set_meta<std::string>("ModelName", U4OpticalSurfaceModel::Name(theModel)) ;
0359 sub->set_meta<std::string>("FinishName", U4OpticalSurfaceFinish::Name(theFinish)) ;
0360
0361 sub->set_meta<int>("Type", theType) ;
0362 sub->set_meta<int>("Model", theModel) ;
0363 sub->set_meta<int>("Finish", theFinish) ;
0364 sub->set_meta<double>("ModelValue", ModelValue ) ;
0365
0366
0367 const G4LogicalBorderSurface* bs = dynamic_cast<const G4LogicalBorderSurface*>(ls) ;
0368 const G4LogicalSkinSurface* ks = dynamic_cast<const G4LogicalSkinSurface*>(ls) ;
0369
0370 if(bs)
0371 {
0372 const G4VPhysicalVolume* _pv1 = bs->GetVolume1();
0373 const G4VPhysicalVolume* _pv2 = bs->GetVolume2();
0374
0375 const char* pv1 = S4::Name<G4VPhysicalVolume>(_pv1) ;
0376 const char* pv2 = S4::Name<G4VPhysicalVolume>(_pv2) ;
0377
0378 sub->set_meta<std::string>("pv1", pv1) ;
0379 sub->set_meta<std::string>("pv2", pv2) ;
0380 sub->set_meta<std::string>("type", "Border" );
0381 }
0382 else if(ks)
0383 {
0384 const G4LogicalVolume* _lv = ks->GetLogicalVolume();
0385 const char* lv = S4::Name<G4LogicalVolume>(_lv);
0386
0387 sub->set_meta<std::string>("lv", lv );
0388 sub->set_meta<std::string>("type", "Skin" );
0389 }
0390 fold->add_subfold( key, sub );
0391 }
0392 return fold ;
0393 }
0394
0395 inline NPFold* U4Surface::MakeFold()
0396 {
0397
0398 std::vector<const G4LogicalSurface*> surfaces ;
0399 Collect(surfaces);
0400
0401 std::vector<std::string> suname_raw ;
0402 U4Surface::CollectRawNames(suname_raw, surfaces);
0403
0404 std::vector<std::string> suname ;
0405 sstr::StripTail_Unique( suname, suname_raw, "0x" );
0406
0407 return MakeFold(surfaces, suname) ;
0408 }
0409
0410
0411
0412
0413
0414
0415
0416
0417
0418
0419
0420
0421
0422 inline G4LogicalSurface* U4Surface::Find( const G4VPhysicalVolume* thePrePV, const G4VPhysicalVolume* thePostPV )
0423 {
0424 if(thePostPV == nullptr || thePrePV == nullptr ) return nullptr ;
0425 G4LogicalSurface* Surface = G4LogicalBorderSurface::GetSurface(thePrePV, thePostPV);
0426 if(Surface == nullptr)
0427 {
0428 G4bool enteredDaughter = thePostPV->GetMotherLogical() == thePrePV->GetLogicalVolume();
0429 if(enteredDaughter)
0430 {
0431 Surface = G4LogicalSkinSurface::GetSurface(thePostPV->GetLogicalVolume());
0432 if(Surface == nullptr)
0433 Surface = G4LogicalSkinSurface::GetSurface(thePrePV->GetLogicalVolume());
0434 }
0435 else
0436 {
0437 Surface = G4LogicalSkinSurface::GetSurface(thePrePV->GetLogicalVolume());
0438 if(Surface == nullptr)
0439 Surface = G4LogicalSkinSurface::GetSurface(thePostPV->GetLogicalVolume());
0440 }
0441 }
0442 return Surface ;
0443 }
0444
0445
0446