Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
O
openstructure
Manage
Activity
Members
Code
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Deploy
Releases
Container Registry
Model registry
Analyze
Contributor analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
schwede
openstructure
Commits
fe622ee2
Commit
fe622ee2
authored
6 months ago
by
Studer Gabriel
Browse files
Options
Downloads
Patches
Plain Diff
bugfix: avoid segfault when aligning empty sequences with parasail
SCHWED-6400
parent
7899f5fb
No related branches found
Branches containing commit
No related tags found
Tags containing commit
No related merge requests found
Changes
1
Hide whitespace changes
Inline
Side-by-side
Showing
1 changed file
modules/seq/alg/src/wrap_parasail.cc
+49
-0
49 additions, 0 deletions
modules/seq/alg/src/wrap_parasail.cc
with
49 additions
and
0 deletions
modules/seq/alg/src/wrap_parasail.cc
+
49
−
0
View file @
fe622ee2
...
...
@@ -17,6 +17,7 @@
// 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
//------------------------------------------------------------------------------
#include
<ost/message.hh>
#include
<ost/seq/alg/wrap_parasail.hh>
#if OST_PARASAIL_ENABLED
...
...
@@ -131,6 +132,45 @@ ost::seq::AlignmentList Align(const ost::seq::ConstSequenceHandle& s1,
return
aln_list
;
}
ost
::
seq
::
AlignmentList
AlignEmptySeq
(
const
ost
::
seq
::
ConstSequenceHandle
&
s1
,
const
ost
::
seq
::
ConstSequenceHandle
&
s2
)
{
ost
::
seq
::
AlignmentHandle
aln
=
ost
::
seq
::
CreateAlignment
();
if
(
s1
.
GetLength
()
==
0
&&
s2
.
GetLength
()
==
0
)
{
ost
::
seq
::
SequenceHandle
new_s1
=
ost
::
seq
::
CreateSequence
(
s1
.
GetName
(),
s1
.
GetString
());
new_s1
.
SetOffset
(
s1
.
GetOffset
());
aln
.
AddSequence
(
new_s1
);
ost
::
seq
::
SequenceHandle
new_s2
=
ost
::
seq
::
CreateSequence
(
s2
.
GetName
(),
s2
.
GetString
());
new_s2
.
SetOffset
(
s2
.
GetOffset
());
aln
.
AddSequence
(
new_s2
);
}
else
if
(
s1
.
GetLength
()
==
0
)
{
String
x
(
s2
.
GetLength
(),
'-'
);
ost
::
seq
::
SequenceHandle
new_s1
=
ost
::
seq
::
CreateSequence
(
s1
.
GetName
(),
x
);
new_s1
.
SetOffset
(
s1
.
GetOffset
());
aln
.
AddSequence
(
new_s1
);
ost
::
seq
::
SequenceHandle
new_s2
=
ost
::
seq
::
CreateSequence
(
s2
.
GetName
(),
s2
.
GetString
());
new_s2
.
SetOffset
(
s2
.
GetOffset
());
aln
.
AddSequence
(
new_s2
);
}
else
if
(
s2
.
GetLength
()
==
0
)
{
ost
::
seq
::
SequenceHandle
new_s1
=
ost
::
seq
::
CreateSequence
(
s1
.
GetName
(),
s1
.
GetString
());
new_s1
.
SetOffset
(
s1
.
GetOffset
());
aln
.
AddSequence
(
new_s1
);
String
x
(
s1
.
GetLength
(),
'-'
);
ost
::
seq
::
SequenceHandle
new_s2
=
ost
::
seq
::
CreateSequence
(
s2
.
GetName
(),
x
);
new_s2
.
SetOffset
(
s2
.
GetOffset
());
aln
.
AddSequence
(
new_s2
);
}
else
{
throw
ost
::
Error
(
"One seq must be empty in AlignEmptySeq"
);
}
if
(
s1
.
HasAttachedView
())
{
aln
.
AttachView
(
0
,
s1
.
GetAttachedView
());
}
if
(
s2
.
HasAttachedView
())
{
aln
.
AttachView
(
1
,
s2
.
GetAttachedView
());
}
ost
::
seq
::
AlignmentList
aln_list
=
{
aln
};
return
aln_list
;
}
}
namespace
ost
{
namespace
seq
{
namespace
alg
{
...
...
@@ -139,6 +179,9 @@ ost::seq::AlignmentList ParaLocalAlign(const ost::seq::ConstSequenceHandle& s1,
const
ost
::
seq
::
ConstSequenceHandle
&
s2
,
ost
::
seq
::
alg
::
SubstWeightMatrixPtr
&
subst
,
int
gap_open
,
int
gap_ext
)
{
if
(
s1
.
GetLength
()
==
0
||
s2
.
GetLength
()
==
0
)
{
return
ost
::
seq
::
AlignmentList
();
}
return
Align
(
s1
,
s2
,
gap_open
,
gap_ext
,
subst
,
parasail_sw_trace_scan_sat
,
true
);
}
...
...
@@ -147,6 +190,9 @@ ost::seq::AlignmentList ParaGlobalAlign(const ost::seq::ConstSequenceHandle& s1,
const
ost
::
seq
::
ConstSequenceHandle
&
s2
,
ost
::
seq
::
alg
::
SubstWeightMatrixPtr
&
subst
,
int
gap_open
,
int
gap_ext
)
{
if
(
s1
.
GetLength
()
==
0
||
s2
.
GetLength
()
==
0
)
{
return
AlignEmptySeq
(
s1
,
s2
);
}
return
Align
(
s1
,
s2
,
gap_open
,
gap_ext
,
subst
,
parasail_nw_trace_scan_sat
,
false
);
}
...
...
@@ -155,6 +201,9 @@ ost::seq::AlignmentList ParaSemiGlobalAlign(const ost::seq::ConstSequenceHandle&
const
ost
::
seq
::
ConstSequenceHandle
&
s2
,
ost
::
seq
::
alg
::
SubstWeightMatrixPtr
&
subst
,
int
gap_open
,
int
gap_ext
)
{
if
(
s1
.
GetLength
()
==
0
||
s2
.
GetLength
()
==
0
)
{
return
AlignEmptySeq
(
s1
,
s2
);
}
return
Align
(
s1
,
s2
,
gap_open
,
gap_ext
,
subst
,
parasail_sg_trace_scan_sat
,
false
);
}
...
...
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment