Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-05-15 07:41:58

0001 #include <cassert>
0002 #include <csignal>
0003 #include <sstream>
0004 
0005 #include "scuda.h"
0006 #include "squad.h"
0007 
0008 #include "NP.hh"
0009 #include "SLOG.hh"
0010 #include "ssys.h"
0011 
0012 #include "QTex.hh"
0013 #include "QU.hh"
0014 #include "QUDA_CHECK.h"
0015 #include "QWls.hh"
0016 
0017 #include "qwls.h"
0018 
0019 const plog::Severity QWls::LEVEL = SLOG::EnvLevel("QWls", "DEBUG");
0020 
0021 const QWls *QWls::INSTANCE = nullptr;
0022 const QWls *QWls::Get()
0023 {
0024     return INSTANCE;
0025 }
0026 
0027 /**
0028 QWls::QWls
0029 ------------
0030 
0031 1. Narrows ICDF from double to float if needed
0032 2. Uploads ICDF into GPU texture
0033 3. Creates qwls instance with device pointers and uploads it
0034 
0035 **/
0036 
0037 QWls::QWls(const NP *wls_icdf, const NP *mat_map, const NP *time_constants, unsigned hd_factor) :
0038     dsrc(wls_icdf->ebyte == 8 ? wls_icdf : nullptr),
0039     src(wls_icdf->ebyte == 4 ? wls_icdf : NP::MakeNarrow(dsrc)),
0040     tex(MakeWlsTex(src, hd_factor)),
0041     wls(MakeInstance(tex, mat_map, time_constants, hd_factor, time_constants->shape[0])),
0042     d_wls(QU::UploadArray<qwls>(wls, 1, "QWls::QWls/d_wls"))
0043 {
0044     INSTANCE = this;
0045 }
0046 
0047 /**
0048 QWls::MakeWlsTex
0049 -------------------
0050 
0051 Creates a 2D CUDA texture from the ICDF array.
0052 Shape: (num_wls*3, 4096, 1) where 3 = HD layers per material.
0053 
0054 **/
0055 
0056 QTex<float> *QWls::MakeWlsTex(const NP *src, unsigned hd_factor)
0057 {
0058     assert(src);
0059     assert(src->shape.size() == 3);
0060 
0061     unsigned ni = src->shape[0]; // height: num_wls * 3
0062     unsigned nj = src->shape[1]; // width: 4096
0063     unsigned nk = src->shape[2]; // 1
0064 
0065     assert(nk == 1);
0066     assert(nj == 4096);
0067     assert(ni % 3 == 0); // must be multiple of 3 (3 HD layers per material)
0068     assert(src->uifc == 'f' && src->ebyte == 4);
0069 
0070     unsigned ny = ni; // height
0071     unsigned nx = nj; // width
0072 
0073     bool normalizedCoords = true;
0074     QTex<float> *tx = new QTex<float>(nx, ny, src->cvalues<float>(), 'L', normalizedCoords, src);
0075 
0076     tx->setHDFactor(hd_factor);
0077     tx->uploadMeta();
0078 
0079     LOG(LEVEL) << " src " << src->desc() << " nx (width) " << nx << " ny (height) " << ny << " tx.HDFactor "
0080                << tx->getHDFactor();
0081 
0082     return tx;
0083 }
0084 
0085 /**
0086 QWls::MakeInstance
0087 ---------------------
0088 
0089 Creates the host-side qwls struct populated with device pointers.
0090 Uploads material_map and time_constants to device memory.
0091 
0092 **/
0093 
0094 qwls *QWls::MakeInstance(const QTex<float> *tex, const NP *mat_map, const NP *time_constants, unsigned hd_factor,
0095                          unsigned num_wls)
0096 {
0097     assert(mat_map);
0098     assert(time_constants);
0099     assert(mat_map->uifc == 'i' && mat_map->ebyte == 4);
0100     assert(time_constants->uifc == 'f' && time_constants->ebyte == 4);
0101 
0102     qwls *w = new qwls;
0103     w->wls_tex = tex->texObj;
0104     w->hd_factor = hd_factor;
0105     w->num_wls = num_wls;
0106     w->tex_height = tex->height;
0107 
0108     // Upload material_map to device
0109     unsigned num_mat = mat_map->shape[0];
0110     int *d_mat_map = nullptr;
0111     size_t mat_map_size = num_mat * sizeof(int);
0112     QUDA_CHECK(cudaMalloc(reinterpret_cast<void **>(&d_mat_map), mat_map_size));
0113     QUDA_CHECK(cudaMemcpy(d_mat_map, mat_map->cvalues<int>(), mat_map_size, cudaMemcpyHostToDevice));
0114     w->material_map = d_mat_map;
0115 
0116     // Upload time_constants to device
0117     float *d_tc = nullptr;
0118     size_t tc_size = num_wls * sizeof(float);
0119     QUDA_CHECK(cudaMalloc(reinterpret_cast<void **>(&d_tc), tc_size));
0120     QUDA_CHECK(cudaMemcpy(d_tc, time_constants->cvalues<float>(), tc_size, cudaMemcpyHostToDevice));
0121     w->time_constants = d_tc;
0122 
0123     return w;
0124 }
0125 
0126 std::string QWls::desc() const
0127 {
0128     std::stringstream ss;
0129     ss << "QWls" << " dsrc " << (dsrc ? dsrc->desc() : "-") << " src " << (src ? src->desc() : "-") << " tex "
0130        << (tex ? tex->desc() : "-");
0131     return ss.str();
0132 }