Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-04-09 07:49:34

0001 #pragma once
0002 /**
0003 sevent.h : host/device communication instance : holder of array pointers
0004 ===========================================================================
0005 
0006 Instantiation of sevent is done by SEvt::SEvt
0007 and the instance is subsequently uploaded to the device after
0008 device buffer allocations hence the *sevent* instance
0009 provides event config and device buffer pointers
0010 both on device and host.
0011 
0012 Note that *num_seed* and *num_photon* will be equal in
0013 normal operation which uses QEvt::setGensteps.
0014 However for clarity separate fields are used to
0015 distinguish photon test running that directly uses
0016 QEvt::setNumPhoton
0017 
0018 QEvt::getDevicePtr returns an instance of sevent.h : as such
0019 sevent.h can be regarded as the GPU side summary of QEvt
0020 related info.
0021 
0022 
0023 Users of sevent.h
0024 -------------------
0025 
0026 qudarap/QU.cc
0027     template instanciations of QU::UploadArray QU::device_alloc QU::copy_host_to_device
0028 
0029 CSGOptiX/CSGOptiX7.cu
0030     simulate : reads params.evt, sets sctx.evt
0031     reads evt.seed evt.genstep evt.max_bounce
0032     writes evt.photon
0033 
0034     simtrace : reads params.evt, evt.num_simtrace, calls evt::add_simtrace
0035 
0036 sysrap/sctx.h
0037     *evt* member of device context *sctx*
0038     sctx::point writes evt.record evt.rec
0039     sctx::trace writes evt.prd
0040     sctx::end writes evt.seq evt.tag evt.flat
0041 
0042 sysrap/SEvt.hh sysrap/SEvt.cc
0043     SEvt::SEvt instanciates hostside sevent evt member
0044     SEvt::init/sevent::init sets maxima
0045 
0046     SEvt::addGenstep/SEvt::setNumPhoton
0047     sets evt.num_photon/seq/tag/flat/record/rec/prd depending on the evt.max
0048 
0049 qudarap/QEvt.hh qudarap/QEvt.cc
0050     **fundamentally integrated with sevent.h**
0051 
0052     QEvt::QEvt
0053 
0054     * gets sevent:evt pointer from SEvt
0055     * allocates space on device for d_evt
0056 
0057     QEvt::setNumPhoton/QEvt::uploadEvt
0058 
0059     * updates device d_evt by copying from evt
0060 
0061 
0062 qudarap/qsim.h
0063     qsim::mock_propagate mocks the CSGOptiX7.cu:simulate
0064     qsim::generate_photon uses evt to read input photons from buffer
0065 
0066 qudarap/QSim.cu
0067     used from several tests for evt.num, evt.seed, evt.genstep
0068 
0069 qudarap/QEvt.cu
0070     sevent.h used for tests of seeding and genstep photon counting
0071 
0072 **/
0073 
0074 #if defined(__CUDACC__) || defined(__CUDABE__)
0075 #    define SEVENT_METHOD __device__ __forceinline__
0076 #else
0077 #    define SEVENT_METHOD inline
0078 #endif
0079 
0080 struct float4 ;
0081 struct float2 ;
0082 struct quad4 ;
0083 struct quad6 ;
0084 
0085 //#if !defined(PRODUCTION)
0086 struct srec ;
0087 struct sseq ;
0088 struct stag ;
0089 struct sflat ;
0090 //#endif
0091 
0092 
0093 struct sphoton ;
0094 struct sphotonlite ;
0095 
0096 #if defined(__CUDACC__) || defined(__CUDABE__)
0097 #else
0098 #include <string>
0099 #include <sstream>
0100 #include <iomanip>
0101 #include "SEventConfig.hh"
0102 #include "NP.hh"
0103 #endif
0104 
0105 struct sevent
0106 {
0107     static constexpr int64_t M = 1000000 ;
0108     static constexpr int64_t genstep_itemsize = 6*4 ;
0109     static constexpr int64_t genstep_numphoton_offset = 3 ;
0110 
0111     static constexpr float w_lo = 60.f ;
0112     static constexpr float w_hi = 820.f ;  // these need to match sdomain.h
0113     static constexpr float w_center = (w_lo+w_hi)/2.f ; // convert wavelength range into center-extent form
0114     static constexpr float w_extent = (w_hi-w_lo)/2.f ;
0115 
0116     float4 center_extent ;
0117     float2 time_domain ;
0118     float2 wavelength_domain ;
0119 
0120     // TODO: make all these below unsigned
0121     // sevent::init sets these max using values from SEventConfig
0122     size_t   max_curand ;
0123     size_t   max_slot   ;
0124     size_t   max_genstep ;  // eg:      100,000
0125     size_t   max_photon  ;  // eg:  100,000,000
0126     size_t   max_simtrace ; // eg: 100,000,000
0127 
0128     int   max_bounce  ; // eg: 0:32  (not including 32)
0129     int   max_record  ; // eg: 10  full step record
0130     int   max_rec     ; // eg: 10  compressed step record
0131     int   max_seq     ; // now restricted to 0 or 1 only
0132     int   max_prd     ; // eg: 16
0133     int   max_tag     ; // 0 or 1 only
0134     int   max_flat    ; // 0 or 1 only
0135     int   max_aux     ; // eg 0, 16, 32 : when greater than zero typically follows max_record
0136     int   max_sup     ;
0137 
0138     int   mode_lite    ; // 0 or 1 (also 2 when debug comparing)
0139 
0140     int   index ;
0141 
0142     //[ counts and pointers, zeroed by sevent::zero
0143     //  only first 4 are always in use, the last 7 are only relevant whilst debugging
0144     //
0145     size_t      num_curand ;
0146     size_t      num_genstep ;
0147     size_t      num_seed ;
0148 
0149     size_t       num_hit ;    // set by QEvt::gatherHit using SU::count_if_sphoton
0150     size_t       num_hitmerged ;
0151     size_t       num_photon ;
0152 
0153     size_t       num_hitlite ;  // set by QEvt::gatherHitLite using SU::count_if_sphotonlite
0154     size_t       num_hitlitemerged ;  // set by QEvt::gatherHitLite using SU::count_if_sphotonlite
0155     size_t       num_photonlite ;
0156 
0157     size_t      num_record ;
0158     size_t      num_rec ;
0159     size_t      num_seq ;
0160     size_t      num_prd ;
0161     size_t      num_tag ;
0162     size_t      num_flat ;
0163     size_t      num_simtrace ;
0164     size_t      num_aux ;
0165     size_t      num_sup ;
0166 
0167     // With QEvt device running the below are pointers to device buffers.
0168     // Most are allocated ONCE ONLY by QEvt::device_alloc_genstep/photon/simtrace
0169     // sized by configured maxima.
0170     // Only the hit buffer is allocated more dynamically depending on num_hit.
0171     // TODO: check for leaking of hit buffers
0172 
0173     quad6*   genstep ;    //QEvt::device_alloc_genstep
0174     int*     seed ;
0175 
0176     sphoton*     photon ;     //QEvt::device_alloc_photon
0177     sphoton*     hit ;        //QEvt::gatherHit_ allocates event by event depending on num_hit
0178     sphoton*     hitmerged ;
0179 
0180     sphotonlite* photonlite ; //QEvt::device_alloc_photon
0181     sphotonlite* hitlite ;
0182     sphotonlite* hitlitemerged ;
0183 
0184     sphoton* record ;
0185     srec*    rec ;
0186     sseq*    seq ;
0187 
0188     quad2*   prd ;
0189     stag*    tag ;
0190     sflat*   flat ;
0191     quad4*   simtrace ;   //QEvt::device_alloc_simtrace
0192     quad4*   aux ;
0193     quad6*   sup ;
0194 
0195     //] counts and pointers
0196 
0197 
0198 #if defined(__CUDACC__) || defined(__CUDABE__)
0199 #else
0200     SEVENT_METHOD void init();
0201     SEVENT_METHOD bool with_photon() const ;
0202     SEVENT_METHOD bool with_photonlite() const ;
0203     SEVENT_METHOD bool no_photon_or_photonlite_alloc() const ;
0204 
0205     SEVENT_METHOD void init_domain(float extent, float time_max);
0206 
0207     SEVENT_METHOD std::string descMax() const ;
0208     SEVENT_METHOD std::string descNum() const ;
0209     SEVENT_METHOD std::string descBuf() const ;
0210     SEVENT_METHOD std::string desc() const ;
0211 
0212     SEVENT_METHOD void get_domain(quad4& dom) const ;
0213     SEVENT_METHOD void get_config(quad4& cfg) const ;
0214     SEVENT_METHOD void get_meta(std::string& meta) const ;
0215 
0216     SEVENT_METHOD void zero();
0217 
0218     template<typename T> SEVENT_METHOD T* get_photon_ptr() const ;
0219     template<typename T> SEVENT_METHOD size_t get_photon_num() const ;
0220 
0221     template<typename T> SEVENT_METHOD T**     get_hitmerged_ptr_ref(); // returning refs so non-const
0222     template<typename T> SEVENT_METHOD size_t* get_hitmerged_num_ref();
0223 
0224 
0225 
0226 #endif
0227 
0228 #ifndef PRODUCTION
0229     SEVENT_METHOD void add_rec( srec& r, unsigned idx, unsigned bounce, const sphoton& p);
0230 #endif
0231     SEVENT_METHOD void add_simtrace( unsigned idx, const quad4& p, const quad2* prd, float tmin );
0232 
0233 };
0234 
0235 
0236 #if defined(__CUDACC__) || defined(__CUDABE__)
0237 #else
0238 
0239 SEVENT_METHOD void sevent::init()
0240 {
0241     // int64_t
0242     max_curand   = SEventConfig::MaxCurand() ;
0243     max_slot     = SEventConfig::MaxSlot() ;
0244     max_genstep  = SEventConfig::MaxGenstep() ;
0245     max_photon   = SEventConfig::MaxPhoton()  ;
0246     max_simtrace = SEventConfig::MaxSimtrace()  ;
0247 
0248     // int
0249     max_bounce   = SEventConfig::MaxBounce()  ;
0250     max_record   = SEventConfig::MaxRecord()  ;  // full step record
0251     max_rec      = SEventConfig::MaxRec()  ;     // compressed step record
0252     max_seq      = SEventConfig::MaxSeq()  ;     // seqhis
0253     max_prd      = SEventConfig::MaxPrd()  ;
0254     max_tag      = SEventConfig::MaxTag()  ;
0255     max_flat     = SEventConfig::MaxFlat()  ;
0256     max_aux      = SEventConfig::MaxAux()  ;
0257     max_sup      = SEventConfig::MaxSup()  ;
0258 
0259     mode_lite     = SEventConfig::ModeLite() ;
0260 
0261     zero();  // pointers and counts
0262 
0263     float max_extent_domain = SEventConfig::MaxExtentDomain() ;
0264     float max_time_domain = SEventConfig::MaxTimeDomain() ;
0265 
0266     init_domain( max_extent_domain, max_time_domain );
0267 }
0268 
0269 
0270 // these used from QEvt::device_alloc_photon
0271 SEVENT_METHOD bool sevent::with_photon() const {     return max_slot > 0 && ( mode_lite == 0 || mode_lite == 2 ) ; }
0272 SEVENT_METHOD bool sevent::with_photonlite() const { return max_slot > 0 && ( mode_lite == 1 || mode_lite == 2 ) ; }
0273 SEVENT_METHOD bool sevent::no_photon_or_photonlite_alloc() const { return photon == nullptr && photonlite == nullptr ; }
0274 
0275 
0276 /**
0277 sevent::init_domain
0278 --------------------
0279 
0280 Domains are mainly(only?) used for the compressed photon "rec",
0281 but now are making more use of the full photon "record"
0282 as thats considerably simpler that "rec".
0283 
0284 **/
0285 SEVENT_METHOD void sevent::init_domain(float max_extent_domain, float max_time_domain)
0286 {
0287     center_extent.x = 0.f ;
0288     center_extent.y = 0.f ;
0289     center_extent.z = 0.f ;
0290     center_extent.w = max_extent_domain ;
0291 
0292     time_domain.x = 0.f ;
0293     time_domain.y = max_time_domain ;
0294 
0295     wavelength_domain.x = w_center ;  // TODO: try to make this constexpr
0296     wavelength_domain.y = w_extent ;
0297 }
0298 
0299 SEVENT_METHOD std::string sevent::descMax() const
0300 {
0301     int w = 5 ;
0302     std::stringstream ss ;
0303     ss
0304         << "sevent::descMax    " << std::endl
0305         << " evt.max_curand    " << std::setw(w) << max_curand   << std::endl
0306         << " evt.max_slot      " << std::setw(w) << max_slot     << std::endl
0307         << " evt.max_genstep   " << std::setw(w) << max_genstep  << std::endl
0308         << " evt.max_photon    " << std::setw(w) << max_photon   << std::endl
0309         << " evt.max_simtrace  " << std::setw(w) << max_simtrace << std::endl
0310         << " evt.max_bounce    " << std::setw(w) << max_bounce   << std::endl
0311         << " evt.max_record    " << std::setw(w) << max_record   << std::endl
0312         << " evt.max_rec       " << std::setw(w) << max_rec      << std::endl
0313         << " evt.max_aux       " << std::setw(w) << max_aux      << std::endl
0314         << " evt.max_sup       " << std::setw(w) << max_sup      << std::endl
0315         << " evt.max_seq       " << std::setw(w) << max_seq      << std::endl
0316         << " evt.max_prd       " << std::setw(w) << max_prd      << std::endl
0317         << " evt.max_tag       " << std::setw(w) << max_tag      << std::endl
0318         << " evt.max_flat      " << std::setw(w) << max_flat     << std::endl
0319         << " evt.mode_lite     " << std::setw(w) << mode_lite    << std::endl
0320         ;
0321 
0322     std::string s = ss.str();
0323     return s ;
0324 }
0325 
0326 SEVENT_METHOD std::string sevent::descNum() const
0327 {
0328     int w = 5 ;
0329     std::stringstream ss ;
0330     ss
0331         << " sevent::descNum          "  << std::endl
0332         << " evt.num_curand           "  << std::setw(w) << num_curand         << std::endl
0333         << " evt.num_genstep          "  << std::setw(w) << num_genstep        << std::endl
0334         << " evt.num_seed             "  << std::setw(w) << num_seed           << std::endl
0335         << " evt.num_photon           "  << std::setw(w) << num_photon         << std::endl
0336         << " evt.num_photonlite       "  << std::setw(w) << num_photonlite     << std::endl
0337         << " evt.num_record           "  << std::setw(w) << num_record         << std::endl
0338         << " evt.num_rec              "  << std::setw(w) << num_rec            << std::endl
0339         << " evt.num_aux              "  << std::setw(w) << num_aux            << std::endl
0340         << " evt.num_sup              "  << std::setw(w) << num_sup            << std::endl
0341         << " evt.num_seq              "  << std::setw(w) << num_seq            << std::endl
0342         << " evt.num_hit              "  << std::setw(w) << num_hit            << std::endl
0343         << " evt.num_hitlite          "  << std::setw(w) << num_hitlite        << std::endl
0344         << " evt.num_hitlitemerged    "  << std::setw(w) << num_hitlitemerged  << std::endl
0345         << " evt.num_simtrace         "  << std::setw(w) << num_simtrace       << std::endl
0346         << " evt.num_prd              "  << std::setw(w) << num_prd            << std::endl
0347         << " evt.num_tag              "  << std::setw(w) << num_tag            << std::endl
0348         << " evt.num_flat             "  << std::setw(w) << num_flat           << std::endl
0349         ;
0350     std::string s = ss.str();
0351     return s ;
0352 }
0353 
0354 SEVENT_METHOD std::string sevent::descBuf() const
0355 {
0356     int w = 5 ;
0357     std::stringstream ss ;
0358     ss
0359         << std::setw(20) << " sevent::descBuf "
0360         << std::endl
0361         << std::setw(20) << " evt.genstep     " << std::setw(w) << ( genstep ? "Y" : "N" ) << " " << std::setw(20) << genstep
0362         << std::setw(20) << " num_genstep "     << std::setw(7) << num_genstep
0363         << std::setw(20) << " max_genstep "     << std::setw(7) << max_genstep
0364         << std::endl
0365         << std::setw(20) << " evt.seed   "      << std::setw(w) << ( seed    ? "Y" : "N" ) << " " << std::setw(20) << seed
0366         << std::setw(20) << " num_seed "        << std::setw(7) << num_seed
0367         << std::setw(20) << " max_photon "      << std::setw(7) << max_photon
0368         << std::endl
0369         << std::setw(20) << " evt.photon "      << std::setw(w) << ( photon  ? "Y" : "N" ) << " " << std::setw(20) << photon
0370         << std::setw(20) << " num_photon "      << std::setw(7) << num_photon
0371         << std::setw(20) << " max_photon "      << std::setw(7) << max_photon
0372         << std::endl
0373         << std::setw(20) << " evt.photonlite "  << std::setw(w) << ( photonlite  ? "Y" : "N" ) << " " << std::setw(20) << photonlite
0374         << std::setw(20) << " num_photonlite "  << std::setw(7) << num_photonlite
0375         << std::setw(20) << " max_photon "      << std::setw(7) << max_photon
0376         << std::endl
0377         << std::setw(20) << " evt.record "      << std::setw(w) << ( record  ? "Y" : "N" ) << " " << std::setw(20) << record
0378         << std::setw(20) << " num_record "      << std::setw(7) << num_record
0379         << std::setw(20) << " max_record "      << std::setw(7) << max_record
0380         << std::endl
0381         << std::setw(20) << " evt.rec "         << std::setw(w) << ( rec  ? "Y" : "N" ) << " " << std::setw(20) << rec
0382         << std::setw(20) << " num_rec "         << std::setw(7) << num_rec
0383         << std::setw(20) << " max_rec "         << std::setw(7) << max_rec
0384         << std::endl
0385         << std::setw(20) << " evt.aux "         << std::setw(w) << ( aux  ? "Y" : "N" ) << " " << std::setw(20) << aux
0386         << std::setw(20) << " num_aux "         << std::setw(7) << num_aux
0387         << std::setw(20) << " max_aux "         << std::setw(7) << max_aux
0388         << std::endl
0389         << std::setw(20) << " evt.sup "         << std::setw(w) << ( sup  ? "Y" : "N" ) << " " << std::setw(20) << sup
0390         << std::setw(20) << " num_sup "         << std::setw(7) << num_sup
0391         << std::setw(20) << " max_sup "         << std::setw(7) << max_sup
0392         << std::endl
0393         << std::setw(20) << " evt.seq "         << std::setw(w) << ( seq     ? "Y" : "N" ) << " " << std::setw(20) << seq
0394         << std::setw(20) << " num_seq "         << std::setw(7) << num_seq
0395         << std::setw(20) << " max_seq "         << std::setw(7) << max_seq
0396         << std::endl
0397         << std::setw(20) << " evt.hit "         << std::setw(w) << ( hit     ? "Y" : "N" ) << " " << std::setw(20) << hit
0398         << std::setw(20) << " num_hit "         << std::setw(7) << num_hit
0399         << std::setw(20) << " max_photon "      << std::setw(7) << max_photon
0400         << std::endl
0401         << std::setw(20) << " evt.simtrace "    << std::setw(w) << ( simtrace  ? "Y" : "N" ) << " " << std::setw(20) << simtrace
0402         << std::setw(20) << " num_simtrace "    << std::setw(7) << num_simtrace
0403         << std::setw(20) << " max_simtrace "    << std::setw(7) << max_simtrace
0404         << std::endl
0405         << std::setw(20) << " evt.prd "         << std::setw(w) << ( prd  ? "Y" : "N" ) << " " << std::setw(20) << prd
0406         << std::setw(20) << " num_prd "         << std::setw(7) << num_prd
0407         << std::setw(20) << " max_prd "         << std::setw(7) << max_prd
0408         << std::endl
0409         << std::setw(20) << " evt.tag "         << std::setw(w) << ( tag  ? "Y" : "N" ) << " " << std::setw(20) << tag
0410         << std::setw(20) << " num_tag "         << std::setw(7) << num_tag
0411         << std::setw(20) << " max_tag "         << std::setw(7) << max_tag
0412         << std::endl
0413         << std::setw(20) << " evt.flat "         << std::setw(w) << ( flat  ? "Y" : "N" ) << " " << std::setw(20) << flat
0414         << std::setw(20) << " num_flat "         << std::setw(7) << num_flat
0415         << std::setw(20) << " max_flat "         << std::setw(7) << max_flat
0416         << std::endl
0417         ;
0418     std::string s = ss.str();
0419     return s ;
0420 }
0421 
0422 
0423 SEVENT_METHOD std::string sevent::desc() const
0424 {
0425     std::stringstream ss ;
0426     ss << descMax() << std::endl ;
0427     ss << descBuf() << std::endl ;
0428     ss << descNum() << std::endl ;
0429     std::string s = ss.str();
0430     return s ;
0431 }
0432 
0433 /**
0434 sevent::get_domain
0435 -------------------
0436 
0437 HMM: could also use metadata (key, value) pairs on the domain NP array
0438 
0439 ana/evt.py::
0440 
0441     2052         post_center, post_extent = self.post_center_extent()  # center and extent are quads, created by combining position and time domain ce
0442     2053         p = self.rx[:,recs,0].astype(np.float32)*post_extent/32767.0 + post_center
0443 
0444 **/
0445 
0446 
0447 SEVENT_METHOD void sevent::get_domain( quad4& dom ) const
0448 {
0449    float4 post_center = make_float4( center_extent.x, center_extent.y, center_extent.z, time_domain.x );
0450    float4 post_extent = make_float4( center_extent.w, center_extent.w, center_extent.w, time_domain.y );
0451    float4 polw_center = make_float4( 0.f, 0.f, 0.f, wavelength_domain.x  );
0452    float4 polw_extent = make_float4( 1.f, 1.f, 1.f, wavelength_domain.y );
0453 
0454    dom.q0.f = post_center ;
0455    dom.q1.f = post_extent ;
0456    dom.q2.f = polw_center ;
0457    dom.q3.f = polw_extent ;
0458 
0459    // xyz duplication allows position and time to be decompressed together
0460    // also polarization and wavelength can be decompressed together using same trick
0461 }
0462 
0463 
0464 SEVENT_METHOD void sevent::get_config( quad4& cfg ) const
0465 {
0466    cfg.q0.u.x = max_genstep ;
0467    cfg.q0.u.y = max_photon ;
0468    cfg.q0.u.z = max_bounce ;
0469    cfg.q0.u.w = 0 ;
0470 
0471    cfg.q1.u.x = max_record ;
0472    cfg.q1.u.y = max_rec ;
0473    cfg.q1.u.z = max_seq ;
0474    cfg.q1.u.w = 0 ;
0475 
0476    cfg.q2.u.x = num_genstep ;
0477    cfg.q2.u.y = num_seed ;
0478    cfg.q2.u.z = num_photon ;
0479    cfg.q2.u.w = num_hit ;
0480 
0481    cfg.q3.u.x = num_record ;
0482    cfg.q3.u.y = num_rec ;
0483    cfg.q3.u.z = num_seq ;
0484    cfg.q3.u.w = 0 ;
0485 }
0486 
0487 
0488 SEVENT_METHOD void sevent::get_meta(std::string& meta) const
0489 {
0490     NP::SetMeta<int64_t>(meta,"evt.max_curand", max_curand);
0491     NP::SetMeta<int64_t>(meta,"evt.max_slot", max_slot);
0492     NP::SetMeta<int64_t>(meta,"evt.max_photon", max_photon);
0493     NP::SetMeta<int64_t>(meta,"evt.num_photon", num_photon);
0494     NP::SetMeta<int64_t>(meta,"evt.num_photonlite", num_photonlite);
0495 
0496     NP::SetMeta<int64_t>(meta,"evt.max_curand/M", max_curand/M);
0497     NP::SetMeta<int64_t>(meta,"evt.max_slot/M", max_slot/M);
0498     NP::SetMeta<int64_t>(meta,"evt.max_photon/M", max_photon/M);
0499     NP::SetMeta<int64_t>(meta,"evt.num_photon/M", num_photon/M);
0500     NP::SetMeta<int64_t>(meta,"evt.num_photonlite/M", num_photonlite/M);
0501 
0502     NP::SetMeta<int64_t>(meta,"evt.max_slot*evt.max_record", max_slot*max_record);
0503     NP::SetMeta<int64_t>(meta,"evt.max_slot*evt.max_record/M", max_slot*max_record/M);
0504 
0505     NP::SetMeta<int64_t>(meta,"evt.max_record", max_record);
0506     NP::SetMeta<int64_t>(meta,"evt.max_record/M", max_record/M);
0507     NP::SetMeta<int64_t>(meta,"evt.max_rec", max_rec);
0508     NP::SetMeta<int64_t>(meta,"evt.max_seq", max_seq);
0509     NP::SetMeta<int64_t>(meta,"evt.max_prd", max_prd);
0510     NP::SetMeta<int64_t>(meta,"evt.max_tag", max_tag);
0511     NP::SetMeta<int64_t>(meta,"evt.max_flat", max_flat);
0512 
0513     NP::SetMeta<int64_t>(meta,"evt.num_record", num_record);
0514     NP::SetMeta<int64_t>(meta,"evt.num_rec", num_rec);
0515     NP::SetMeta<int64_t>(meta,"evt.num_seq", num_seq);
0516     NP::SetMeta<int64_t>(meta,"evt.num_prd", num_prd);
0517     NP::SetMeta<int64_t>(meta,"evt.num_tag", num_tag);
0518     NP::SetMeta<int64_t>(meta,"evt.num_flat", num_flat);
0519 }
0520 
0521 
0522 
0523 /**
0524 sevent::zero
0525 --------------
0526 
0527 CAUTION with QEvt device side running the pointers
0528 are to device buffers that are allocated once based on
0529 configured maxima and then reused for each launch
0530 
0531 **/
0532 SEVENT_METHOD void sevent::zero()
0533 {
0534     index = 0 ;
0535 
0536     //num_slot = 0 ;
0537     num_curand = 0 ;
0538     num_genstep = 0 ;
0539     num_seed  = 0 ;
0540     num_hit = 0 ;
0541     num_photon = 0 ;
0542 
0543     num_hitlite = 0 ;
0544     num_hitlitemerged = 0 ;
0545     num_photonlite = 0 ;
0546 
0547     num_record = 0 ;
0548     num_rec = 0 ;
0549     num_aux = 0 ;
0550     num_sup = 0 ;
0551     num_seq = 0 ;
0552     num_prd = 0 ;
0553     num_tag = 0 ;
0554     num_flat = 0 ;
0555     num_simtrace = 0 ;
0556 
0557 
0558     genstep = nullptr ;
0559     seed = nullptr ;
0560     hit = nullptr ;
0561     hitlite = nullptr ;
0562     hitlitemerged = nullptr ;
0563     photon = nullptr ;
0564     photonlite = nullptr ;
0565 
0566     record = nullptr ;
0567     rec = nullptr ;
0568     aux = nullptr ;
0569     sup = nullptr ;
0570     seq = nullptr ;
0571     prd = nullptr ;
0572     tag = nullptr ;
0573     flat = nullptr ;
0574     simtrace = nullptr ;
0575 }
0576 
0577 
0578 
0579 
0580 
0581 template<typename T> SEVENT_METHOD T*           sevent::get_photon_ptr() const { return nullptr; }
0582 template<>           SEVENT_METHOD sphoton*     sevent::get_photon_ptr<sphoton>() const { return photon; }
0583 template<>           SEVENT_METHOD sphotonlite* sevent::get_photon_ptr<sphotonlite>() const { return photonlite; }
0584 
0585 template<typename T> SEVENT_METHOD size_t sevent::get_photon_num() const { return 0; }
0586 template<>           SEVENT_METHOD size_t sevent::get_photon_num<sphoton>() const { return num_photon; }
0587 template<>           SEVENT_METHOD size_t sevent::get_photon_num<sphotonlite>() const { return num_photonlite; }
0588 
0589 template<typename T> SEVENT_METHOD T**           sevent::get_hitmerged_ptr_ref() {              return nullptr; }
0590 template<>           SEVENT_METHOD sphoton**     sevent::get_hitmerged_ptr_ref<sphoton>() {     return &hitmerged; }
0591 template<>           SEVENT_METHOD sphotonlite** sevent::get_hitmerged_ptr_ref<sphotonlite>() { return &hitlitemerged; }
0592 
0593 template<typename T> SEVENT_METHOD size_t* sevent::get_hitmerged_num_ref() {              return nullptr ; }
0594 template<>           SEVENT_METHOD size_t* sevent::get_hitmerged_num_ref<sphoton>() {     return &num_hitmerged; }
0595 template<>           SEVENT_METHOD size_t* sevent::get_hitmerged_num_ref<sphotonlite>() { return &num_hitlitemerged; }
0596 
0597 
0598 
0599 #endif     // ends host only block
0600 
0601 
0602 #ifndef PRODUCTION
0603 /**
0604 sevent::add_rec
0605 ----------------
0606 
0607 Populates compressed "srec& r" from "const sphoton& p" using the domains
0608 and copies into evt->rec array using the (idx,bounce) slot.
0609 
0610 Note no flags. Have moved to using the full uncompressed "record"
0611 for simplicity instead of this domain compressed one.
0612 
0613 **/
0614 SEVENT_METHOD void  sevent::add_rec( srec& r, unsigned idx, unsigned bounce, const sphoton& p )
0615 {
0616     r.set_position(     p.pos,  center_extent );
0617     r.set_time(         p.time, time_domain );
0618     r.set_polarization( p.pol );
0619     r.set_wavelength(   p.wavelength, wavelength_domain );
0620 
0621     rec[max_rec*idx+bounce] = r ;
0622 }
0623 #endif
0624 
0625 /**
0626 sevent::add_simtrace : fills simtrace[idx] with values from p, prd
0627 --------------------------------------------------------------------
0628 
0629 NB simtrace "photon" *a* is very different from real sphoton.
0630 
0631 a.q0.f
0632     prd.q0.f normal, distance, aka "isect"
0633 a.q1
0634     intersect position from pos+t*dir, 0.
0635 a.q2
0636     initial pos, prd.globalPrimIdx_boundary (formerly boundary, and before that tmin)
0637 a.q3
0638     initial dir, prd.identity
0639 
0640 
0641 +----+----------------+----------------+----------------+-------------------------+--------------------------+
0642 | q  |      x         |      y         |     z          |      w                  |  notes                   |
0643 +====+================+================+================+=========================+==========================+
0644 |    |  nrm.x         |  nrm.y         |  nrm.z         |  t/distance             | surface normal at isect, |
0645 | q0 |                |                |                |                         | ray trace distance       |
0646 |    |                |                |                |                         |                          |
0647 +----+----------------+----------------+----------------+-------------------------+--------------------------+
0648 |    |  pos.x         |  pos.y         | pos.z          |  tmin                   | intersect_position, tmin |
0649 | q1 |                |                |                |                         |                          |
0650 |    |                |                |                |                         |                          |
0651 +----+----------------+----------------+----------------+-------------------------+--------------------------+
0652 |    |  pos0.x        |  pos0.y        |  pos0.z        | (unsigned)              | origin_position,boundary |
0653 | q2 |                |                |                | globalPrimIdx_boundary  |                          |
0654 |    |                |                |                | (2,3)                   |                          |
0655 +----+----------------+----------------+----------------+-------------------------+--------------------------+
0656 |    | dir0.x         |  dir0.y        |  dir0.z        | (unsigned)              | initial direction, id    |
0657 | q3 | (3,0)          |                |                | identity                |                          |
0658 |    |                |                |                | (3,3)                   |                          |
0659 +----+----------------+----------------+----------------+-------------------------+--------------------------+
0660 
0661 NB : simtrace has almost totally different content and layout compared with photon and record
0662 
0663 SEE qsim::generate_photon_simtrace
0664 
0665 HUH: INPUT SIMTRACE p HAS DIFFERENT LAYOUT FROM THE OUTPUT
0666 
0667 
0668 **/
0669 
0670 SEVENT_METHOD void sevent::add_simtrace( unsigned idx, const quad4& p, const quad2* prd, float tmin )
0671 {
0672     float t = prd->distance() ;
0673     quad4 a ;
0674 
0675     a.q0.f  = prd->q0.f ;                // float4 : intersect normal and distance
0676 
0677     a.q1.f.x = p.q0.f.x + t*p.q1.f.x ;   // intersect position
0678     a.q1.f.y = p.q0.f.y + t*p.q1.f.y ;
0679     a.q1.f.z = p.q0.f.z + t*p.q1.f.z ;
0680     a.q1.f.w = tmin ;   // was 0.f
0681 
0682     a.q2.f.x = p.q0.f.x ;                //  initial pos
0683     a.q2.f.y = p.q0.f.y ;
0684     a.q2.f.z = p.q0.f.z ;
0685     a.q2.u.w = prd->globalPrimIdx_boundary() ;
0686 
0687     a.q3.f.x = p.q1.f.x ;                // initial dir
0688     a.q3.f.y = p.q1.f.y ;
0689     a.q3.f.z = p.q1.f.z ;
0690     //a.q3.u.w = prd->identity() ;
0691     a.q3.u.w = prd->iindex_identity() ;
0692 
0693     simtrace[idx] = a ;
0694     // NB this is the from zero slot idx not the offset "photon_idx"
0695     // for multi-launch simtrace with separate buffers for each launch
0696 }
0697