6 #include "../CombBLAS.h"
24 TwitterEdge(
short mycount,
bool myfollow, time_t mylatest):count(mycount), follower(myfollow), latest(mylatest) {};
27 bool TweetWithinInterval (time_t begin, time_t end)
const {
return ((count > 0) && (begin <= latest && latest <= end)); };
28 bool TweetSince (time_t begin)
const {
return ((count > 0) && (begin <= latest)); };
29 bool LastTweetBy (time_t end)
const {
return ((count > 0) && (latest <= end)); };
31 operator bool ()
const {
return true; } ;
35 cout <<
"Error: TwitterEdge::operator+=() shouldn't be executed" << endl;
37 follower |= rhs.follower;
39 latest = max(latest, rhs.latest);
44 return ((follower == b.follower) && (latest == b.latest) && (count == b.count));
55 template <
typename IT>
63 if( twe.follower == 0 && twe.latest == 0 && twe.count == 0)
85 MPI_Datatype getMPIType()
87 return MPIType<TwitterEdge>();
90 void binaryfill(FILE * rFile, IT & row, IT & col,
TwitterEdge & val)
92 TwitterInteraction twi;
93 size_t entryLength = fread (&twi,
sizeof(TwitterInteraction),1,rFile);
96 val =
TwitterEdge(twi.retweets, twi.follow, twi.twtime);
98 cout <<
"Not enough bytes read in binaryfill " << endl;
102 template <
typename c,
typename t>
116 int year, month, day, hour, min, sec;
117 sscanf (date.c_str(),
"%d-%d-%d",&year, &month, &day);
118 sscanf (time.c_str(),
"%d:%d:%d",&hour, &min, &sec);
120 memset(&timeinfo, 0,
sizeof(
struct tm));
121 timeinfo.tm_year = year - 1900;
122 timeinfo.tm_mon = month - 1 ;
123 timeinfo.tm_mday = day;
124 timeinfo.tm_hour = hour;
125 timeinfo.tm_min = min;
126 timeinfo.tm_sec = sec;
127 tw.latest = timegm(&timeinfo);
128 if(tw.latest == -1) { cout <<
"Can not parse time date" << endl; exit(-1);}
139 template <
typename c,
typename t>
140 void save(std::basic_ostream<c,t>& os,
const TwitterEdge & tw, IT row, IT col)
142 os << row <<
"\t" << col <<
"\t";
143 os << tw.follower <<
"\t";
144 os << tw.count <<
"\t";
145 os << tw.latest << endl;
148 struct TwitterInteraction
166 return (
id == rhs.
id);
170 return (
id != rhs.
id);
174 cout <<
"Adding parent with id: " << rhs.
id <<
" to this one with id " <<
id << endl;
185 template <
typename IT>
191 os <<
"Parent=" << twe.
id;
200 template <
typename IT>
207 template <
typename SR,
typename T>
208 void select2nd(
void * invec,
void * inoutvec,
int * len, MPI_Datatype *datatype);
211 template <
typename SR,
typename VECTYPE>
212 static VECTYPE filtered_select2nd(
const TwitterEdge & arg1,
const VECTYPE & arg2, time_t & sincedate)
217 memset(&timeinfo, 0,
sizeof(
struct tm));
218 int year, month, day, hour, min, sec;
219 year = 2009; month = 7; day = 1;
220 hour = 0; min = 0; sec = 0;
222 timeinfo.tm_year = year - 1900;
223 timeinfo.tm_mon = month - 1 ;
224 timeinfo.tm_mday = day;
225 timeinfo.tm_hour = hour;
226 timeinfo.tm_min = min;
227 timeinfo.tm_sec = sec;
228 sincedate = timegm(&timeinfo);
231 outs <<
"Initializing since date (only once) to " << sincedate << endl;
263 static bool returnedSAID(
bool setFlagTo =
false)
265 static bool flag =
false;
277 static MPI_Op mpi_op()
279 MPI_Op_create(select2nd<LatestRetwitterBFS,ParentType>,
false, &MPI_BFSADD);
285 return filtered_select2nd<LatestRetwitterBFS>(arg1, arg2, sincedate);
289 y = add(y, multiply(a, x));
296 template <
typename SR,
typename T>
297 void select2nd(
void * invec,
void * inoutvec,
int * len, MPI_Datatype *datatype)
299 T * pinvec =
static_cast<T*
>(invec);
300 T * pinoutvec =
static_cast<T*
>(inoutvec);
301 for (
int i = 0; i < *len; i++)
303 pinoutvec[i] =
SR::add(pinvec[i], pinoutvec[i]);
310 struct getfringe:
public std::binary_function<ParentType, ParentType, ParentType>
321 struct seldegree:
public std::binary_function<ParentType, int64_t, int64_t>
331 struct passifthere:
public std::binary_function<ParentType, int64_t, bool>
343 bool operator()(
double m,
double c)
const
352 uint8_t operator() (
double t1,
double t2)
385 static double id() {
return 0.0; }
390 static bool returnedSAID(
bool setFlagTo =
false)
392 static bool flag =
false;
399 static double add(
const double & arg1,
const double & arg2)
401 return std::min(arg1, arg2);
404 static MPI_Op mpi_op()
409 static double multiply(
const TwitterEdge & arg1,
const double & arg2)
411 return filtered_select2nd<LatestRetwitterMIS>(arg1, arg2, sincedate);
415 y = add(y, multiply(a, x));
424 static double id() {
return 0.0; }
429 static bool returnedSAID(
bool setFlagTo =
false)
431 static bool flag =
false;
438 static double add(
const double & arg1,
const double & arg2)
443 static MPI_Op mpi_op()
445 MPI_Op_create(select2nd<LatestRetwitterSelect2nd,double>,
false, &MPI_SEL2NDADD);
446 return MPI_SEL2NDADD;
449 static double multiply(
const TwitterEdge & arg1,
const double & arg2)
451 return filtered_select2nd<LatestRetwitterSelect2nd>(arg1, arg2, sincedate);
455 y = add(y, multiply(a, x));