ocl_kernels_ell.cl 5.34 KB
Newer Older
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
// **************************************************************************
//
//    PARALUTION   www.paralution.com
//
//    Copyright (C) 2015  PARALUTION Labs UG (haftungsbeschränkt) & Co. KG
//                        Am Hasensprung 6, 76571 Gaggenau
//                        Handelsregister: Amtsgericht Mannheim, HRA 706051
//                        Vertreten durch:
//                        PARALUTION Labs Verwaltungs UG (haftungsbeschränkt)
//                        Am Hasensprung 6, 76571 Gaggenau
//                        Handelsregister: Amtsgericht Mannheim, HRB 721277
//                        Geschäftsführer: Dimitar Lukarski, Nico Trost
//
//    This program is free software: you can redistribute it and/or modify
//    it under the terms of the GNU General Public License as published by
//    the Free Software Foundation, either version 3 of the License, or
//    (at your option) any later version.
//
//    This program is distributed in the hope that it will be useful,
//    but WITHOUT ANY WARRANTY; without even the implied warranty of
//    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
//    GNU General Public License for more details.
//
//    You should have received a copy of the GNU General Public License
//    along with this program.  If not, see <http://www.gnu.org/licenses/>.
//
// **************************************************************************



// PARALUTION version 1.1.0 


// Nathan Bell and Michael Garland
// Efficient Sparse Matrix-Vector Multiplication on {CUDA}
// NVR-2008-004 / NVIDIA Technical Report
__kernel void kernel_ell_spmv(         const IndexType num_rows,
                                       const IndexType num_cols,
                                       const IndexType num_cols_per_row,
                              __global const IndexType *Acol,
                              __global const ValueType *Aval,
                              __global const ValueType *x,
                              __global       ValueType *y) {

  IndexType row = get_global_id(0);

  if (row < num_rows) {

    ValueType sum = (ValueType) 0;

    for (IndexType n=0; n<num_cols_per_row; ++n) {

      const IndexType ind = n * num_rows + row;
      const IndexType col = Acol[ind];

      if ((col >= 0) && (col < num_cols))
        sum += Aval[ind] * x[col];

    }

    y[row] = sum;

  }

}

// Nathan Bell and Michael Garland
// Efficient Sparse Matrix-Vector Multiplication on {CUDA}
// NVR-2008-004 / NVIDIA Technical Report
__kernel void kernel_ell_add_spmv(         const IndexType num_rows,
                                           const IndexType num_cols,
                                           const IndexType num_cols_per_row,
                                  __global const IndexType *Acol,
                                  __global const ValueType *Aval,
                                           const ValueType scalar,
                                  __global const ValueType *x,
                                  __global       ValueType *y) {

  IndexType row = get_global_id(0);

  if (row < num_rows) {

    ValueType sum = (ValueType) 0;

    for (IndexType n=0; n<num_cols_per_row; ++n) {

      const IndexType ind = n * num_rows + row;
      const IndexType col = Acol[ind];
      
      if ((col >= 0) && (col < num_cols))
        sum += Aval[ind] * x[col];

    }
        
    y[row] += scalar * sum;

  }

}

__kernel void kernel_ell_max_row(         const IndexType nrow,
                                 __global const IndexType *data,
                                 __global       IndexType *out,
                                          const IndexType  GROUP_SIZE,
                                          const IndexType  LOCAL_SIZE) {

    IndexType tid = get_local_id(0);

    __local IndexType sdata[BLOCK_SIZE];

    sdata[tid] = 0;

    IndexType max;

    IndexType gid = GROUP_SIZE * get_group_id(0) + tid;

    for (IndexType i = 0; i < LOCAL_SIZE; ++i, gid += BLOCK_SIZE) {

      if (gid < nrow) {
        max = data[gid+1] - data[gid];
        if (max > sdata[tid])
          sdata[tid] = max;
      }

    }

    barrier(CLK_LOCAL_MEM_FENCE);

    for (IndexType i = BLOCK_SIZE/2; i > 0; i /= 2) {

      if (tid < i)
        if (sdata[tid+i] > sdata[tid]) sdata[tid] = sdata[tid+i];

      barrier(CLK_LOCAL_MEM_FENCE);

    }

    if (tid == 0)
      out[get_group_id(0)] = sdata[tid];

}

__kernel void kernel_ell_csr_to_ell(         const IndexType nrow,
                                             const IndexType max_row,
                                    __global const IndexType *src_row_offset,
                                    __global const IndexType *src_col,
                                    __global const ValueType *src_val,
                                    __global       IndexType *ell_col,
                                    __global       ValueType *ell_val) {

  IndexType ai = get_global_id(0);
  IndexType aj;
  IndexType n = 0;
  IndexType ell_ind;

  if (ai < nrow) {

    for (aj=src_row_offset[ai]; aj<src_row_offset[ai+1]; ++aj) {

      ell_ind = n * nrow + ai;

      ell_col[ell_ind] = src_col[aj];
      ell_val[ell_ind] = src_val[aj];

      ++n;

    }

    for (aj=src_row_offset[ai+1]-src_row_offset[ai]; aj<max_row; ++aj) {

      ell_ind = n * nrow + ai;

      ell_col[ell_ind] = (int) -1;
      ell_val[ell_ind] = (ValueType) 0;

      ++n;

    }

  }

}