Skip to content

Conversation

@GomezGab
Copy link
Collaborator

@GomezGab GomezGab commented May 11, 2025

LAGraph_RichClubCoefficient.c, rcc_demo.c, test_RichClubCoefficient.c, and LG_check_RCC.c

  • Get the rich club coefficient of a graph.
  • Given a Symmetric Graph with no self-edges, LAGraph_RichClubCoefficient will calculate the rich club coefficients of the graph.
  • The values will be output as a sparse GrB_Vector, and the rich club coefficient of k will be found at the closest entry at or above k.
  • check_RCC is the pure C version of the algorithm

LAGraph_FastAssign.c and FastAssign_demo.c

  • A GBv10 algorithm
  • Essentially moves CSC trick found in FastSV7 and in SwapEdges (WIP) to its own function.
  • This method is used for cases where you want to build a vector equivalent to the following for loop:
    for (j = 0 ; j < n ; j++)
      {
          uint64_t i = I_vec [j] ;
          c [i] += X_vec [j] ;
      }
    
  • Where += is defined by the accum and dup operators.
  • See more details in function comments.

LG_CC_FastSV7_FA.c, cc_demo.c, and test_ConnectedComponents.c

  • Demonstrates a use case for FastAssign by using FastAssign in the FastSV7 algorithm for GraphBLASv10.
  • Passes LG_BRUTAL tests
  • Code is simpler but performs as well as FastSV7. (In fact, it is consistently ~5% faster on GAP-twitter in my testing)
  • Piggybacked off of cc_demo.c and test_ConnectedComponents.c for testing and benchmarking.

LG_CC_FastSV7.c and LAGraph.h

  • Moved out the definition for USING_GRAPHBLAS_V10

LAGr_SwapEdges.c, LAGraph_SwapEdges.c, test_SwapEdges.c, and SwapEdges_demo.c

  • Randomizes undirected, unweighted Graphs.
  • LAGr allows the user to customize attempted swaps per edge, minimum swap count per edge, and the number of swaps to do.
  • LAGraph_SwapEdges only allows the user to provide the number of desired swaps pre-edge to perform on a Graph, with the other values predetermined.

GomezGab and others added 30 commits April 17, 2024 19:04
)) ;
LAGRAPH_TRY (LAGraph_New (G_new, &A_new, LAGraph_ADJACENCY_DIRECTED, msg)) ;
LG_FREE_WORK ;
return (num_swaps >= totSwaps)? GrB_SUCCESS : LAGRAPH_INSUFFICIENT_SWAPS ;
Copy link
Collaborator Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I defined LAGRAPH_INSUFFICIENT_SWAPS as 2100 in LAGraphX.h. Was this correct?
Should I return the number of swaps I made? As its own output variable, or do I make totSwaps input/output?

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I'm reluctant to add a new enum value to the LAGraph return values that is unique to just one algorithm. How about returning LAGRAPH_CONVERGENCE_FAILURE instead? The method iterates, trying to find swaps, and if it doesn't succeed, it does need to return an error.

The description of LAGRAPH_CONVERGENCE_FAILURE is quite similar, "An iterative process failed to converge to a good solution".

What do you think?

Copy link
Collaborator Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I think returning an error is reasonable. Although currently, I return the error, but still return a valid G_new graph with however many swaps I managed to get before going under the per-loop threshold. Is this a reasonable thing to do, or does returning an error imply the output graph should be freed?

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I see. Maybe we should talk it over. In your current design, LAGRAPH_INSUFFICIENT_SWAPS is a warning, not an error, right?

@GomezGab GomezGab marked this pull request as ready for review May 18, 2025 22:55
#if GxB_IMPLEMENTATION >= GxB_VERSION (10,0,0)
#define USING_GRAPHBLAS_V10 1
#else
#define USING_GRAPHBLAS_V10 0

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Since this #define is internal, but accessible to any LAGraph algorithm, let's rename it to "LG_SUITESPARSE_GRAPHBLAS_V10". That way, it mimics the LAGRAPH_SUITESPARSE definition (which is in LAGraph.h), but it has an LG_* prefix which means it's for internal use only.

@DrTimothyAldenDavis
Copy link
Member

Is this ready to merge?

@GomezGab
Copy link
Collaborator Author

Ready to go

@DrTimothyAldenDavis DrTimothyAldenDavis merged commit 5eee3cd into GraphBLAS:with_GraphBLAS_v10 May 23, 2025
4 checks passed
@DrTimothyAldenDavis
Copy link
Member

Looks great. Merged into the with_GraphBLAS_v10 branch.

@GomezGab GomezGab deleted the RCC_GB10 branch May 23, 2025 22:36
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Labels

None yet

Projects

None yet

Development

Successfully merging this pull request may close these issues.

2 participants