COMBINATORIAL_BLAS 1.6
 
Loading...
Searching...
No Matches
harmonic.cpp
Go to the documentation of this file.
1/*
2//@HEADER
3// *****************************************************************************
4//
5// HPCGraph: Graph Computation on High Performance Computing Systems
6// Copyright (2016) Sandia Corporation
7//
8// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9// the U.S. Government retains certain rights in this software.
10//
11// Redistribution and use in source and binary forms, with or without
12// modification, are permitted provided that the following conditions are
13// met:
14//
15// 1. Redistributions of source code must retain the above copyright
16// notice, this list of conditions and the following disclaimer.
17//
18// 2. Redistributions in binary form must reproduce the above copyright
19// notice, this list of conditions and the following disclaimer in the
20// documentation and/or other materials provided with the distribution.
21//
22// 3. Neither the name of the Corporation nor the names of the
23// contributors may be used to endorse or promote products derived from
24// this software without specific prior written permission.
25//
26// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact George M. Slota (gmslota@sandia.gov)
39// Siva Rajamanickam (srajama@sandia.gov)
40// Kamesh Madduri (madduri@cse.psu.edu)
41//
42// *****************************************************************************
43//@HEADER
44*/
45
46
47#include <mpi.h>
48#include <omp.h>
49#include <stdio.h>
50#include <stdlib.h>
51#include <stdint.h>
52#include <fstream>
53
54#include "dist_graph.h"
55#include "comms.h"
56#include "util.h"
57#include "harmonic.h"
58
59#define HC_NOT_VISITED 18446744073709551615
60#define HC_VISITED 18446744073709551614
61
62extern int procid, nprocs;
63extern bool verbose, debug, verify;
64
65
67 uint64_t* distances, uint64_t root)
68{
69 if (debug) { printf("Task %d run_harmonic() start\n", procid); }
70 double elt = 0.0;
71 if (verbose) {
72 MPI_Barrier(MPI_COMM_WORLD);
73 elt = omp_get_wtime();
74 }
75
76 q->queue_size = 0;
77 q->next_size = 0;
78 q->send_size = 0;
79
80 uint64_t root_index = get_value(&g->map, root);
81 if (root_index != NULL_KEY && root_index < g->n_local)
82 {
83 q->queue[0] = root;
84 q->queue_size = 1;
85 }
86
87 uint64_t distance = 0;
88 comm->global_queue_size = 1;
89#pragma omp parallel default(shared)
90{
93
94#pragma omp for
95 for (uint64_t i = 0; i < g->n_total; ++i)
96 distances[i] = HC_NOT_VISITED;
97
98 while (comm->global_queue_size)
99 {
100 if (debug && tq.tid == 0) {
101 printf("Task %d Distance %u run_harmonic() GQ: %li, TQ: %u\n",
102 procid, distance, comm->global_queue_size, q->queue_size);
103 }
104
105#pragma omp for schedule(guided) nowait
106 for (int64_t i = 0; i < q->queue_size; ++i)
107 {
108 uint64_t vert = q->queue[i];
109 uint64_t vert_index = get_value(&g->map, vert);
110 if (distances[vert_index] != HC_NOT_VISITED &&
111 distances[vert_index] != HC_VISITED)
112 continue;
113 distances[vert_index] = distance;
114
115 uint64_t in_degree = in_degree(g, vert_index);
116 uint64_t* ins = in_vertices(g, vert_index);
117 for (uint64_t j = 0; j < in_degree; ++j)
118 {
119 uint64_t in_index = ins[j];
120 if (distances[in_index] == HC_NOT_VISITED)
121 {
122 distances[in_index] = HC_VISITED;
123
124 if (in_index < g->n_local)
125 add_vid_to_queue(&tq, q, g->local_unmap[in_index]);
126 else
127 add_vid_to_send(&tq, q, in_index);
128 }
129 }
130 }
131
132 empty_queue(&tq, q);
133 empty_send(&tq, q);
134#pragma omp barrier
135
136#pragma omp single
137 {
138 exchange_verts(g, comm, q);
139 ++distance;
140 }
141 } // end while
142
144} // end parallel
145
146 if (verbose) {
147 elt = omp_get_wtime() - elt;
148 printf("Task %d run_harmonic() time %9.6f (s)\n", procid, elt);
149 }
150 if (debug) { printf("Task %d run_harmonic() success\n", procid); }
151
152 return 0;
153}
154
156 uint64_t* input_list, double* centralities,
157 char* output_file)
158{
159 if (verbose) printf("Task %d centralities to %s\n", procid, output_file);
160
161 if (procid == 0)
162 {
163 std::ofstream outfile;
164 outfile.open(output_file);
165
166 outfile << "Vertex, Centrality" << std::endl;
167
168 for (uint64_t i = 0; i < num_to_output; ++i)
169 outfile << input_list[i] << "," << centralities[i] << std::endl;
170
171 outfile.close();
172 }
173
174 if (verbose) printf("Task %d done writing centralities\n", procid);
175
176 return 0;
177}
178
179double harmonic_calc(dist_graph_t* g, uint64_t* distances, uint64_t root)
180{
181 double my_hc = 0.0;
182 uint64_t count = 0;
183 for (uint64_t i = 0; i < g->n_local; ++i)
184 if (distances[i] > 0 &&
185 distances[i] != HC_VISITED && distances[i] != HC_NOT_VISITED)
186 {
187 my_hc += (1.0 / (double)distances[i]);
188 ++count;
189 }
190
191 double global_hc = 0.0;
192 MPI_Allreduce(&my_hc, &global_hc, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
193
194 if (procid == 0)
195 printf("Vertex %lu - Harmonic Centrality %lf\n", root, global_hc);
196
197 return global_hc;
198}
199
201 char* output_file,
202 uint64_t num_to_output, uint64_t* input_list)
203{
204 if (debug) { printf("Task %d harmonic_dist() start\n", procid); }
205
206 MPI_Barrier(MPI_COMM_WORLD);
207 double elt = omp_get_wtime();
208
209 uint64_t* distances = (uint64_t*)malloc(g->n_total*sizeof(uint64_t));
210 double* centralities = (double*)malloc(num_to_output*sizeof(double));
211
212 for (uint64_t i = 0; i < num_to_output; ++i)
213 {
214 uint64_t root = input_list[i];
215 run_harmonic(g, comm, q, distances, root);
216 centralities[i] = harmonic_calc(g, distances, root);
217 }
218
219 MPI_Barrier(MPI_COMM_WORLD);
220 elt = omp_get_wtime() - elt;
221 if (procid == 0) printf("Harmonic time %9.6f (s)\n", elt);
222
223 if (output) {
224 harmonic_output(g, num_to_output, input_list, centralities, output_file);
225 }
226
227 free(distances);
228
229 if (debug) printf("Task %d harmonic_dist() success\n", procid);
230 return 0;
231}
232
Mac OS X ATTR com apple quarantine q
Definition ._remapper.cpp:1
void init_thread_queue(thread_queue_t *tq)
Definition comms.cpp:136
void clear_thread_queue(thread_queue_t *tq)
Definition comms.cpp:157
bool output
Definition comms.cpp:56
void empty_queue(thread_queue_t *tq, queue_data_t *q)
Definition comms.h:730
void add_vid_to_queue(thread_queue_t *tq, queue_data_t *q, uint64_t vertex_id)
Definition comms.h:721
void exchange_verts(dist_graph_t *g, mpi_data_t *comm, queue_data_t *q)
Definition comms.h:184
void add_vid_to_send(thread_queue_t *tq, queue_data_t *q, uint64_t vertex_id)
Definition comms.h:743
void empty_send(thread_queue_t *tq, queue_data_t *q)
Definition comms.h:752
long int64_t
Definition compat.h:21
#define in_degree(g, n)
Definition dist_graph.h:55
#define in_vertices(g, n)
Definition dist_graph.h:57
#define NULL_KEY
Definition fast_map.h:58
uint64_t get_value(fast_map *map, uint64_t key)
Definition fast_map.h:132
int harmonic_dist(dist_graph_t *g, mpi_data_t *comm, queue_data_t *q, char *output_file, uint64_t num_to_output, uint64_t *input_list)
Definition harmonic.cpp:200
int procid
Definition main.cpp:55
bool debug
Definition harmonic.cpp:63
#define HC_NOT_VISITED
Definition harmonic.cpp:59
bool verify
Definition harmonic.cpp:63
int harmonic_output(dist_graph_t *g, uint64_t num_to_output, uint64_t *input_list, double *centralities, char *output_file)
Definition harmonic.cpp:155
double harmonic_calc(dist_graph_t *g, uint64_t *distances, uint64_t root)
Definition harmonic.cpp:179
#define HC_VISITED
Definition harmonic.cpp:60
bool verbose
Definition main.cpp:56
int nprocs
Definition harmonic.cpp:62
int run_harmonic(dist_graph_t *g, mpi_data_t *comm, queue_data_t *q, uint64_t *distances, uint64_t root)
Definition harmonic.cpp:66
unsigned __int64 uint64_t
Definition stdint.h:90
uint64_t n_total
Definition dist_graph.h:69
uint64_t * local_unmap
Definition dist_graph.h:80
uint64_t n_local
Definition dist_graph.h:66
fast_map map
Definition dist_graph.h:83
uint64_t global_queue_size
Definition comms.h:87
int32_t tid
Definition comms.h:101