1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568
|
/*
* The MIT License
*
* Wavefront Alignment Algorithms
* Copyright (c) 2017 by Santiago Marco-Sola <santiagomsola@gmail.com>
*
* This file is part of Wavefront Alignment Algorithms.
*
* Permission is hereby granted, free of charge, to any person obtaining a copy
* of this software and associated documentation files (the "Software"), to deal
* in the Software without restriction, including without limitation the rights
* to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
* copies of the Software, and to permit persons to whom the Software is
* furnished to do so, subject to the following conditions:
*
* The above copyright notice and this permission notice shall be included in all
* copies or substantial portions of the Software.
*
* THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
* IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
* FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
* AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
* LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
* OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
* SOFTWARE.
*
* PROJECT: Wavefront Alignment Algorithms
* AUTHOR(S): Santiago Marco-Sola <santiagomsola@gmail.com>
* DESCRIPTION: Wavefront Alignment Algorithms benchmarking tool
*/
#include <omp.h>
#include "align_benchmark_params.h"
#include "utils/commons.h"
#include "utils/sequence_buffer.h"
#include "system/profiler_timer.h"
#include "alignment/score_matrix.h"
#include "edit/edit_dp.h"
#include "gap_linear/nw.h"
#include "gap_affine/swg.h"
#include "gap_affine2p/affine2p_matrix.h"
#include "gap_affine2p/affine2p_dp.h"
#include "wavefront/wavefront_align.h"
#include "benchmark/benchmark_indel.h"
#include "benchmark/benchmark_edit.h"
#include "benchmark/benchmark_gap_linear.h"
#include "benchmark/benchmark_gap_affine.h"
#include "benchmark/benchmark_gap_affine2p.h"
/*
* Algorithms
*/
bool align_benchmark_is_wavefront(
const alignment_algorithm_type algorithm) {
return algorithm == alignment_indel_wavefront ||
algorithm == alignment_edit_wavefront ||
algorithm == alignment_gap_linear_wavefront ||
algorithm == alignment_gap_affine_wavefront ||
algorithm == alignment_gap_affine2p_wavefront;
}
/*
* Benchmark UTest
*/
void align_pairwise_test() {
// Patters & Texts
char * pattern = "GATTACA";
char * text = "GATCACTA";
// Penalties
linear_penalties_t linear_penalties = {
.match = 0,
.mismatch = 4,
.indel = 2,
};
affine_penalties_t affine_penalties = {
.match = 0,
.mismatch = 4, //9,
.gap_opening = 6, //13,
.gap_extension = 2,
};
// Ends
const int pattern_begin_free = 0;
const int pattern_end_free = 0;
const int text_begin_free = 0;
const int text_end_free = 0;
const bool endsfree =
pattern_begin_free>0 || pattern_end_free>0 ||
text_begin_free>0 || text_end_free>0;
/*
* Gap-Affine
*/
// Allocate
wavefront_aligner_attr_t attributes = wavefront_aligner_attr_default;
attributes.distance_metric = gap_affine;
attributes.linear_penalties = linear_penalties;
attributes.affine_penalties = affine_penalties;
attributes.heuristic.strategy = wf_heuristic_none;
attributes.heuristic.min_wavefront_length = 256;
attributes.heuristic.max_distance_threshold = 4096;
attributes.heuristic.steps_between_cutoffs = 10;
attributes.alignment_scope = compute_alignment; // compute_score
attributes.memory_mode = wavefront_memory_med;
attributes.alignment_form.span = (endsfree) ? alignment_endsfree : alignment_end2end;
attributes.alignment_form.pattern_begin_free = pattern_begin_free;
attributes.alignment_form.pattern_end_free = pattern_end_free;
attributes.alignment_form.text_begin_free = text_begin_free;
attributes.alignment_form.text_end_free = text_end_free;
attributes.plot.enabled = false;
wavefront_aligner_t* const wf_aligner = wavefront_aligner_new(&attributes);
// Align
wavefront_align(wf_aligner,
pattern,strlen(pattern),text,strlen(text));
// CIGAR
fprintf(stderr,">> WFA2");
cigar_print_pretty(stderr,
pattern,strlen(pattern),text,strlen(text),
wf_aligner->cigar,wf_aligner->mm_allocator);
fprintf(stderr,"SCORE: %d \n",cigar_score_gap_affine(
wf_aligner->cigar,&affine_penalties));
// Plot
if (attributes.plot.enabled) {
FILE* const wf_plot = fopen("test.wfa","w");
wavefront_plot_print(wf_plot,wf_aligner);
fclose(wf_plot);
}
// Free
wavefront_aligner_delete(wf_aligner);
}
/*
* Configuration
*/
wavefront_aligner_t* align_input_configure_wavefront(
align_input_t* const align_input) {
// Set attributes
wavefront_aligner_attr_t attributes = wavefront_aligner_attr_default;
attributes.memory_mode = parameters.wfa_memory_mode;
if (parameters.wfa_score_only) {
attributes.alignment_scope = compute_score;
}
// WF-Heuristic
switch (parameters.wfa_heuristic) {
case wf_heuristic_none:
attributes.heuristic.strategy = wf_heuristic_none;
break;
case wf_heuristic_banded_static:
attributes.heuristic.strategy = wf_heuristic_banded_static;
attributes.heuristic.min_k = parameters.wfa_heuristic_p1;
attributes.heuristic.max_k = parameters.wfa_heuristic_p2;
break;
case wf_heuristic_banded_adaptive:
attributes.heuristic.strategy = wf_heuristic_banded_adaptive;
attributes.heuristic.min_k = parameters.wfa_heuristic_p1;
attributes.heuristic.max_k = parameters.wfa_heuristic_p2;
attributes.heuristic.steps_between_cutoffs = parameters.wfa_heuristic_p3;
break;
case wf_heuristic_wfadaptive:
attributes.heuristic.strategy = wf_heuristic_wfadaptive;
attributes.heuristic.min_wavefront_length = parameters.wfa_heuristic_p1;
attributes.heuristic.max_distance_threshold = parameters.wfa_heuristic_p2;
attributes.heuristic.steps_between_cutoffs = parameters.wfa_heuristic_p3;
break;
case wf_heuristic_xdrop:
attributes.heuristic.strategy = wf_heuristic_xdrop;
attributes.heuristic.xdrop = parameters.wfa_heuristic_p1;
attributes.heuristic.steps_between_cutoffs = parameters.wfa_heuristic_p2;
break;
case wf_heuristic_zdrop:
attributes.heuristic.strategy = wf_heuristic_zdrop;
attributes.heuristic.zdrop = parameters.wfa_heuristic_p1;
attributes.heuristic.steps_between_cutoffs = parameters.wfa_heuristic_p2;
break;
default:
break;
}
// Select flavor
switch (parameters.algorithm) {
case alignment_indel_wavefront:
attributes.distance_metric = indel;
break;
case alignment_edit_wavefront:
attributes.distance_metric = edit;
break;
case alignment_gap_linear_wavefront:
attributes.distance_metric = gap_linear;
attributes.linear_penalties = parameters.linear_penalties;
break;
case alignment_gap_affine_wavefront:
attributes.distance_metric = gap_affine;
attributes.affine_penalties = parameters.affine_penalties;
break;
case alignment_gap_affine2p_wavefront:
attributes.distance_metric = gap_affine_2p;
attributes.affine2p_penalties = parameters.affine2p_penalties;
break;
default:
return NULL; // No WF selected
break;
}
// Select alignment form
attributes.alignment_form.span = (parameters.endsfree) ? alignment_endsfree : alignment_end2end;
// Misc
if (parameters.wfa_match_funct_arguments != NULL) {
attributes.match_funct = parameters.wfa_match_funct;
attributes.match_funct_arguments = parameters.wfa_match_funct_arguments;
}
attributes.plot.enabled = (parameters.plot != 0);
attributes.plot.align_level = (parameters.plot < 0) ? -1 : parameters.plot - 1;
attributes.system.verbose = parameters.verbose;
attributes.system.max_memory_abort = parameters.wfa_max_memory;
attributes.system.max_alignment_score = parameters.wfa_max_score;
attributes.system.max_num_threads = parameters.wfa_max_threads;
// Allocate
return wavefront_aligner_new(&attributes);
}
void align_input_configure_global(
align_input_t* const align_input) {
// Clear
benchmark_align_input_clear(align_input);
// Penalties
align_input->linear_penalties = parameters.linear_penalties;
align_input->affine_penalties = parameters.affine_penalties;
align_input->affine2p_penalties = parameters.affine2p_penalties;
// Alignment form
align_input->ends_free = parameters.endsfree;
// Output
align_input->output_file = parameters.output_file;
align_input->output_full = parameters.output_full;
// MM
align_input->mm_allocator = mm_allocator_new(BUFFER_SIZE_1M);
// WFA
if (align_benchmark_is_wavefront(parameters.algorithm)) {
align_input->wf_aligner = align_input_configure_wavefront(align_input);
} else {
align_input->wf_aligner = NULL;
}
// PROFILE/STATS
timer_reset(&align_input->timer);
// DEBUG
align_input->debug_flags = 0;
align_input->debug_flags |= parameters.check_metric;
if (parameters.check_display) align_input->debug_flags |= ALIGN_DEBUG_DISPLAY_INFO;
if (parameters.check_correct) align_input->debug_flags |= ALIGN_DEBUG_CHECK_CORRECT;
if (parameters.check_score) align_input->debug_flags |= ALIGN_DEBUG_CHECK_SCORE;
if (parameters.check_alignments) align_input->debug_flags |= ALIGN_DEBUG_CHECK_ALIGNMENT;
align_input->check_linear_penalties = ¶meters.linear_penalties;
align_input->check_affine_penalties = ¶meters.affine_penalties;
align_input->check_affine2p_penalties = ¶meters.affine2p_penalties;
align_input->check_bandwidth = parameters.check_bandwidth;
align_input->verbose = parameters.verbose;
}
void align_input_configure_local(
align_input_t* const align_input) {
// Ends-free configuration
if (parameters.endsfree) {
const int plen = align_input->pattern_length;
const int tlen = align_input->text_length;
align_input->pattern_begin_free = nominal_prop_u32(plen,parameters.pattern_begin_free);
align_input->pattern_end_free = nominal_prop_u32(plen,parameters.pattern_end_free);
align_input->text_begin_free = nominal_prop_u32(tlen,parameters.text_begin_free);
align_input->text_end_free = nominal_prop_u32(tlen,parameters.text_end_free);
if (align_benchmark_is_wavefront(parameters.algorithm)) {
wavefront_aligner_set_alignment_free_ends(align_input->wf_aligner,
align_input->pattern_begin_free,align_input->pattern_end_free,
align_input->text_begin_free,align_input->text_end_free);
}
}
// Custom extend-match function
if (parameters.wfa_match_funct != NULL) {
match_function_params.pattern = align_input->pattern;
match_function_params.pattern_length = align_input->pattern_length;
match_function_params.text = align_input->text;
match_function_params.text_length = align_input->text_length;
}
}
void align_benchmark_free(
align_input_t* const align_input) {
if (align_input->wf_aligner) wavefront_aligner_delete(align_input->wf_aligner);
mm_allocator_delete(align_input->mm_allocator);
}
/*
* I/O
*/
bool align_benchmark_read_input(
FILE* input_file,
char** line1,
char** line2,
size_t* line1_allocated,
size_t* line2_allocated,
const int seqs_processed,
align_input_t* const align_input) {
// Parameters
int line1_length=0, line2_length=0;
// Read queries
line1_length = getline(line1,line1_allocated,input_file);
if (line1_length==-1) return false;
line2_length = getline(line2,line2_allocated,input_file);
if (line1_length==-1) return false;
// Configure input
align_input->sequence_id = seqs_processed;
align_input->pattern = *line1 + 1;
align_input->pattern_length = line1_length - 2;
align_input->pattern[align_input->pattern_length] = '\0';
align_input->text = *line2 + 1;
align_input->text_length = line2_length - 2;
align_input->text[align_input->text_length] = '\0';
return true;
}
/*
* Display
*/
void align_benchmark_print_progress(
const int seqs_processed) {
const uint64_t time_elapsed_alg = timer_get_current_total_ns(¶meters.timer_global);
const float rate_alg = (float)seqs_processed/(float)TIMER_CONVERT_NS_TO_S(time_elapsed_alg);
fprintf(stderr,"...processed %d reads (alignment = %2.3f seq/s)\n",seqs_processed,rate_alg);
}
void align_benchmark_print_results(
align_input_t* const align_input,
const int seqs_processed,
const bool print_stats) {
// Print benchmark results
fprintf(stderr,"[Benchmark]\n");
fprintf(stderr,"=> Total.reads %d\n",seqs_processed);
fprintf(stderr,"=> Time.Benchmark ");
timer_print(stderr,¶meters.timer_global,NULL);
if (parameters.num_threads == 1) {
fprintf(stderr," => Time.Alignment ");
timer_print(stderr,&align_input->timer,¶meters.timer_global);
} else {
for (int i=0;i<parameters.num_threads;++i) {
fprintf(stderr," => Time.Alignment.Thread.%0d ",i);
timer_print(stderr,&align_input[i].timer,¶meters.timer_global);
}
}
// Print Stats
const bool checks_enabled =
parameters.check_display || parameters.check_correct ||
parameters.check_score || parameters.check_alignments;
if (checks_enabled && parameters.num_threads==1) {
const bool print_wf_stats = (parameters.algorithm == alignment_gap_affine_wavefront);
benchmark_print_stats(stderr,align_input,print_wf_stats);
}
}
void align_benchmark_plot_wf(
align_input_t* const align_input,
const int seq_id) {
// Setup filename
char filename[500];
if (parameters.output_filename != NULL) {
sprintf(filename,"%s.%03d.plot",parameters.output_filename,seq_id);
} else {
sprintf(filename,"%s.%03d.plot",parameters.input_filename,seq_id);
}
// Open file
FILE* const wf_plot = fopen(filename,"w");
wavefront_plot_print(wf_plot,align_input->wf_aligner);
fclose(wf_plot);
}
/*
* Benchmark
*/
void align_benchmark_run_algorithm(
align_input_t* const align_input) {
// Sequence-dependent configuration
align_input_configure_local(align_input);
// Select algorithm
switch (parameters.algorithm) {
// Indel
case alignment_indel_wavefront:
benchmark_indel_wavefront(align_input);
break;
// Edit
case alignment_edit_bpm:
benchmark_edit_bpm(align_input);
break;
case alignment_edit_dp:
benchmark_edit_dp(align_input);
break;
case alignment_edit_dp_banded:
benchmark_edit_dp_banded(align_input,parameters.bandwidth);
break;
case alignment_edit_wavefront:
benchmark_edit_wavefront(align_input);
break;
// Gap-linear
case alignment_gap_linear_nw:
benchmark_gap_linear_nw(align_input,¶meters.linear_penalties);
break;
case alignment_gap_linear_wavefront:
benchmark_gap_linear_wavefront(align_input,¶meters.linear_penalties);
break;
// Gap-affine
case alignment_gap_affine_swg:
benchmark_gap_affine_swg(align_input,¶meters.affine_penalties);
break;
case alignment_gap_affine_swg_endsfree:
benchmark_gap_affine_swg_endsfree(
align_input,¶meters.affine_penalties);
break;
case alignment_gap_affine_swg_banded:
benchmark_gap_affine_swg_banded(align_input,
¶meters.affine_penalties,parameters.bandwidth);
break;
case alignment_gap_affine_wavefront:
benchmark_gap_affine_wavefront(align_input,¶meters.affine_penalties);
break;
// Gap-affine 2p
case alignment_gap_affine2p_dp:
benchmark_gap_affine2p_dp(align_input,¶meters.affine2p_penalties);
break;
case alignment_gap_affine2p_wavefront:
benchmark_gap_affine2p_wavefront(align_input,¶meters.affine2p_penalties);
break;
default:
fprintf(stderr,"Algorithm not implemented\n");
exit(1);
break;
}
}
void align_benchmark_sequential() {
// PROFILE
timer_reset(¶meters.timer_global);
timer_start(¶meters.timer_global);
// I/O files
parameters.input_file = fopen(parameters.input_filename, "r");
if (parameters.input_file == NULL) {
fprintf(stderr,"Input file '%s' couldn't be opened\n",parameters.input_filename);
exit(1);
}
if (parameters.output_filename != NULL) {
parameters.output_file = fopen(parameters.output_filename, "w");
}
// Global configuration
align_input_t align_input;
align_input_configure_global(&align_input);
// Read-align loop
int seqs_processed = 0, progress = 0;
while (true) {
// Read input sequence-pair
const bool input_read = align_benchmark_read_input(
parameters.input_file,¶meters.line1,¶meters.line2,
¶meters.line1_allocated,¶meters.line2_allocated,
seqs_processed,&align_input);
if (!input_read) break;
// Execute the selected algorithm
align_benchmark_run_algorithm(&align_input);
// Update progress
++seqs_processed;
if (++progress == parameters.progress) {
progress = 0;
if (parameters.verbose >= 0) align_benchmark_print_progress(seqs_processed);
}
// DEBUG
// mm_allocator_print(stderr,align_input.wf_aligner->mm_allocator,false);
// mm_allocator_print(stderr,align_input.wf_aligner->bialigner->alg_forward->mm_allocator,false);
// mm_allocator_print(stderr,align_input.wf_aligner->bialigner->alg_reverse->mm_allocator,false);
// mm_allocator_print(stderr,align_input.wf_aligner->bialigner->alg_subsidiary->mm_allocator,false);
// Plot
if (parameters.plot != 0) align_benchmark_plot_wf(&align_input,seqs_processed);
}
// Print benchmark results
timer_stop(¶meters.timer_global);
if (parameters.verbose >= 0) align_benchmark_print_results(&align_input,seqs_processed,true);
// Free
align_benchmark_free(&align_input);
fclose(parameters.input_file);
if (parameters.output_file) fclose(parameters.output_file);
free(parameters.line1);
free(parameters.line2);
}
void align_benchmark_parallel() {
// PROFILE
timer_reset(¶meters.timer_global);
timer_start(¶meters.timer_global);
// Open input file
parameters.input_file = fopen(parameters.input_filename, "r");
if (parameters.input_file == NULL) {
fprintf(stderr,"Input file '%s' couldn't be opened\n",parameters.input_filename);
exit(1);
}
if (parameters.output_filename != NULL) {
parameters.output_file = fopen(parameters.output_filename, "w");
}
// Global configuration
align_input_t align_input[parameters.num_threads];
for (int tid=0;tid<parameters.num_threads;++tid) {
align_input_configure_global(align_input+tid);
}
// Read-align loop
sequence_buffer_t* const sequence_buffer = sequence_buffer_new(2*parameters.batch_size,100);
int seqs_processed = 0, progress = 0, seqs_batch = 0;
while (true) {
// Read batch-input sequence-pair
sequence_buffer_clear(sequence_buffer);
for (seqs_batch=0;seqs_batch<parameters.batch_size;++seqs_batch) {
const bool seqs_pending = align_benchmark_read_input(
parameters.input_file,¶meters.line1,¶meters.line2,
¶meters.line1_allocated,¶meters.line2_allocated,
seqs_processed,align_input);
if (!seqs_pending) break;
// Add pair pattern-text
sequence_buffer_add_pair(sequence_buffer,
align_input->pattern,align_input->pattern_length,
align_input->text,align_input->text_length);
}
if (seqs_batch == 0) break;
// Parallel processing of the sequences batch
#pragma omp parallel num_threads(parameters.num_threads)
{
int tid = omp_get_thread_num();
#pragma omp for
for (int seq_idx=0;seq_idx<seqs_batch;++seq_idx) {
// Configure sequence
sequence_offset_t* const offset = sequence_buffer->offsets + seq_idx;
align_input[tid].sequence_id = seqs_processed;
align_input[tid].pattern = sequence_buffer->buffer + offset->pattern_offset;
align_input[tid].pattern_length = offset->pattern_length;
align_input[tid].text = sequence_buffer->buffer + offset->text_offset;
align_input[tid].text_length = offset->text_length;
// Execute the selected algorithm
align_benchmark_run_algorithm(align_input+tid);
}
}
// Update progress
seqs_processed += seqs_batch;
progress += seqs_batch;
if (progress >= parameters.progress) {
progress -= parameters.progress;
if (parameters.verbose >= 0) align_benchmark_print_progress(seqs_processed);
}
// DEBUG
//mm_allocator_print(stderr,align_input->wf_aligner->mm_allocator,false);
}
// Print benchmark results
timer_stop(¶meters.timer_global);
if (parameters.verbose >= 0) align_benchmark_print_results(align_input,seqs_processed,true);
// Free
for (int tid=0;tid<parameters.num_threads;++tid) {
align_benchmark_free(align_input+tid);
}
sequence_buffer_delete(sequence_buffer);
fclose(parameters.input_file);
if (parameters.output_file) fclose(parameters.output_file);
free(parameters.line1);
free(parameters.line2);
}
/*
* Main
*/
int main(int argc,char* argv[]) {
// Parsing command-line options
parse_arguments(argc,argv);
// Select option
if (parameters.algorithm == alignment_test) {
align_pairwise_test();
} else {
// Execute benchmark
if (parameters.num_threads == 1) {
align_benchmark_sequential();
} else {
align_benchmark_parallel();
}
}
}
|