Bioplib
Protein Structure C Library
|
Perform Needleman & Wunsch sequence alignment. More...
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <ctype.h>
#include "SysDefs.h"
#include "macros.h"
#include "array.h"
#include "general.h"
#include "seq.h"
Go to the source code of this file.
Data Structures | |
struct | XY |
Macros | |
#define | MAX3(c, d, e) (MAX(MAX((c),(d)),(e))) |
#define | DATAENV "DATADIR" /* Environment variable or assign */ |
#define | MAXBUFF 400 |
#define | MAXWORD 16 |
Functions | |
int | blAlign (char *seq1, int length1, char *seq2, int length2, BOOL verbose, BOOL identity, int penalty, char *align1, char *align2, int *align_len) |
int | blAffinealign (char *seq1, int length1, char *seq2, int length2, BOOL verbose, BOOL identity, int penalty, int penext, char *align1, char *align2, int *align_len) |
int | blAffinealignuc (char *seq1, int length1, char *seq2, int length2, BOOL verbose, BOOL identity, int penalty, int penext, char *align1, char *align2, int *align_len) |
BOOL | blReadMDM (char *mdmfile) |
int | blCalcMDMScore (char resa, char resb) |
int | blCalcMDMScoreUC (char resa, char resb) |
int | blZeroMDM (void) |
void | blSetMDMScoreWeight (char resa, char resb, REAL weight) |
Perform Needleman & Wunsch sequence alignment.
This code is NOT IN THE PUBLIC DOMAIN, but it may be copied according to the conditions laid out in the accompanying file COPYING.DOC.
The code may be modified as required, but any modifications must be documented so that the person responsible can be identified.
The code may not be sold commercially or included as part of a commercial product except as described in the file COPYING.DOC.
A simple Needleman & Wunsch Dynamic Programming alignment of 2 sequences.
A window is not used so the routine may be a bit slow on long sequences.
First call ReadMDM() to read the mutation data matrix, then call align() to align the sequences.
Definition in file align.c.
#define DATAENV "DATADIR" /* Environment variable or assign */ |
int blAffinealign | ( | char * | seq1, |
int | length1, | ||
char * | seq2, | ||
int | length2, | ||
BOOL | verbose, | ||
BOOL | identity, | ||
int | penalty, | ||
int | penext, | ||
char * | align1, | ||
char * | align2, | ||
int * | align_len | ||
) |
[in] | *seq1 | First sequence |
[in] | length1 | First sequence length |
[in] | *seq2 | Second sequence |
[in] | length2 | Second sequence length |
[in] | verbose | Display N&W matrix |
[in] | identity | Use identity matrix |
[in] | penalty | Gap insertion penalty value |
[in] | penext | Extension penalty |
[out] | *align1 | Sequence 1 aligned |
[out] | *align2 | Sequence 2 aligned |
[out] | *align_len | Alignment length |
Perform simple N&W alignment of seq1 and seq2. No window is used, so will be slow for long sequences.
Note that you must allocate sufficient memory for the aligned sequences. The easy way to do this is to ensure that align1 and align2 are of length (length1+length2).
NOTE AND CHANGES SHOULD BE PROPAGATED TO affinealignuc() ******
int blAffinealignuc | ( | char * | seq1, |
int | length1, | ||
char * | seq2, | ||
int | length2, | ||
BOOL | verbose, | ||
BOOL | identity, | ||
int | penalty, | ||
int | penext, | ||
char * | align1, | ||
char * | align2, | ||
int * | align_len | ||
) |
[in] | *seq1 | First sequence |
[in] | length1 | First sequence length |
[in] | *seq2 | Second sequence |
[in] | length2 | Second sequence length |
[in] | verbose | Display N&W matrix |
[in] | identity | Use identity matrix |
[in] | penalty | Gap insertion penalty value |
[in] | penext | Extension penalty |
[out] | *align1 | Sequence 1 aligned |
[out] | *align2 | Sequence 2 aligned |
[out] | *align_len | Alignment length |
Perform simple N&W alignment of seq1 and seq2. No window is used, so will be slow for long sequences.
Note that you must allocate sufficient memory for the aligned sequences. The easy way to do this is to ensure that align1 and align2 are of length (length1+length2).
NOTE AND CHANGES SHOULD BE PROPAGATED TO affinealign() ******
int blAlign | ( | char * | seq1, |
int | length1, | ||
char * | seq2, | ||
int | length2, | ||
BOOL | verbose, | ||
BOOL | identity, | ||
int | penalty, | ||
char * | align1, | ||
char * | align2, | ||
int * | align_len | ||
) |
[in] | *seq1 | First sequence |
[in] | length1 | First sequence length |
[in] | *seq2 | Second sequence |
[in] | length2 | Second sequence length |
[in] | verbose | Display N&W matrix |
[in] | identity | Use identity matrix |
[in] | penalty | Gap insertion penalty value |
[out] | *align1 | Sequence 1 aligned |
[out] | *align2 | Sequence 2 aligned |
[out] | *align_len | Alignment length |
Perform simple N&W alignment of seq1 and seq2. No window is used, so will be slow for long sequences.
A single gap penalty is used, so gap extension incurrs no further penalty.
Note that you must allocate sufficient memory for the aligned sequences. The easy way to do this is to ensure that align1 and align2 are of length (length1+length2).
int blCalcMDMScore | ( | char | resa, |
char | resb | ||
) |
[in] | resa | First residue |
[in] | resb | Second residue |
Calculate score from static globally stored mutation data matrix
If both residues are set as '\0' it will simply silence all warnings
int blCalcMDMScoreUC | ( | char | resa, |
char | resb | ||
) |
[in] | resa | First residue |
[in] | resb | Second residue |
Calculate score from static globally stored mutation data matrix
BOOL blReadMDM | ( | char * | mdmfile | ) |
[in] | *mdmfile | Mutation data matrix filename |
Read mutation data matrix into static global arrays. The matrix may have comments at the start introduced with a ! in the first column. The matrix must be complete (i.e. a triangular matrix will not work). A line describing the residue types must appear, and may be placed before or after the matrix itself
void blSetMDMScoreWeight | ( | char | resa, |
char | resb, | ||
REAL | weight | ||
) |