File: WFAligner.cpp

package info (click to toggle)
libwfa2 2.3.3-4
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 10,072 kB
  • sloc: ansic: 13,812; python: 540; cpp: 500; makefile: 268; sh: 176; lisp: 41
file content (378 lines) | stat: -rw-r--r-- 12,704 bytes parent folder | download | duplicates (2)
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
/*
 *                             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: C++ bindings for the WaveFront Alignment modules
 */

#include "WFAligner.hpp"

extern "C" {
  #include "../../wavefront/wavefront_align.h"
}

/*
 * Namespace
 */
namespace wfa {

/*
 * General Wavefront Aligner
 */
WFAligner::WFAligner(
    const AlignmentScope alignmentScope,
    const MemoryModel memoryModel) {
  this->attributes = wavefront_aligner_attr_default;
  switch (memoryModel) {
    case MemoryHigh: this->attributes.memory_mode = wavefront_memory_high; break;
    case MemoryMed: this->attributes.memory_mode = wavefront_memory_med; break;
    case MemoryLow: this->attributes.memory_mode = wavefront_memory_low; break;
    case MemoryUltralow: this->attributes.memory_mode = wavefront_memory_ultralow; break;
    default: this->attributes.memory_mode = wavefront_memory_high; break;
  }
  this->attributes.alignment_scope = (alignmentScope==Score) ? compute_score : compute_alignment;
  //this->attributes.system.verbose = 2;
  this->wfAligner = nullptr;
}
WFAligner::~WFAligner() {
  wavefront_aligner_delete(wfAligner);
}
/*
 * Align End-to-end
 */
WFAligner::AlignmentStatus WFAligner::alignEnd2EndLambda(
    const int patternLength,
    const int textLength) {
  // Configure
  wavefront_aligner_set_alignment_end_to_end(wfAligner);
  // Align (using custom matching function)
  return (WFAligner::AlignmentStatus) wavefront_align(wfAligner,NULL,patternLength,NULL,textLength);
}
WFAligner::AlignmentStatus WFAligner::alignEnd2End(
    const char* const pattern,
    const int patternLength,
    const char* const text,
    const int textLength) {
  // Configure
  wavefront_aligner_set_alignment_end_to_end(wfAligner);
  // Align
  return (WFAligner::AlignmentStatus) wavefront_align(wfAligner,pattern,patternLength,text,textLength);
}
WFAligner::AlignmentStatus WFAligner::alignEnd2End(
    std::string& pattern,
    std::string& text) {
  // Delegate
  return alignEnd2End(pattern.c_str(),pattern.length(),text.c_str(),text.length());
}
/*
 * Align Ends-free
 */
WFAligner::AlignmentStatus WFAligner::alignEndsFreeLambda(
    const int patternLength,
    const int patternBeginFree,
    const int patternEndFree,
    const int textLength,
    const int textBeginFree,
    const int textEndFree) {
  // Configure
  wavefront_aligner_set_alignment_free_ends(wfAligner,
      patternBeginFree,patternEndFree,
      textBeginFree,textEndFree);
  // Align (using custom matching function)
  return (WFAligner::AlignmentStatus) wavefront_align(wfAligner,NULL,patternLength,NULL,textLength);
}
WFAligner::AlignmentStatus WFAligner::alignEndsFree(
    const char* const pattern,
    const int patternLength,
    const int patternBeginFree,
    const int patternEndFree,
    const char* const text,
    const int textLength,
    const int textBeginFree,
    const int textEndFree) {
  // Configure
  wavefront_aligner_set_alignment_free_ends(wfAligner,
      patternBeginFree,patternEndFree,
      textBeginFree,textEndFree);
  // Align
  return (WFAligner::AlignmentStatus) wavefront_align(wfAligner,pattern,patternLength,text,textLength);
}
WFAligner::AlignmentStatus WFAligner::alignEndsFree(
    std::string& pattern,
    const int patternBeginFree,
    const int patternEndFree,
    std::string& text,
    const int textBeginFree,
    const int textEndFree) {
  // Delegate
  return alignEndsFree(
      pattern.c_str(),pattern.length(),
      patternBeginFree,patternEndFree,
      text.c_str(),text.length(),
      textBeginFree,textEndFree);
}
/*
 * Alignment resume
 */
WFAligner::AlignmentStatus WFAligner::alignResume() {
  // Resume alignment
  return (WFAligner::AlignmentStatus) wavefront_align_resume(wfAligner);
}
/*
 * Heuristics
 */
void WFAligner::setHeuristicNone() {
  wavefront_aligner_set_heuristic_none(wfAligner);
}
void WFAligner::setHeuristicBandedStatic(
    const int band_min_k,
    const int band_max_k) {
  wavefront_aligner_set_heuristic_banded_static(
      wfAligner,band_min_k,band_max_k);
}
void WFAligner::setHeuristicBandedAdaptive(
    const int band_min_k,
    const int band_max_k,
    const int steps_between_cutoffs) {
  wavefront_aligner_set_heuristic_banded_adaptive(
      wfAligner,band_min_k,band_max_k,steps_between_cutoffs);
}
void WFAligner::setHeuristicWFadaptive(
    const int min_wavefront_length,
    const int max_distance_threshold,
    const int steps_between_cutoffs) {
  wavefront_aligner_set_heuristic_wfadaptive(
      wfAligner,min_wavefront_length,
      max_distance_threshold,steps_between_cutoffs);
}
void WFAligner::setHeuristicWFmash(
    const int min_wavefront_length,
    const int max_distance_threshold,
    const int steps_between_cutoffs) {
  wavefront_aligner_set_heuristic_wfmash(
      wfAligner,min_wavefront_length,
      max_distance_threshold,steps_between_cutoffs);
}
void WFAligner::setHeuristicXDrop(
    const int xdrop,
    const int steps_between_cutoffs) {
  wavefront_aligner_set_heuristic_xdrop(
      wfAligner,xdrop,steps_between_cutoffs);
}
void WFAligner::setHeuristicZDrop(
    const int zdrop,
    const int steps_between_cutoffs) {
  wavefront_aligner_set_heuristic_zdrop(
      wfAligner,zdrop,steps_between_cutoffs);
}
/*
 * Custom extend-match function (lambda)
 */
void WFAligner::setMatchFunct(
    int (*matchFunct)(int,int,void*),
    void* matchFunctArguments) {
  wavefront_aligner_set_match_funct(wfAligner,matchFunct,matchFunctArguments);
}
/*
 * Limits
 */
void WFAligner::setMaxAlignmentScore(
    const int maxAlignmentScore) {
  wavefront_aligner_set_max_alignment_score(
      wfAligner,maxAlignmentScore);
}
void WFAligner::setMaxMemory(
    const uint64_t maxMemoryResident,
    const uint64_t maxMemoryAbort) {
  wavefront_aligner_set_max_memory(wfAligner,maxMemoryResident,maxMemoryAbort);
}
// Parallelization
void WFAligner::setMaxNumThreads(
        const int maxNumThreads) {
    wavefront_aligner_set_max_num_threads(wfAligner, maxNumThreads);
}
void WFAligner::setMinOffsetsPerThread(
        const int minOffsetsPerThread) {
    wavefront_aligner_set_min_offsets_per_thread(wfAligner, minOffsetsPerThread);
}
/*
 * Accessors
 */
int WFAligner::getAlignmentScore() {
  return wfAligner->cigar->score;
}
int WFAligner::getAlignmentStatus() {
  return wfAligner->align_status.status;
}
void WFAligner::getAlignmentCigar(
    char** const cigarOperations,
    int* cigarLength) {
 *cigarOperations = wfAligner->cigar->operations + wfAligner->cigar->begin_offset;
 *cigarLength = wfAligner->cigar->end_offset - wfAligner->cigar->begin_offset;
}
std::string WFAligner::getAlignmentCigar() {
  // Fetch CIGAR
  char* buffer;
  int bufferLength;
  getAlignmentCigar(&buffer,&bufferLength);
  // Create string and return
  return std::string(buffer,bufferLength);
}
/*
 * Misc
 */
char* WFAligner::strError(
    const int wfErrorCode) {
  return wavefront_align_strerror(wfErrorCode);
}
void WFAligner::setVerbose(
    const int verbose) {
  wfAligner->system.verbose = verbose;
}
/*
 * Indel Aligner (a.k.a Longest Common Subsequence - LCS)
 */
WFAlignerIndel::WFAlignerIndel(
    const AlignmentScope alignmentScope,
    const MemoryModel memoryModel) :
        WFAligner(alignmentScope,memoryModel) {
  attributes.distance_metric = indel;
  wfAligner = wavefront_aligner_new(&attributes);
}
/*
 * Edit Aligner (a.k.a Levenshtein)
 */
WFAlignerEdit::WFAlignerEdit(
    const AlignmentScope alignmentScope,
    const MemoryModel memoryModel) :
        WFAligner(alignmentScope,memoryModel) {
  attributes.distance_metric = edit;
  wfAligner = wavefront_aligner_new(&attributes);
}
/*
 * Gap-Linear Aligner (a.k.a Needleman-Wunsch)
 */
WFAlignerGapLinear::WFAlignerGapLinear(
    const int mismatch,
    const int indel,
    const AlignmentScope alignmentScope,
    const MemoryModel memoryModel) :
        WFAligner(alignmentScope,memoryModel) {
  attributes.distance_metric = gap_linear;
  attributes.linear_penalties.match = 0;
  attributes.linear_penalties.mismatch = mismatch;
  attributes.linear_penalties.indel = indel;
  wfAligner = wavefront_aligner_new(&attributes);
}
WFAlignerGapLinear::WFAlignerGapLinear(
    const int match,
    const int mismatch,
    const int indel,
    const AlignmentScope alignmentScope,
    const MemoryModel memoryModel) :
        WFAligner(alignmentScope,memoryModel) {
  attributes.distance_metric = gap_linear;
  attributes.linear_penalties.match = match;
  attributes.linear_penalties.mismatch = mismatch;
  attributes.linear_penalties.indel = indel;
  wfAligner = wavefront_aligner_new(&attributes);
}
/*
 * Gap-Affine Aligner (a.k.a Smith-Waterman-Gotoh)
 */
WFAlignerGapAffine::WFAlignerGapAffine(
    const int mismatch,
    const int gapOpening,
    const int gapExtension,
    const AlignmentScope alignmentScope,
    const MemoryModel memoryModel) :
        WFAligner(alignmentScope,memoryModel) {
  attributes.distance_metric = gap_affine;
  attributes.affine_penalties.match = 0;
  attributes.affine_penalties.mismatch = mismatch;
  attributes.affine_penalties.gap_opening = gapOpening;
  attributes.affine_penalties.gap_extension = gapExtension;
  wfAligner = wavefront_aligner_new(&attributes);
}
WFAlignerGapAffine::WFAlignerGapAffine(
    const int match,
    const int mismatch,
    const int gapOpening,
    const int gapExtension,
    const AlignmentScope alignmentScope,
    const MemoryModel memoryModel) :
        WFAligner(alignmentScope,memoryModel) {
  attributes.distance_metric = gap_affine;
  attributes.affine_penalties.match = match;
  attributes.affine_penalties.mismatch = mismatch;
  attributes.affine_penalties.gap_opening = gapOpening;
  attributes.affine_penalties.gap_extension = gapExtension;
  wfAligner = wavefront_aligner_new(&attributes);
}
/*
 * Gap-Affine Dual-Cost Aligner (a.k.a. concave 2-pieces)
 */
WFAlignerGapAffine2Pieces::WFAlignerGapAffine2Pieces(
    const int mismatch,
    const int gapOpening1,
    const int gapExtension1,
    const int gapOpening2,
    const int gapExtension2,
    const AlignmentScope alignmentScope,
    const MemoryModel memoryModel) :
        WFAligner(alignmentScope,memoryModel) {
  attributes.distance_metric = gap_affine_2p;
  attributes.affine2p_penalties.match = 0;
  attributes.affine2p_penalties.mismatch = mismatch;
  attributes.affine2p_penalties.gap_opening1 = gapOpening1;
  attributes.affine2p_penalties.gap_extension1 = gapExtension1;
  attributes.affine2p_penalties.gap_opening2 = gapOpening2;
  attributes.affine2p_penalties.gap_extension2 = gapExtension2;
  wfAligner = wavefront_aligner_new(&attributes);
}
WFAlignerGapAffine2Pieces::WFAlignerGapAffine2Pieces(
    const int match,
    const int mismatch,
    const int gapOpening1,
    const int gapExtension1,
    const int gapOpening2,
    const int gapExtension2,
    const AlignmentScope alignmentScope,
    const MemoryModel memoryModel) :
        WFAligner(alignmentScope,memoryModel) {
  attributes.distance_metric = gap_affine_2p;
  attributes.affine2p_penalties.match = match;
  attributes.affine2p_penalties.mismatch = mismatch;
  attributes.affine2p_penalties.gap_opening1 = gapOpening1;
  attributes.affine2p_penalties.gap_extension1 = gapExtension1;
  attributes.affine2p_penalties.gap_opening2 = gapOpening2;
  attributes.affine2p_penalties.gap_extension2 = gapExtension2;
  wfAligner = wavefront_aligner_new(&attributes);
}

} /* namespace wfa */