Skip to content

Commit 8165243

Browse files
Add highlight feature for specified nodes in visualization
1 parent 497d046 commit 8165243

File tree

1 file changed

+62
-2
lines changed

1 file changed

+62
-2
lines changed

src/subcommand/viz_main.cpp

Lines changed: 62 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -12,6 +12,8 @@
1212
#include "utils.hpp"
1313
#include "colorbrewer.hpp"
1414
#include "split.hpp"
15+
#include <fstream>
16+
#include <unordered_set>
1517

1618
//#define debug_odgi_viz
1719

@@ -72,6 +74,7 @@ namespace odgi {
7274
{'N', "color-by-uncalled-bases"});
7375
args::ValueFlag<char> color_by_prefix(viz_opts, "CHAR", "Color paths by their names looking at the prefix before the given"
7476
" character CHAR.",{'s', "color-by-prefix"});
77+
args::ValueFlag<std::string> highlight_node_ids(viz_opts, "FILE", "Color nodes listed in FILE (one id per row) in red and all other nodes in grey.", {'J', "highlight-node-ids"});
7578
// TODO
7679
args::ValueFlag<std::string> _name_prefixes(viz_opts, "FILE", "Merge paths beginning with prefixes listed (one per line) in *FILE*.", {'M', "prefix-merges"});
7780
args::ValueFlag<std::string> _ignore_prefix(viz_opts, "PREFIX", "Ignore paths starting with the given *PREFIX*.", {'I', "ignore-prefix"});
@@ -237,6 +240,33 @@ namespace odgi {
237240
return 1;
238241
}
239242

243+
std::unordered_set<int64_t> highlighted_nodes;
244+
const bool highlight_nodes_requested = highlight_node_ids;
245+
if (highlight_nodes_requested) {
246+
std::ifstream ids_in(args::get(highlight_node_ids));
247+
if (!ids_in) {
248+
std::cerr << "[odgi::viz] error: unable to open file "
249+
<< args::get(highlight_node_ids) << " with node identifiers." << std::endl;
250+
return 1;
251+
}
252+
std::string line;
253+
while (std::getline(ids_in, line)) {
254+
const auto first = line.find_first_not_of(" \t\r\n");
255+
if (first == std::string::npos) {
256+
continue;
257+
}
258+
const auto last = line.find_last_not_of(" \t\r\n");
259+
const std::string trimmed = line.substr(first, last - first + 1);
260+
try {
261+
highlighted_nodes.insert(std::stoll(trimmed));
262+
} catch (const std::exception &) {
263+
std::cerr << "[odgi::viz] error: invalid node id '" << trimmed
264+
<< "' in " << args::get(highlight_node_ids) << "." << std::endl;
265+
return 1;
266+
}
267+
}
268+
}
269+
240270
const uint64_t num_threads = args::get(nthreads) ? args::get(nthreads) : 1;
241271

242272
//NOTE: this sample will overwrite the file or test.png without warning!
@@ -1259,6 +1289,11 @@ namespace odgi {
12591289

12601290
graph.for_each_step_in_path(path, [&](const step_handle_t &occ) {
12611291
h = graph.get_handle_of_step(occ);
1292+
const int64_t node_id = graph.get_id(h);
1293+
const bool node_highlighted = highlight_nodes_requested && highlighted_nodes.count(node_id);
1294+
const uint8_t highlight_r = node_highlighted ? 255 : 180;
1295+
const uint8_t highlight_g = node_highlighted ? 0 : 180;
1296+
const uint8_t highlight_b = node_highlighted ? 0 : 180;
12621297
p = position_map[number_bool_packing::unpack_number(h) - shift];
12631298
hl = graph.get_length(h);
12641299

@@ -1306,8 +1341,18 @@ namespace odgi {
13061341
}
13071342

13081343
if (curr_bin - 1 >= pangenomic_start_pos && curr_bin - 1 <= pangenomic_end_pos) {
1344+
const double draw_scale = highlight_nodes_requested ? 1.0 : x;
1345+
const uint8_t draw_r = highlight_nodes_requested
1346+
? highlight_r
1347+
: static_cast<uint8_t>(std::max(0.0, std::min(255.0, (double) path_r * draw_scale)));
1348+
const uint8_t draw_g = highlight_nodes_requested
1349+
? highlight_g
1350+
: static_cast<uint8_t>(std::max(0.0, std::min(255.0, (double) path_g * draw_scale)));
1351+
const uint8_t draw_b = highlight_nodes_requested
1352+
? highlight_b
1353+
: static_cast<uint8_t>(std::max(0.0, std::min(255.0, (double) path_b * draw_scale)));
13091354
add_path_step(image, width, curr_bin - 1 - pangenomic_start_pos, path_y,
1310-
(float) path_r * x, (float) path_g * x, (float) path_b * x);
1355+
draw_r, draw_g, draw_b);
13111356
}
13121357

13131358
}
@@ -1322,6 +1367,11 @@ namespace odgi {
13221367
/// Loop over all the steps along a path, from first through last and draw them
13231368
graph.for_each_step_in_path(path, [&](const step_handle_t &occ) {
13241369
handle_t h = graph.get_handle_of_step(occ);
1370+
const int64_t node_id = graph.get_id(h);
1371+
const bool node_highlighted = highlight_nodes_requested && highlighted_nodes.count(node_id);
1372+
const uint8_t highlight_r = node_highlighted ? 255 : 180;
1373+
const uint8_t highlight_g = node_highlighted ? 0 : 180;
1374+
const uint8_t highlight_b = node_highlighted ? 0 : 180;
13251375
uint64_t p = position_map[number_bool_packing::unpack_number(h) - shift];
13261376
uint64_t hl = graph.get_length(h);
13271377
// make contects for the bases in the node
@@ -1342,8 +1392,18 @@ namespace odgi {
13421392
}
13431393

13441394
if ((p + i) >= pangenomic_start_pos && (p + i) <= pangenomic_end_pos) {
1395+
const double draw_scale = highlight_nodes_requested ? 1.0 : x;
1396+
const uint8_t draw_r = highlight_nodes_requested
1397+
? highlight_r
1398+
: static_cast<uint8_t>(std::max(0.0, std::min(255.0, (double) path_r * draw_scale)));
1399+
const uint8_t draw_g = highlight_nodes_requested
1400+
? highlight_g
1401+
: static_cast<uint8_t>(std::max(0.0, std::min(255.0, (double) path_g * draw_scale)));
1402+
const uint8_t draw_b = highlight_nodes_requested
1403+
? highlight_b
1404+
: static_cast<uint8_t>(std::max(0.0, std::min(255.0, (double) path_b * draw_scale)));
13451405
add_path_step(image, width, p + i - pangenomic_start_pos, path_y,
1346-
(float) path_r * x, (float) path_g * x, (float) path_b * x);
1406+
draw_r, draw_g, draw_b);
13471407
}
13481408
}
13491409

0 commit comments

Comments
 (0)