33#ifndef __STDC_CONSTANT_MACROS
34#define __STDC_CONSTANT_MACROS
36#ifndef __STDC_LIMIT_MACROS
37#define __STDC_LIMIT_MACROS
48#include "CombBLAS/CombBLAS.h"
56template <
typename T1,
typename T2>
60 static T_promote id(){
return std::numeric_limits<T_promote>::max(); };
82template <
typename IT,
typename NT,
typename DER>
94 [](
IT pv,
IT gpv) ->
short {
return static_cast<short>(
pv ==
gpv); },
103 false,
static_cast<short>(0));
112 false,
static_cast<IT>(0));
130template <
typename IT,
typename NT,
typename DER>
140 [](std::pair<IT, IT>
x,
short isStar){
return std::make_pair(0,0);},
148 [](std::pair<IT, IT>
x,
IT f){
return std::make_pair(
f,0);},
149 [](std::pair<IT, IT>
x,
IT f){
return true;},
154 [](std::pair<IT, IT>
x,
IT mnf){
return std::make_pair(std::get<0>(
x),
mnf);},
155 [](std::pair<IT, IT>
x,
IT mnf){
return std::get<0>(
x) >
mnf;},
159 [](std::pair<IT, IT> val,
IT ind){
return std::get<0>(val);},
160 [](std::pair<IT, IT> val,
IT ind){
return std::make_pair(ind, std::get<1>(val));},
161 [](std::pair<IT, IT>
val1, std::pair<IT, IT>
val2){
return val2;} );
168 [](std::pair<IT, IT>
x,
IT f){
return std::get<1>(
x);},
169 [](std::pair<IT, IT>
x,
IT f){
return true;},
174template <
typename IT,
typename NT,
typename DER>
185 [](std::pair<IT, IT>
x,
short isStar){
return std::make_pair(0,0);},
193 [](std::pair<IT, IT>
x,
IT f){
return std::make_pair(
f,0);},
194 [](std::pair<IT, IT>
x,
IT f){
return true;},
199 [](std::pair<IT, IT>
x,
IT mnf){
return std::make_pair(std::get<0>(
x),
mnf);},
200 [](std::pair<IT, IT>
x,
IT mnf){
return std::get<0>(
x) !=
mnf;},
204 [](std::pair<IT, IT> val,
IT ind){
return std::get<0>(val);},
205 [](std::pair<IT, IT> val,
IT ind){
return std::make_pair(ind, std::get<1>(val));},
206 [](std::pair<IT, IT>
val1, std::pair<IT, IT>
val2){
return val1;} );
211 [](std::pair<IT, IT>
x,
IT f){
return std::get<1>(
x);},
212 [](std::pair<IT, IT>
x,
IT f){
return true;},
217template <
typename IT>
224template <
typename IT,
typename NT,
typename DER>
234template <
typename IT,
typename NT,
typename DER>
237 DER* spSeq =
A.seqptr();
261template <
typename IT>
265 cclabel.ApplyInd([](
IT val,
IT ind){
return val==ind ? -1 : val;});
270 return roots.getnnz();
275template <
typename IT,
typename NT,
typename DER>
283 std::ostringstream
outs;
322template <
typename IT>
325 for(
IT i=0; i<
nCC; i++)
333template <
typename IT>
341 std::ostringstream
outs;
344 outs <<
"Size of the first component: " <<
cc1.getnnz() << std::endl;
345 outs <<
"Size of the second component: " <<
cc2.getnnz() << std::endl;
346 outs <<
"Size of the third component: " <<
cc3.getnnz() << std::endl;
347 outs <<
"Size of the fourth component: " <<
cc4.getnnz() << std::endl;
351template <
typename IT>
355 for(
IT i=0; i<
nCC; i++)
367 for(
IT i=0; i<
ccSizes.LocArrSize(); i++)
382 std::cout <<
"The largest component size: " <<
largestCCSise << std::endl;
384 output.open(
"hist.txt", std::ios_base::app );
static void Print(const std::string &s)
FullyDistVec< IT, short > StarCheck(const SpParMat< IT, NT, DER > &A, FullyDistVec< IT, IT > &father)
IT LabelCC(FullyDistVec< IT, IT > &father, FullyDistVec< IT, IT > &cclabel)
void Correctness(const SpParMat< IT, NT, DER > &A, FullyDistVec< IT, IT > &cclabel, IT nCC, FullyDistVec< IT, IT > father)
void UnconditionalHook(const SpParMat< IT, NT, DER > &A, FullyDistVec< IT, IT > &father)
bool neigborsInSameCC(const SpParMat< IT, NT, DER > &A, FullyDistVec< IT, IT > &cclabel)
void First4Clust(FullyDistVec< IT, IT > &cc)
void Shortcut(FullyDistVec< IT, IT > &father)
FullyDistVec< IT, IT > CC(SpParMat< IT, NT, DER > &A, IT &nCC)
void PrintCC(FullyDistVec< IT, IT > CC, IT nCC)
void HistCC(FullyDistVec< IT, IT > CC, IT nCC)
void ConditionalHook(const SpParMat< IT, NT, DER > &A, FullyDistVec< IT, IT > &father)
static bool returnedSAID()
promote_trait< T1, T2 >::T_promote T_promote
static T_promote add(const T_promote &arg1, const T_promote &arg2)
static T_promote multiply(const T1 &arg1, const T2 &arg2)
static T2 add(const T2 &arg1, const T2 &arg2)
static T2 multiply(const T1 &arg1, const T2 &arg2)
static void axpy(const T1 a, const T2 &x, T_promote &y)