scirs2-core 0.4.3

Core utilities and common functionality for SciRS2 (scirs2-core)
Documentation
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
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
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
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
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
//! General matrix-vector multiplication (GEMV) kernels for GPU
//!
//! Implements y = alpha * A * x + beta * y where:
//! - A is an M x N matrix
//! - x is an N-dimensional vector
//! - y is an M-dimensional vector
//! - alpha and beta are scalar values

use std::collections::HashMap;

use crate::gpu::kernels::{
    BaseKernel, DataType, GpuKernel, KernelMetadata, KernelParams, OperationType,
};
use crate::gpu::{GpuBackend, GpuError};

/// General matrix-vector multiplication kernel
pub struct GemvKernel {
    base: BaseKernel,
}

impl Default for GemvKernel {
    fn default() -> Self {
        Self::new()
    }
}

impl GemvKernel {
    /// Create a new GEMV kernel
    pub fn new() -> Self {
        let metadata = KernelMetadata {
            workgroup_size: [256, 1, 1],
            local_memory_usage: 1024, // 1 KB local memory for reduction
            supports_tensor_cores: false,
            operationtype: OperationType::ComputeIntensive,
            backend_metadata: HashMap::new(),
        };

        let cuda_source = r#"
extern "C" __global__ void gemv(
    const float* __restrict__ matrix,  // M x N matrix (row-major)
    const float* __restrict__ vector,  // N-dimensional vector
    float* __restrict__ result,        // M-dimensional result vector
    float alpha,
    float beta,
    int M,  // Number of rows
    int N   // Number of columns
) {
    int row = blockIdx.x * blockDim.x + threadIdx.x;

    if (row < M) {
        float sum = 0.0f;

        // Compute dot product of matrix row with vector
        for (int col = 0; col < N; col++) {
            sum += matrix[row * N + col] * vector[col];
        }

        // Apply alpha and beta coefficients
        result[row] = alpha * sum + beta * result[row];
    }
}

// Optimized version using shared memory for larger matrices
extern "C" __global__ void gemv_shared(
    const float* __restrict__ matrix,
    const float* __restrict__ vector,
    float* __restrict__ result,
    float alpha,
    float beta,
    int M,
    int N
) {
    extern __shared__ float shared_vector[];

    int row = blockIdx.x * blockDim.x + threadIdx.x;
    int tid = threadIdx.x;

    // Load vector into shared memory in chunks
    for (int i = tid; i < N; i += blockDim.x) {
        if (i < N) {
            shared_vector[i] = vector[i];
        }
    }
    __syncthreads();

    if (row < M) {
        float sum = 0.0f;

        // Compute dot product using shared memory vector
        for (int col = 0; col < N; col++) {
            sum += matrix[row * N + col] * shared_vector[col];
        }

        result[row] = alpha * sum + beta * result[row];
    }
}
"#
        .to_string();

        let rocm_source = cuda_source.clone();

        let wgpu_source = r#"
struct Uniforms {
    alpha: f32,
    beta: f32,
    M: u32,  // Number of rows
    N: u32,  // Number of columns
};

@group(0) @binding(0) var<uniform> uniforms: Uniforms;
@group(0) @binding(1) var<storage, read> matrix: array<f32>;   // M x N matrix
@group(0) @binding(2) var<storage, read> vector: array<f32>;   // N-dimensional vector
@group(0) @binding(3) var<storage, write> result: array<f32>;  // M-dimensional result

@compute @workgroup_size(256)
fn gemv(@builtin(global_invocation_id) global_id: vec3<u32>) {
    let row = global_id.x;

    if (row < uniforms.M) {
        var sum = 0.0;

        // Compute dot product of matrix row with vector
        for (var col = 0u; col < uniforms.N; col = col + 1u) {
            let matrix_idx = row * uniforms.N + col;
            sum = sum + matrix[matrix_idx] * vector[col];
        }

        // Apply alpha and beta coefficients
        result[row] = uniforms.alpha * sum + uniforms.beta * result[row];
    }
}
"#
        .to_string();

        let metal_source = r#"
#include <metal_stdlib>
using namespace metal;

kernel void gemv(
    const device float* matrix [[buffer(0)]],    // M x N matrix
    const device float* vector [[buffer(1)]],    // N-dimensional vector
    device float* result [[buffer(2)]],          // M-dimensional result
    constant float& alpha [[buffer(3)]],
    constant float& beta [[buffer(4)]],
    constant uint& M [[buffer(5)]],              // Number of rows
    constant uint& N [[buffer(6)]],              // Number of columns
    uint gid [[thread_position_in_grid]])
{
    if (gid < M) {
        float sum = 0.0f;

        // Compute dot product of matrix row with vector
        for (uint col = 0; col < N; col++) {
            sum += matrix[gid * N + col] * vector[col];
        }

        // Apply alpha and beta coefficients
        result[gid] = alpha * sum + beta * result[gid];
    }
}

// Optimized version using threadgroup memory
kernel void gemv_tiled(
    const device float* matrix [[buffer(0)]],
    const device float* vector [[buffer(1)]],
    device float* result [[buffer(2)]],
    constant float& alpha [[buffer(3)]],
    constant float& beta [[buffer(4)]],
    constant uint& M [[buffer(5)]],
    constant uint& N [[buffer(6)]],
    uint gid [[thread_position_in_grid]],
    uint lid [[thread_position_in_threadgroup]],
    uint blockSize [[threads_per_threadgroup]])
{
    threadgroup float shared_vector[256];  // Shared vector storage

    // Load vector into threadgroup memory
    for (uint i = lid; i < N; i += blockSize) {
        if (i < N) {
            shared_vector[i] = vector[i];
        }
    }
    threadgroup_barrier(mem_flags::mem_threadgroup);

    if (gid < M) {
        float sum = 0.0f;

        // Compute using shared vector
        for (uint col = 0; col < N; col++) {
            sum += matrix[gid * N + col] * shared_vector[col];
        }

        result[gid] = alpha * sum + beta * result[gid];
    }
}
"#
        .to_string();

        let opencl_source = r#"
__kernel void gemv(
    __global const float* matrix,   // M x N matrix
    __global const float* vector,   // N-dimensional vector
    __global float* result,         // M-dimensional result
    const float alpha,
    const float beta,
    const int M,                    // Number of rows
    const int N)                    // Number of columns
{
    int row = get_global_id(0);

    if (row < M) {
        float sum = 0.0f;

        // Compute dot product of matrix row with vector
        for (int col = 0; col < N; col++) {
            sum += matrix[row * N + col] * vector[col];
        }

        // Apply alpha and beta coefficients
        result[row] = alpha * sum + beta * result[row];
    }
}

// Version with local memory optimization
__kernel void gemv_local(
    __global const float* matrix,
    __global const float* vector,
    __global float* result,
    const float alpha,
    const float beta,
    const int M,
    const int N,
    __local float* local_vector)
{
    int row = get_global_id(0);
    int lid = get_local_id(0);
    int local_size = get_local_size(0);

    // Load vector into local memory
    for (int i = lid; i < N; i += local_size) {
        if (i < N) {
            local_vector[i] = vector[i];
        }
    }
    barrier(CLK_LOCAL_MEM_FENCE);

    if (row < M) {
        float sum = 0.0f;

        // Compute using local vector
        for (int col = 0; col < N; col++) {
            sum += matrix[row * N + col] * local_vector[col];
        }

        result[row] = alpha * sum + beta * result[row];
    }
}
"#
        .to_string();

        Self {
            base: BaseKernel::new(
                "gemv",
                &cuda_source,
                &rocm_source,
                &wgpu_source,
                &metal_source,
                &opencl_source,
                metadata,
            ),
        }
    }
}

impl GpuKernel for GemvKernel {
    fn name(&self) -> &str {
        self.base.name()
    }

    fn source_for_backend(&self, backend: GpuBackend) -> Result<String, GpuError> {
        self.base.source_for_backend(backend)
    }

    fn metadata(&self) -> KernelMetadata {
        self.base.metadata()
    }

    fn can_specialize(&self, params: &KernelParams) -> bool {
        matches!(params.datatype, DataType::Float32 | DataType::Float64)
    }

    fn specialize(&self, params: &KernelParams) -> Result<Box<dyn GpuKernel>, GpuError> {
        if !self.can_specialize(params) {
            return Err(GpuError::SpecializationNotSupported);
        }

        // For now, return the same kernel (no type specialization implemented yet)
        Ok(Box::new(Self::new()))
    }
}

/// Batched GEMV kernel for processing multiple matrix-vector multiplications
pub struct BatchGemvKernel {
    base: BaseKernel,
}

impl Default for BatchGemvKernel {
    fn default() -> Self {
        Self::new()
    }
}

impl BatchGemvKernel {
    /// Create a new batched GEMV kernel
    pub fn new() -> Self {
        let metadata = KernelMetadata {
            workgroup_size: [16, 16, 1],
            local_memory_usage: 2048,
            supports_tensor_cores: false,
            operationtype: OperationType::ComputeIntensive,
            backend_metadata: HashMap::new(),
        };

        let cuda_source = r#"
extern "C" __global__ void batch_gemv(
    const float* __restrict__ matrices,  // Batch of M x N matrices
    const float* __restrict__ vectors,   // Batch of N-dimensional vectors
    float* __restrict__ results,         // Batch of M-dimensional results
    float alpha,
    float beta,
    int batch_size,
    int M,  // Number of rows per matrix
    int N   // Number of columns per matrix
) {
    int batch_idx = blockIdx.z;
    int row = blockIdx.x * blockDim.x + threadIdx.x;

    if (batch_idx < batch_size && row < M) {
        // Calculate offsets for this batch
        int matrix_offset = batch_idx * M * N;
        int vector_offset = batch_idx * N;
        int result_offset = batch_idx * M;

        float sum = 0.0f;

        // Compute dot product of matrix row with vector
        for (int col = 0; col < N; col++) {
            sum += matrices[matrix_offset + row * N + col] *
                   vectors[vector_offset + col];
        }

        // Apply alpha and beta coefficients
        results[result_offset + row] = alpha * sum + beta * results[result_offset + row];
    }
}
"#
        .to_string();

        let rocm_source = cuda_source.clone();

        let wgpu_source = r#"
struct Uniforms {
    alpha: f32,
    beta: f32,
    batch_size: u32,
    M: u32,  // Number of rows per matrix
    N: u32,  // Number of columns per matrix
};

@group(0) @binding(0) var<uniform> uniforms: Uniforms;
@group(0) @binding(1) var<storage, read> matrices: array<f32>;  // Batch of matrices
@group(0) @binding(2) var<storage, read> vectors: array<f32>;   // Batch of vectors
@group(0) @binding(3) var<storage, write> results: array<f32>;  // Batch of results

@compute @workgroup_size(16, 16, 1)
fn batch_gemv(@builtin(global_invocation_id) global_id: vec3<u32>) {
    let batch_idx = global_id.z;
    let row = global_id.x;

    if (batch_idx < uniforms.batch_size && row < uniforms.M) {
        // Calculate offsets for this batch
        let matrix_offset = batch_idx * uniforms.M * uniforms.N;
        let vector_offset = batch_idx * uniforms.N;
        let result_offset = batch_idx * uniforms.M;

        var sum = 0.0;

        // Compute dot product
        for (var col = 0u; col < uniforms.N; col = col + 1u) {
            let matrix_idx = matrix_offset + row * uniforms.N + col;
            let vector_idx = vector_offset + col;
            sum = sum + matrices[matrix_idx] * vectors[vector_idx];
        }

        // Apply coefficients
        let result_idx = result_offset + row;
        results[result_idx] = uniforms.alpha * sum + uniforms.beta * results[result_idx];
    }
}
"#
        .to_string();

        let metal_source = r#"
#include <metal_stdlib>
using namespace metal;

kernel void batch_gemv(
    const device float* matrices [[buffer(0)]],   // Batch of matrices
    const device float* vectors [[buffer(1)]],    // Batch of vectors
    device float* results [[buffer(2)]],          // Batch of results
    constant float& alpha [[buffer(3)]],
    constant float& beta [[buffer(4)]],
    constant uint& batch_size [[buffer(5)]],
    constant uint& M [[buffer(6)]],               // Rows per matrix
    constant uint& N [[buffer(7)]],               // Columns per matrix
    uint3 gid [[thread_position_in_grid]])
{
    uint batch_idx = gid.z;
    uint row = gid.x;

    if (batch_idx < batch_size && row < M) {
        // Calculate offsets
        uint matrix_offset = batch_idx * M * N;
        uint vector_offset = batch_idx * N;
        uint result_offset = batch_idx * M;

        float sum = 0.0f;

        // Compute dot product
        for (uint col = 0; col < N; col++) {
            sum += matrices[matrix_offset + row * N + col] *
                   vectors[vector_offset + col];
        }

        // Apply coefficients
        results[result_offset + row] = alpha * sum + beta * results[result_offset + row];
    }
}
"#
        .to_string();

        let opencl_source = r#"
__kernel void batch_gemv(
    __global const float* matrices,
    __global const float* vectors,
    __global float* results,
    const float alpha,
    const float beta,
    const int batch_size,
    const int M,
    const int N)
{
    int batch_idx = get_global_id(2);
    int row = get_global_id(0);

    if (batch_idx < batch_size && row < M) {
        // Calculate offsets
        int matrix_offset = batch_idx * M * N;
        int vector_offset = batch_idx * N;
        int result_offset = batch_idx * M;

        float sum = 0.0f;

        // Compute dot product
        for (int col = 0; col < N; col++) {
            sum += matrices[matrix_offset + row * N + col] *
                   vectors[vector_offset + col];
        }

        // Apply coefficients
        results[result_offset + row] = alpha * sum + beta * results[result_offset + row];
    }
}
"#
        .to_string();

        Self {
            base: BaseKernel::new(
                "batch_gemv",
                &cuda_source,
                &rocm_source,
                &wgpu_source,
                &metal_source,
                &opencl_source,
                metadata,
            ),
        }
    }
}

impl GpuKernel for BatchGemvKernel {
    fn name(&self) -> &str {
        self.base.name()
    }

    fn source_for_backend(&self, backend: GpuBackend) -> Result<String, GpuError> {
        self.base.source_for_backend(backend)
    }

    fn metadata(&self) -> KernelMetadata {
        self.base.metadata()
    }

    fn can_specialize(&self, params: &KernelParams) -> bool {
        matches!(params.datatype, DataType::Float32 | DataType::Float64)
    }

    fn specialize(&self, params: &KernelParams) -> Result<Box<dyn GpuKernel>, GpuError> {
        if !self.can_specialize(params) {
            return Err(GpuError::SpecializationNotSupported);
        }

        Ok(Box::new(Self::new()))
    }
}