TMM : Transfer Matrix Method
multi-layer thin films, coherent calc:
header-only GPU/CPU : C4MultiLayrStack.h
https://github.com/simoncblyth/customgeant4/
Natural 2 volume PMT geometry : HAMA and NNVT
Custom Boundary POM
+---------------pmt-Pyrex----------------+ | | | | | | | +~inner~Vacuum~~~~~~~~~~~+ | | ! ! | | ! ! | | ! ! | | ! ! | | ! ! | | + + | | | | | | | | | | | | | | | | | | | | | | +------------------------+ | | | | | | | +----------------------------------------+
OpSurfaceName[0] == '@'
Custom Boundary Process : Advantages
Old FastSim POM | 4 Solid, 4 LV, 4 PV |
Custom Boundary POM | 2 Solid, 2 LV, 2 PV |
Disadvantages
https://github.com/simoncblyth/customgeant4/
PostStepDoIt type switch
if( m_custom_status == 'Y' ) // CustomART handling { G4double rand = G4UniformRand(); if ( rand < theAbsorption ) { DoAbsorption(); // A } else { DielectricDielectric(); // R or T } } else if (type == dielectric_metal) { DielectricMetal(); } ...
Standard DoAbsorption DielectricDielectric reused
G4OpBoundaryProcess unless OpticalSurfaceNam[0] == '@'/'#'
if( OpticalSurfaceName0 == '@' || OpticalSurfaceName0 == '#' ) { if( m_custom_art->local_z(aTrack) < 0. ) // Lower hemi : Standard { m_custom_status = 'Z' ; } else if( OpticalSurfaceName0 == '@') // MultiFilm ART POM { m_custom_status = 'Y' ; m_custom_art->doIt(aTrack, aStep) ; type = dielectric_dielectric ; theModel = glisur ; theFinish = polished ; } else if( OpticalSurfaceName0 == '#' ) // Traditional POM { m_custom_status = '-' ; type = dielectric_metal ; theModel = glisur ; theReflectivity = 0. ; theTransmittance = 0. ; theEfficiency = 1. ; } }
DielectricDielectric expects 2-way (theReflectivity,theTransmittance) => so C4CustomART.h rescales 3-way (A,R,T)
(theAbsorption,theReflectivity,theTransmittance) <= (A, R/(1-A), T/(1-A))
C4CustomART customizations
doIt C4MultiLayrStack.h TMM calc, only for:
reuse standard G4OpBoundaryProcess
#include "C4IPMTAccessor.h" #include "C4MultiLayrStack.h" #include "C4Touchable.h" struct C4CustomART { const C4IPMTAccessor* accessor ; G4double& theAbsorption ; // doIt sets these G4double& theReflectivity ; G4double& theTransmittance ; G4double& theEfficiency ; const G4ThreeVector& theGlobalPoint ; const G4ThreeVector& OldMomentum ; const G4ThreeVector& OldPolarization ; const G4ThreeVector& theRecoveredNormal ; const G4double& thePhotonMomentum ; C4CustomART( const C4IPMTAccessor* accessor, G4double& theAbsorption, G4double& theReflectivity, G4double& theTransmittance, G4double& theEfficiency, const G4ThreeVector& theGlobalPoint, const G4ThreeVector& OldMomentum, const G4ThreeVector& OldPolarization, const G4ThreeVector& theRecoveredNormal, const G4double& thePhotonMomentum ); double local_z( const G4Track& aTrack ); void doIt(const G4Track& aTrack, const G4Step& aStep ); };
https://github.com/simoncblyth/customgeant4/blob/main/C4CustomART.h
Why Low Dependency Access ?
QE scan over 100 energy points, all PMTs
IPMTSimParamSvc | 7.70 s |
PMTSimParamData | 0.30 s |
IPMTSimParamSvc/PMTSimParamData | 26.10 |
MR 126 MERGED Feb 3, 2023 |
Full PMT access from Custom4 + Opticks without Svc : How ?
SimSvc/PMTSimParamSvc/PMTSimParamSvc/PMTAccessor.h
#include "C4IPMTAccessor.h" struct PMTAccessor : public C4IPMTAccessor { PMTAccessor(const PMTSimParamData* data); int get_num_lpmt() const ; double get_pmtid_qe( int pmtid, double energy ) const ; ...
Pure virtual protocol base : Custom4/C4IPMTAccessor.h
struct C4IPMTAccessor // only standard types in API { virtual int get_num_lpmt() const = 0 ; virtual double get_pmtid_qe( int pmtid, double energy ) const = 0 ; virtual double get_qescale( int pmtid ) const = 0 ; ...
Consistent PMT data access : JUNOSW + Opticks + Custom4
qsim::propagate_at_surface_CustomART
int qsim::propagate_at_surface_CustomART( unsigned& flag, curandStateXORWOW& rng, sctx& ctx) const { const sphoton& p = ctx.p ; const float3* nrm = (float3*)&ctx.prd->q0.f.x ; int lpmtid = ctx.prd->identity(); //0->17611 float minus_cos_theta = dot(p.mom, *nrm); float dot_pol_x_mom_nrm = dot(p.pol,cross(p.mom,*nrm)); // enables calc of S/P-polarization fraction float ARTE[4] ; pmt->get_lpmtid_ARTE(ARTE, lpmtid, p.wavelength, minus_cos_theta, dot_pol_x_mom_nrm ); const float& A = ARTE[0]; // Absorb const float& T = ARTE[2]; // Transmit (R+T=1) const float& E = ARTE[3]; // Effic: _qe/A_norm float u_A = curand_uniform(&rng); int action = u_A < A ? BREAK : CONTINUE ; if( action == BREAK ){ float u_E = curand_uniform(&rng) ; flag = u_E < E ? SURFACE_DETECT : SURFACE_ABSORB ; } else { propagate_at_boundary( flag, rng, ctx, T ); } return action ; }
https://bitbucket.org/simoncblyth/opticks/src/master/qudarap/qpmt.h
Shapes
In [6]: h.f.art.shape Out[6]: (9, 900, 4, 4) In [7]: h.f.comp.shape Out[7]: (9, 900, 1, 4, 4, 2)
CFLayrTest a : LayrTest : scan__R12860__cpu_thr_double b : LayrTest : scan__R12860__cpu_thr_float c : LayrTest : scan__R12860__gpu_thr_double d : LayrTest : scan__R12860__gpu_thr_float e : SPMT_test : sscan f : QPMTTest : qscan g : QPMT_MockTest : qscan h : QPMT_Test : qscan ## tab, rst = ts.cf_table(tt, pmtcat, excl=excl) # excl 0.05 ## rst +------------------------------+----------+----------+----------+----------+----------+----------+----------+----------+ | R12860 art\comp 0.05| a:ctd| b:ctf| c:gtd| d:gtf| e:| f:| g:| h:| +==============================+==========+==========+==========+==========+==========+==========+==========+==========+ | a:ctd| 0| 4.829e-05| 1.066e-14| 4.829e-05| 8.644e-05| 3.422e-05| 4.255e-05| 3.685e-05| +------------------------------+----------+----------+----------+----------+----------+----------+----------+----------+ | b:ctf| 9.317e-07| 0| 4.829e-05| 5.722e-06| 4.578e-05| 5.15e-05| 5.15e-05| 5.531e-05| +------------------------------+----------+----------+----------+----------+----------+----------+----------+----------+ | c:gtd| 1.582e-15| 9.317e-07| 0| 4.829e-05| 8.644e-05| 3.422e-05| 4.255e-05| 3.685e-05| +------------------------------+----------+----------+----------+----------+----------+----------+----------+----------+ | d:gtf| 7.958e-07| 8.792e-07| 7.958e-07| 0| 4.196e-05| 5.341e-05| 5.722e-05| 5.913e-05| +------------------------------+----------+----------+----------+----------+----------+----------+----------+----------+ | e:| 2.956e-06| 3.159e-06| 2.956e-06| 3.07e-06| 0| 5.341e-05| 5.15e-05| 4.959e-05| +------------------------------+----------+----------+----------+----------+----------+----------+----------+----------+ | f:| 3.094e-06| 3.338e-06| 3.094e-06| 3.159e-06| 1.192e-06| 0| 1.335e-05| 9.537e-06| +------------------------------+----------+----------+----------+----------+----------+----------+----------+----------+ | g:| 2.956e-06| 3.159e-06| 2.956e-06| 3.07e-06| 7.451e-07| 1.103e-06| 0| 7.629e-06| +------------------------------+----------+----------+----------+----------+----------+----------+----------+----------+ | h:| 2.907e-06| 3.248e-06| 2.907e-06| 2.921e-06| 7.451e-07| 8.941e-07| 8.419e-07| 0| +------------------------------+----------+----------+----------+----------+----------+----------+----------+----------+
Largest differences between ART (9, 900, 4, 4) from various tests (float,double,GPU,CPU,mock-CPU,standalone...)