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
|
//------------------------------------------------------------------------------
// GB_subassign_06n: C(I,J)<M> = A ; no S
//------------------------------------------------------------------------------
// SuiteSparse:GraphBLAS, Timothy A. Davis, (c) 2017-2022, All Rights Reserved.
// SPDX-License-Identifier: Apache-2.0
//------------------------------------------------------------------------------
// Method 06n: C(I,J)<M> = A ; no S
// M: present
// Mask_comp: false
// C_replace: false
// accum: NULL
// A: matrix
// S: none (see also GB_subassign_06s)
// FULL: if A and C are dense, then C remains dense.
// If A is sparse and C dense, C will likely become sparse, except if M(i,j)=0
// wherever A(i,j) is not present. So if M==A is aliased and A is sparse, then
// C remains dense. Need C(I,J)<A,struct>=A kernel. Then in that case, if C
// is dense it remains dense, even if A is sparse. If that change is made,
// this kernel can start with converting C to sparse if A is sparse.
// C is not bitmap: GB_bitmap_assign is used if C is bitmap.
// M and A are not bitmap: 06s is used instead, if M or A are bitmap.
#include "GB_subassign_methods.h"
GrB_Info GB_subassign_06n
(
GrB_Matrix C,
// input:
const GrB_Index *I,
const int64_t nI,
const int Ikind,
const int64_t Icolon [3],
const GrB_Index *J,
const int64_t nJ,
const int Jkind,
const int64_t Jcolon [3],
const GrB_Matrix M,
const bool Mask_struct,
const GrB_Matrix A,
GB_Context Context
)
{
//--------------------------------------------------------------------------
// check inputs
//--------------------------------------------------------------------------
ASSERT (!GB_IS_BITMAP (C)) ; ASSERT (!GB_IS_FULL (C)) ;
ASSERT (!GB_IS_BITMAP (M)) ; // Method 06n is not used for M bitmap
ASSERT (!GB_IS_BITMAP (A)) ; // Method 06n is not used for A bitmap
ASSERT (!GB_aliased (C, M)) ; // NO ALIAS of C==M
ASSERT (!GB_aliased (C, A)) ; // NO ALIAS of C==A
ASSERT_MATRIX_OK (C, "C input for 06n", GB0) ;
ASSERT_MATRIX_OK (M, "M input for 06n", GB0) ;
ASSERT_MATRIX_OK (A, "A input for 06n", GB0) ;
//--------------------------------------------------------------------------
// get inputs
//--------------------------------------------------------------------------
GB_EMPTY_TASKLIST ;
GB_MATRIX_WAIT_IF_JUMBLED (C) ;
GB_MATRIX_WAIT_IF_JUMBLED (M) ;
GB_MATRIX_WAIT_IF_JUMBLED (A) ;
GB_GET_C ; // C must not be bitmap
int64_t zorig = C->nzombies ;
const int64_t Cnvec = C->nvec ;
const int64_t *restrict Ch = C->h ;
const int64_t *restrict Cp = C->p ;
const bool C_is_hyper = (Ch != NULL) ;
GB_GET_C_HYPER_HASH ;
GB_GET_MASK ;
GB_GET_A ;
const int64_t *restrict Ah = A->h ;
const int64_t Anvec = A->nvec ;
const bool A_is_hyper = (Ah != NULL) ;
GrB_BinaryOp accum = NULL ;
GB_OK (GB_hyper_hash_build (A, Context)) ;
const int64_t *restrict A_Yp = (A_is_hyper) ? A->Y->p : NULL ;
const int64_t *restrict A_Yi = (A_is_hyper) ? A->Y->i : NULL ;
const int64_t *restrict A_Yx = (A_is_hyper) ? A->Y->x : NULL ;
const int64_t A_hash_bits = (A_is_hyper) ? (A->Y->vdim - 1) : 0 ;
//--------------------------------------------------------------------------
// Method 06n: C(I,J)<M> = A ; no S
//--------------------------------------------------------------------------
// Time: O(nnz(M)*(log(a)+log(c)), where a and c are the # of entries in a
// vector of A and C, respectively. The entries in the intersection of M
// (where the entries are true) and the matrix addition C(I,J)+A must be
// examined. This method scans M, and searches for entries in A and C(I,J)
// using two binary searches. If M is very dense, this method can be
// slower than Method 06s. This method is selected if nnz (A) >= nnz (M).
// Compare with Methods 05 and 07, which use a similar algorithmic outline
// and parallelization strategy.
//--------------------------------------------------------------------------
// Parallel: slice M into coarse/fine tasks (Method 05, 06n, 07)
//--------------------------------------------------------------------------
GB_SUBASSIGN_ONE_SLICE (M) ; // M cannot be jumbled
//--------------------------------------------------------------------------
// phase 1: create zombies, update entries, and count pending tuples
//--------------------------------------------------------------------------
#pragma omp parallel for num_threads(nthreads) schedule(dynamic,1) \
reduction(+:nzombies)
for (taskid = 0 ; taskid < ntasks ; taskid++)
{
//----------------------------------------------------------------------
// get the task descriptor
//----------------------------------------------------------------------
GB_GET_TASK_DESCRIPTOR_PHASE1 ;
//----------------------------------------------------------------------
// compute all vectors in this task
//----------------------------------------------------------------------
for (int64_t k = kfirst ; k <= klast ; k++)
{
//------------------------------------------------------------------
// get j, the kth vector of M
//------------------------------------------------------------------
int64_t j = GBH (Mh, k) ;
GB_GET_VECTOR (pM, pM_end, pA, pA_end, Mp, k, Mvlen) ;
int64_t mjnz = pM_end - pM ;
if (mjnz == 0) continue ;
//------------------------------------------------------------------
// get A(:,j)
//------------------------------------------------------------------
int64_t pA, pA_end ;
GB_LOOKUP_VECTOR (pA, pA_end, A, j) ;
int64_t ajnz = pA_end - pA ;
bool ajdense = (ajnz == Avlen) ;
int64_t pA_start = pA ;
//------------------------------------------------------------------
// get jC, the corresponding vector of C
//------------------------------------------------------------------
GB_LOOKUP_VECTOR_jC (fine_task, taskid) ;
int64_t cjnz = pC_end - pC_start ;
if (cjnz == 0 && ajnz == 0) continue ;
bool cjdense = (cjnz == Cvlen) ;
//------------------------------------------------------------------
// C(I,jC)<M(:,j)> = A(:,j) ; no S
//------------------------------------------------------------------
if (cjdense && ajdense)
{
//--------------------------------------------------------------
// C(:,jC) and A(:,j) are both dense
//--------------------------------------------------------------
for ( ; pM < pM_end ; pM++)
{
//----------------------------------------------------------
// update C(iC,jC), but only if M(iA,j) allows it
//----------------------------------------------------------
if (GB_mcast (Mx, pM, msize))
{
int64_t iA = GBI (Mi, pM, Mvlen) ;
GB_iC_DENSE_LOOKUP ;
// find iA in A(:,j)
// A(:,j) is dense; no need for binary search
pA = pA_start + iA ;
ASSERT (GBI (Ai, pA, Avlen) == iA) ;
// ----[C A 1] or [X A 1]-----------------------
// [C A 1]: action: ( =A ): copy A to C, no acc
// [X A 1]: action: ( undelete ): zombie lives
GB_noaccum_C_A_1_matrix ;
}
}
}
else if (cjdense)
{
//--------------------------------------------------------------
// C(:,jC) is dense, A(:,j) is sparse
//--------------------------------------------------------------
for ( ; pM < pM_end ; pM++)
{
//----------------------------------------------------------
// update C(iC,jC), but only if M(iA,j) allows it
//----------------------------------------------------------
if (GB_mcast (Mx, pM, msize))
{
int64_t iA = GBI (Mi, pM, Mvlen) ;
GB_iC_DENSE_LOOKUP ;
// find iA in A(:,j)
bool aij_found ;
int64_t apright = pA_end - 1 ;
GB_BINARY_SEARCH (iA, Ai, pA, apright, aij_found) ;
if (!aij_found)
{
// C (iC,jC) is present but A (i,j) is not
// ----[C . 1] or [X . 1]---------------------------
// [C . 1]: action: ( delete ): becomes zombie
// [X . 1]: action: ( X ): still zombie
GB_DELETE_ENTRY ;
}
else
{
// ----[C A 1] or [X A 1]---------------------------
// [C A 1]: action: ( =A ): copy A to C, no accum
// [X A 1]: action: ( undelete ): zombie lives
GB_noaccum_C_A_1_matrix ;
}
}
}
}
else if (ajdense)
{
//--------------------------------------------------------------
// C(:,jC) is sparse, A(:,j) is dense
//--------------------------------------------------------------
for ( ; pM < pM_end ; pM++)
{
//----------------------------------------------------------
// update C(iC,jC), but only if M(iA,j) allows it
//----------------------------------------------------------
if (GB_mcast (Mx, pM, msize))
{
int64_t iA = GBI (Mi, pM, Mvlen) ;
// find C(iC,jC) in C(:,jC)
GB_iC_BINARY_SEARCH ;
// lookup iA in A(:,j)
pA = pA_start + iA ;
ASSERT (GBI (Ai, pA, Avlen) == iA) ;
if (cij_found)
{
// ----[C A 1] or [X A 1]---------------------------
// [C A 1]: action: ( =A ): copy A into C, no accum
// [X A 1]: action: ( undelete ): zombie lives
GB_noaccum_C_A_1_matrix ;
}
else
{
// C (iC,jC) is not present, A (i,j) is present
// ----[. A 1]--------------------------------------
// [. A 1]: action: ( insert )
task_pending++ ;
}
}
}
}
else
{
//--------------------------------------------------------------
// C(:,jC) and A(:,j) are both sparse
//--------------------------------------------------------------
for ( ; pM < pM_end ; pM++)
{
//----------------------------------------------------------
// update C(iC,jC), but only if M(iA,j) allows it
//----------------------------------------------------------
if (GB_mcast (Mx, pM, msize))
{
int64_t iA = GBI (Mi, pM, Mvlen) ;
// find C(iC,jC) in C(:,jC)
GB_iC_BINARY_SEARCH ;
// find iA in A(:,j)
bool aij_found ;
int64_t apright = pA_end - 1 ;
GB_BINARY_SEARCH (iA, Ai, pA, apright, aij_found) ;
if (cij_found && aij_found)
{
// ----[C A 1] or [X A 1]---------------------------
// [C A 1]: action: ( =A ): copy A into C, no accum
// [X A 1]: action: ( undelete ): zombie lives
GB_noaccum_C_A_1_matrix ;
}
else if (!cij_found && aij_found)
{
// C (iC,jC) is not present, A (i,j) is present
// ----[. A 1]--------------------------------------
// [. A 1]: action: ( insert )
task_pending++ ;
}
else if (cij_found && !aij_found)
{
// C (iC,jC) is present but A (i,j) is not
// ----[C . 1] or [X . 1]---------------------------
// [C . 1]: action: ( delete ): becomes zombie
// [X . 1]: action: ( X ): still zombie
GB_DELETE_ENTRY ;
}
}
}
}
}
GB_PHASE1_TASK_WRAPUP ;
}
//--------------------------------------------------------------------------
// phase 2: insert pending tuples
//--------------------------------------------------------------------------
GB_PENDING_CUMSUM ;
zorig = C->nzombies ;
#pragma omp parallel for num_threads(nthreads) schedule(dynamic,1) \
reduction(&&:pending_sorted)
for (taskid = 0 ; taskid < ntasks ; taskid++)
{
//----------------------------------------------------------------------
// get the task descriptor
//----------------------------------------------------------------------
GB_GET_TASK_DESCRIPTOR_PHASE2 ;
//----------------------------------------------------------------------
// compute all vectors in this task
//----------------------------------------------------------------------
for (int64_t k = kfirst ; k <= klast ; k++)
{
//------------------------------------------------------------------
// get j, the kth vector of M
//------------------------------------------------------------------
int64_t j = GBH (Mh, k) ;
GB_GET_VECTOR (pM, pM_end, pA, pA_end, Mp, k, Mvlen) ;
int64_t mjnz = pM_end - pM ;
if (mjnz == 0) continue ;
//------------------------------------------------------------------
// get A(:,j)
//------------------------------------------------------------------
int64_t pA, pA_end ;
GB_LOOKUP_VECTOR (pA, pA_end, A, j) ;
int64_t ajnz = pA_end - pA ;
if (ajnz == 0) continue ;
bool ajdense = (ajnz == Avlen) ;
int64_t pA_start = pA ;
//------------------------------------------------------------------
// get jC, the corresponding vector of C
//------------------------------------------------------------------
GB_LOOKUP_VECTOR_jC (fine_task, taskid) ;
bool cjdense = ((pC_end - pC_start) == Cvlen) ;
//------------------------------------------------------------------
// C(I,jC)<M(:,j)> = A(:,j)
//------------------------------------------------------------------
if (!cjdense)
{
//--------------------------------------------------------------
// C(:,jC) is sparse; use binary search for C
//--------------------------------------------------------------
for ( ; pM < pM_end ; pM++)
{
//----------------------------------------------------------
// update C(iC,jC), but only if M(iA,j) allows it
//----------------------------------------------------------
if (GB_mcast (Mx, pM, msize))
{
int64_t iA = GBI (Mi, pM, Mvlen) ;
// find iA in A(:,j)
if (ajdense)
{
// A(:,j) is dense; no need for binary search
pA = pA_start + iA ;
ASSERT (GBI (Ai, pA, Avlen) == iA) ;
}
else
{
// A(:,j) is sparse; use binary search
int64_t apright = pA_end - 1 ;
bool aij_found ;
GB_BINARY_SEARCH (iA, Ai, pA, apright, aij_found) ;
if (!aij_found) continue ;
}
// find C(iC,jC) in C(:,jC)
GB_iC_BINARY_SEARCH ;
if (!cij_found)
{
// C (iC,jC) is not present, A (i,j) is present
// ----[. A 1]--------------------------------------
// [. A 1]: action: ( insert )
GB_PENDING_INSERT_aij ;
}
}
}
}
}
GB_PHASE2_TASK_WRAPUP ;
}
//--------------------------------------------------------------------------
// finalize the matrix and return result
//--------------------------------------------------------------------------
GB_SUBASSIGN_WRAPUP ;
}
|