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
|
//------------------------------------------------------------------------------
// GB_extract_vector_list: extract vector indices for all entries in a matrix
//------------------------------------------------------------------------------
// SuiteSparse:GraphBLAS, Timothy A. Davis, (c) 2017-2022, All Rights Reserved.
// SPDX-License-Identifier: Apache-2.0
//------------------------------------------------------------------------------
// Constructs a list of vector indices for each entry in a matrix. Creates
// the output J for GB_extractTuples, and I for GB_transpose when the qsort
// method is used.
// TODO: use #include "GB_positional_op_ijp.c" here
#include "GB_ek_slice.h"
#define GB_FREE_ALL \
{ \
GB_WERK_POP (A_ek_slicing, int64_t) ; \
}
GrB_Info GB_extract_vector_list // extract vector list from a matrix
(
// output:
int64_t *restrict J, // size nnz(A) or more
// input:
const GrB_Matrix A,
GB_Context Context
)
{
//--------------------------------------------------------------------------
// check inputs
//--------------------------------------------------------------------------
ASSERT (J != NULL) ;
ASSERT (A != NULL) ;
ASSERT (GB_JUMBLED_OK (A)) ; // pattern not accessed
ASSERT (GB_ZOMBIES_OK (A)) ; // pattern not accessed
ASSERT (!GB_IS_BITMAP (A)) ;
//--------------------------------------------------------------------------
// get A
//--------------------------------------------------------------------------
const int64_t *restrict Ap = A->p ;
const int64_t *restrict Ah = A->h ;
const int64_t avlen = A->vlen ;
//--------------------------------------------------------------------------
// determine the max number of threads to use
//--------------------------------------------------------------------------
GB_GET_NTHREADS_MAX (nthreads_max, chunk, Context) ;
//--------------------------------------------------------------------------
// slice the entries for each task
//--------------------------------------------------------------------------
GB_WERK_DECLARE (A_ek_slicing, int64_t) ;
int A_ntasks, A_nthreads ;
GB_SLICE_MATRIX (A, 2, chunk) ;
//--------------------------------------------------------------------------
// extract the vector index for each entry
//--------------------------------------------------------------------------
int tid ;
#pragma omp parallel for num_threads(A_nthreads) schedule(dynamic,1)
for (tid = 0 ; tid < A_ntasks ; tid++)
{
// if kfirst > klast then task tid does no work at all
int64_t kfirst = kfirst_Aslice [tid] ;
int64_t klast = klast_Aslice [tid] ;
for (int64_t k = kfirst ; k <= klast ; k++)
{
//------------------------------------------------------------------
// find the part of A(:,k) to be operated on by this task
//------------------------------------------------------------------
int64_t j = GBH (Ah, k) ;
int64_t pA_start, pA_end ;
GB_get_pA (&pA_start, &pA_end, tid, k,
kfirst, klast, pstart_Aslice, Ap, avlen) ;
//------------------------------------------------------------------
// extract vector indices of A(:,j)
//------------------------------------------------------------------
for (int64_t p = pA_start ; p < pA_end ; p++)
{
J [p] = j ;
}
}
}
//--------------------------------------------------------------------------
// free workspace and return result
//--------------------------------------------------------------------------
GB_FREE_ALL ;
return (GrB_SUCCESS) ;
}
|