Class containing an alignment
This class stores the alignment sequences with it's names, residues and extra information.
It contains multiple methods regarding the sequences.
It also contains submodules that provide methods for Calculating statistics, Trim the alignment and Printing sequences information.
More...
#include <Alignment.h>
Classes | |
class | sequencesMatrix |
Internal Alignment Class that represents a sequences matrix. More... | |
Public Member Functions | |
bool | fillMatrices (bool aligned, bool checkInvalidChars=true) |
Method to initialize matrices and check if the sequences have been correctly loaded and are free of errors It checks if sequences contain unknown characters, sets the isAligned flag and initializes 'saveResidues' and 'saveSequences', depending on the sizes of the sequences (whether or not they have the same length). More... | |
Alignment () | |
Constructor. More... | |
Alignment (Alignment &) | |
Copy Constructor. The copy constructor does not make a copy of all structures present on the original alignment. Instead, all alignments that share a common origin, share some information, and also, a counter (Alignment::SeqRef) of how many alignments share the same information. This allows us to keep the information on memory while any alignment is using it. When an alignment is destroyed, the counter is substracted by 1, and, if it is 0 after the substraction, the shared data is destroyed. More... | |
~Alignment () | |
Destructor. The destructor takes care of removing the non-shared information of each alignment when is destroyed. It also substracts 1 from the shared counter Alignment::SeqRef. When Alignment::SeqRef arrives to 0, the shared data is destroyed, including Alignment::SeqRef. More... | |
int | getNumSpecies () |
Number of sequences getter. More... | |
void | getSequences (std::string *names) |
Getter for the sequences names. More... | |
void | getSequences (std::string *names, int *lenghts) |
Getter for the sequences names and its lenghts. More... | |
void | getSequences (std::string *names, std::string *sequences, int *lenghts) |
Getter for the sequences, its names and lenghts. More... | |
bool | getSequenceNameOrder (std::string *names, int *orderVector) |
Method to map two sets of sequences, own and external. The method accepts a vector of names, and test them against the alignment own sequences' names. If the sequence is present, its index will be stored in orderVector. . More... | |
int | getNumAminos () |
Residues number getter. It counts gaps as residue. More... | |
void | setWindowsSize (int ghWindow, int shWindow) |
Windows setter. More... | |
void | setBlockSize (int blockSize) |
BlockSize Setter. More... | |
void | calculateSeqOverlap () |
Calculates overlap between sequences. More... | |
void | printSeqIdentity () |
Method to print different identity values computed from the alignment. In this method we assess the identity values matrix, as well as different average values. Moreover, the method computes which one is the most similar sequence. More... | |
void | printSeqOverlap () |
Prints the overlap between sequences. More... | |
int | getAlignmentType () const |
Alignment type getter. See SequenceTypes. More... | |
bool | isFileAligned () |
isAligned getter. More... | |
Alignment * | getTranslationCDS (Alignment *proteinAlignment) |
Method to back translate a protein alignment using the sequences present on the current alignment. Current alignment should contain the original sequences which were translated to protein sequences, without any gap. More... | |
bool | checkCorrespondence (std::string *names, int *lenghts, int totalInputSequences, int multiple) |
Function to check CDS file. . More... | |
void | calculateColIdentity (float *columnIdentity) |
Method that calculates the columns identity value. This is, the frequency of the most present element in the column, being residue, indetermination and gap allowed. More... | |
void | setKeepSequencesFlag (bool newFlagValue) |
Keep Sequences setter. More... | |
void | printAlignmentInfo (std::ostream &output) |
Print information about sequences number, average sequence length, maximum and minimum sequences length, etc. More... | |
bool | prepareCodingSequence (bool splitByStopCodon, bool ignoreStopCodon, Alignment *proteinAlignment) |
Method to check if the CDS file is correct. Based on nature of residues: DNA/RNA (Most of the residues) There is no gaps on the whole dataset. Each sequence is a multiple of 3. It will also remove or split stop codons depending on the flags passed. More... | |
bool | alignmentSummaryHTML (const Alignment &trimmedAlig, const char *const destFile) |
Method to report the trimming results in HTML. Outputs an HTML file that shows visually what has been done to the alignment (removed sequences or residues) and the stats used to do so. More... | |
bool | statSVG (const char *const destFile) |
bool | alignmentSummarySVG (Alignment &trimmedAlig, const char *destFile, int blocks) |
Method to report the trimming results in SVG. Outputs a SVG file that shows visually what has been done to the alignment (removed sequences or residues) and the stats used to do so. More... | |
void | updateSequencesAndResiduesNums (bool countSequences=true, bool countResidues=true) |
Updates the sequence number and residue number based on saveResidues and saveSequences. More... | |
Public Attributes | |
Cleaner * | Cleaning = nullptr |
Trimming submodule. It contains methods and variables related to trimming. More... | |
statistics::Manager * | Statistics = nullptr |
Statistics submodule. It contains methods and variables related to statistics calculation and reporting. More... | |
sequencesMatrix * | SequencesMatrix = nullptr |
SequencesMatrix submodule This helper module stores the sequences without gaps. Useful when comparing MSA. More... | |
int * | SeqRef = nullptr |
Pointer to keep a count of how many alignments depend on the same shared-sequences. Increased on each trim step, and automatically reduced on each alignment destroy. If it arrives to 0, the last alignment being destroyed will also take care of releasing memory for shared information. More... | |
int | originalNumberOfSequences = 0 |
Number of sequences the alignment had when it was loaded. More... | |
int | numberOfSequences = 0 |
Number of sequences present on the alignment. More... | |
int | originalNumberOfResidues = 0 |
Number of residues the alignment had when it was loaded. More... | |
int | numberOfResidues = 0 |
Number of residues present on the alignment if it is aligned. More... | |
bool | isAligned = false |
Flag that indicates if all sequences on the alignment have the same length (Including gaps). More... | |
std::string * | sequences = nullptr |
Vector containing the sequences. More... | |
std::string * | seqsName = nullptr |
String vector containing the sequences names. More... | |
std::string * | seqsInfo = nullptr |
Vector containing extra information about each sequence. More... | |
std::string | filename |
Filename where this alignment was loaded from. More... | |
std::string | alignmentInfo |
String containing information of the alignment as a whole. More... | |
float ** | identities = nullptr |
2D Matrix that represents how much each sequence resembles other sequences in the same MSA More... | |
float ** | overlaps = nullptr |
Matrix that stores the overlap between sequences of an alignment. See Alignment::calculateSeqOverlap. More... | |
int * | saveResidues = nullptr |
Vector containing which residues will be kept/removed. -1 Indicates a 'removed/toRemove'. Otherwise residue/column will be kept. More... | |
int * | saveSequences = nullptr |
Vector containing which sequences will be kept/removed. -1 Indicates a 'removed/toRemove'. Otherwise sequence will be kept. More... | |
Private Attributes | |
int | dataType = SequenceTypes::NotDefined |
Flag variable that represents the type of the alignment. More... | |
Class containing an alignment
This class stores the alignment sequences with it's names, residues and extra information.
It contains multiple methods regarding the sequences.
It also contains submodules that provide methods for Calculating statistics, Trim the alignment and Printing sequences information.
Definition at line 49 of file Alignment.h.
Alignment::Alignment | ( | void | ) |
Constructor.
Definition at line 46 of file Alignment.cpp.
References Cleaner::Cleaner(), Cleaning, dataType, identities, isAligned, statistics::Manager::Manager(), NotDefined, numberOfResidues, numberOfSequences, originalNumberOfResidues, originalNumberOfSequences, saveResidues, saveSequences, SeqRef, seqsInfo, seqsName, sequences, SequencesMatrix, and Statistics.
Referenced by getTranslationCDS(), and FormatHandling::FormatManager::splitAlignmentKeeping().
Alignment::Alignment | ( | Alignment & | originalAlignment | ) |
Copy Constructor.
The copy constructor does not make a copy of all structures present on the original alignment.
Instead, all alignments that share a common origin, share some information, and also, a counter (Alignment::SeqRef) of how many alignments share the same information.
This allows us to keep the information on memory while any alignment is using it.
When an alignment is destroyed, the counter is substracted by 1, and, if it is 0 after the substraction, the shared data is destroyed.
Definition at line 87 of file Alignment.cpp.
References alignmentInfo, Cleaner::Cleaner(), Cleaning, dataType, filename, identities, isAligned, statistics::Manager::Manager(), numberOfResidues, numberOfSequences, originalNumberOfResidues, originalNumberOfSequences, saveResidues, saveSequences, SeqRef, seqsInfo, seqsName, sequences, SequencesMatrix, and Statistics.
Referenced by Cleaner::cleanByCutValueFallBehind(), Cleaner::cleanByCutValueOverpass(), Cleaner::cleanByCutValueOverpassOrEquals(), Cleaner::cleanOverlapSeq(), Cleaner::cleanStrict(), Cleaner::getClustering(), statistics::Consistency::perform(), Cleaner::removeColumns(), and Cleaner::removeSequences().
Alignment::~Alignment | ( | void | ) |
Destructor. The destructor takes care of removing the non-shared information of each alignment when is destroyed.
It also substracts 1 from the shared counter Alignment::SeqRef.
When Alignment::SeqRef arrives to 0, the shared data is destroyed, including Alignment::SeqRef.
Definition at line 134 of file Alignment.cpp.
References Cleaning, identities, numberOfSequences, saveResidues, saveSequences, SeqRef, seqsInfo, seqsName, sequences, SequencesMatrix, and Statistics.
Alignment::alignmentSummaryHTML | ( | const Alignment & | trimmedAlig, |
const char *const | destFile | ||
) |
Method to report the trimming results in HTML.
Outputs an HTML file that shows visually what has been done to the alignment (removed sequences or residues) and the stats used to do so.
trimmedAlig | Trimmed alignment to compare with. |
destFile | Filename where to save the results. |
Definition at line 1167 of file Alignment.cpp.
References statistics::Manager::consistency, debug, utils::determineColor(), filename, statistics::Manager::gaps, statistics::Gaps::getGapsWindow(), statistics::Similarity::getMdkWindowedVector(), statistics::Consistency::getValues(), isAligned, utils::max(), NotAligned, numberOfResidues, numberOfSequences, originalNumberOfResidues, originalNumberOfSequences, reporting::reportManager::report(), saveResidues, saveSequences, seqsName, sequences, statistics::Manager::similarity, and Statistics.
Referenced by trimAlManager::output_reports().
Alignment::alignmentSummarySVG | ( | Alignment & | trimmedAlig, |
const char * | destFile, | ||
int | blocks | ||
) |
Method to report the trimming results in SVG.
Outputs a SVG file that shows visually what has been done to the alignment (removed sequences or residues) and the stats used to do so.
trimmedAlig | Trimmed alignment obtained from the object calling this function. |
destFile | Filename where to save the results. |
blocks | Size in residues to report. 120 by default. |
Definition at line 1708 of file Alignment.cpp.
References statistics::Manager::consistency, debug, utils::determineColor(), filename, statistics::Manager::gaps, statistics::Gaps::getGapsWindow(), statistics::Similarity::getMdkWindowedVector(), statistics::Consistency::getValues(), isAligned, utils::max(), NotAligned, numberOfResidues, numberOfSequences, originalNumberOfResidues, originalNumberOfSequences, reporting::reportManager::report(), saveResidues, saveSequences, seqsName, sequences, statistics::Manager::similarity, and Statistics.
Referenced by trimAlManager::output_reports().
Alignment::calculateColIdentity | ( | float * | columnIdentity | ) |
Method that calculates the columns identity value.
This is, the frequency of the most present element in the column, being residue, indetermination and gap allowed.
[out] | columnIdentity | Vector to fill with identities for each column. |
Definition at line 1074 of file Alignment.cpp.
References AA, getAlignmentType(), numberOfResidues, numberOfSequences, and sequences.
Alignment::calculateSeqOverlap | ( | ) |
Calculates overlap between sequences.
Definition at line 378 of file Alignment.cpp.
References AA, getAlignmentType(), numberOfResidues, numberOfSequences, overlaps, and sequences.
Referenced by printSeqOverlap().
Alignment::checkCorrespondence | ( | std::string * | names, |
int * | lenghts, | ||
int | totalInputSequences, | ||
int | multiple = 1 |
||
) |
Function to check CDS file.
.
names | Vector containing the names to check. |
lenghts | Vector containing the length of each sequence. |
totalInputSequences | Number of sequences present. |
multiple | Multiplies the length of each sequence by this number. |
Definition at line 716 of file Alignment.cpp.
References debug, IncludingIndeterminationSymbols, LessNucleotidesThanExpected, utils::min(), numberOfSequences, utils::removeCharacter(), reporting::reportManager::report(), seqsName, SequenceNotPresentInCDS, sequences, and SequenceWillBeCut.
Referenced by trimAlManager::check_backtranslation_infile_names_correspondence().
Alignment::fillMatrices | ( | bool | aligned, |
bool | checkInvalidChars = true |
||
) |
Method to initialize matrices and check if the sequences have been correctly loaded and are free of errors
It checks if sequences contain unknown characters, sets the isAligned flag and initializes 'saveResidues' and 'saveSequences', depending on the sizes of the sequences (whether or not they have the same length).
aligned | Flag to make the method check if the alignment is aligned or not. |
Definition at line 794 of file Alignment.cpp.
References debug, filename, isAligned, NotAligned, numberOfResidues, numberOfSequences, reporting::reportManager::report(), saveResidues, saveSequences, seqsName, sequences, SequencesNotSameSize, and UnknownCharacter.
Alignment::getAlignmentType | ( | void | ) | const |
Alignment type getter.
See SequenceTypes.
Definition at line 469 of file Alignment.cpp.
References utils::checkAlignmentType(), dataType, NotDefined, numberOfResidues, numberOfSequences, and sequences.
Referenced by calculateColIdentity(), statistics::Similarity::calculateMatrixIdentity(), Cleaner::calculateSeqIdentity(), calculateSeqOverlap(), Cleaner::calculateSpuriousVector(), statistics::Similarity::calculateVectors(), trimAlManager::check_coding_sequences_type(), trimAlManager::create_or_use_similarity_matrix(), main(), statistics::Consistency::perform(), and prepareCodingSequence().
Alignment::getNumAminos | ( | void | ) |
Residues number getter. It counts gaps as residue.
Definition at line 252 of file Alignment.cpp.
References numberOfResidues.
Referenced by trimAlManager::check_block_size(), trimAlManager::CleanResiduesNonAuto(), statistics::Consistency::compareAndChoose(), statistics::Consistency::forceComparison(), and statistics::Consistency::perform().
Alignment::getNumSpecies | ( | void | ) |
Number of sequences getter.
Definition at line 245 of file Alignment.cpp.
References numberOfSequences.
Referenced by trimAlManager::check_backtranslation_infile_names_correspondence(), trimAlManager::check_clusters_incompatibilities(), statistics::Consistency::compareAndChoose(), statistics::Consistency::forceComparison(), and prepareCodingSequence().
Alignment::getSequenceNameOrder | ( | std::string * | names, |
int * | orderVector | ||
) |
Method to map two sets of sequences, own and external.
The method accepts a vector of names, and test them against the alignment own sequences' names.
If the sequence is present, its index will be stored in orderVector.
.
names | Vector of sequences names to map. | |
[out] | orderVector | Vector of orders to fill. |
Definition at line 445 of file Alignment.cpp.
References numberOfSequences, and seqsName.
Referenced by statistics::Consistency::compareAndChoose(), and statistics::Consistency::forceComparison().
Alignment::getSequences | ( | std::string * | names | ) |
Getter for the sequences names.
[out] | names | Vector of sequences names to fill. |
Definition at line 416 of file Alignment.cpp.
References numberOfSequences, and seqsName.
Referenced by statistics::Consistency::compareAndChoose(), and statistics::Consistency::forceComparison().
Alignment::getSequences | ( | std::string * | names, |
int * | lenghts | ||
) |
Getter for the sequences names and its lenghts.
[out] | names | Vector of sequences names to fill. |
[out] | lenghts | Vector of lenghts to fill. |
Definition at line 424 of file Alignment.cpp.
References numberOfSequences, utils::removeCharacter(), seqsName, and sequences.
Referenced by trimAlManager::check_backtranslation_infile_names_correspondence().
Alignment::getSequences | ( | std::string * | names, |
std::string * | sequences, | ||
int * | lenghts | ||
) |
Getter for the sequences, its names and lenghts.
[out] | names | Vector of sequences names to fill. |
[out] | sequences | Vector of sequences to fill. |
[out] | lenghts | Vector of lenghts to fill. |
Definition at line 434 of file Alignment.cpp.
References numberOfSequences, utils::removeCharacter(), seqsName, and sequences.
Referenced by prepareCodingSequence().
Alignment::getTranslationCDS | ( | Alignment * | proteinAlignment | ) |
Method to back translate a protein alignment using the sequences present on the current alignment.
Current alignment should contain the original sequences which were translated to protein sequences, without any gap.
proteinAlignment | Protein Alignment to use as origin. |
Definition at line 169 of file Alignment.cpp.
References Alignment(), numberOfResidues, numberOfSequences, originalNumberOfResidues, originalNumberOfSequences, saveResidues, saveSequences, seqsInfo, seqsName, and sequences.
Referenced by trimAlManager::postprocess_alignment().
Alignment::isFileAligned | ( | void | ) |
isAligned getter.
Definition at line 495 of file Alignment.cpp.
References isAligned.
Referenced by trimAlManager::check_backtranslations(), trimAlManager::check_file_aligned(), trimAlManager::clean_alignment(), statistics::Consistency::perform(), and printAlignmentInfo().
Alignment::prepareCodingSequence | ( | bool | splitByStopCodon, |
bool | ignoreStopCodon, | ||
Alignment * | proteinAlignment | ||
) |
Method to check if the CDS file is correct.
Based on nature of residues: DNA/RNA (Most of the residues)
There is no gaps on the whole dataset.
Each sequence is a multiple of 3.
It will also remove or split stop codons depending on the flags passed.
splitByStopCodon | Flag that informs the method to split sequences if it finds any stop codon. |
ignoreStopCodon | Flag that informs the method to stop reading sequence if it finds any stop codon. |
proteinAlignment | Alignment containing protein sequences which names contains all names in the current alignment. |
Definition at line 528 of file Alignment.cpp.
References AA, CDScontainsProteinSequences, CuttingSequence, debug, getAlignmentType(), getNumSpecies(), getSequences(), numberOfSequences, reporting::reportManager::report(), seqsName, SequenceContainsGap, SequenceHasStopCodon, SequenceNotMultipleOfThree, and sequences.
Referenced by trimAlManager::check_and_prepare_coding_sequence().
Alignment::printAlignmentInfo | ( | std::ostream & | output | ) |
Print information about sequences number, average sequence length, maximum and minimum sequences length, etc.
output | Output stream. |
Definition at line 868 of file Alignment.cpp.
References isFileAligned(), numberOfResidues, numberOfSequences, seqsName, and sequences.
Referenced by main().
Alignment::printSeqIdentity | ( | void | ) |
Method to print different identity values computed from the alignment.
In this method we assess the identity values matrix, as well as different average values.
Moreover, the method computes which one is the most similar sequence.
Definition at line 912 of file Alignment.cpp.
References Cleaner::calculateSeqIdentity(), Cleaning, identities, utils::max(), numberOfSequences, originalNumberOfSequences, and seqsName.
Referenced by trimAlManager::print_statistics().
Alignment::printSeqOverlap | ( | void | ) |
Prints the overlap between sequences.
Definition at line 1001 of file Alignment.cpp.
References calculateSeqOverlap(), utils::max(), numberOfSequences, overlaps, and seqsName.
Referenced by trimAlManager::print_statistics().
Alignment::setBlockSize | ( | int | blockSize | ) |
BlockSize Setter.
blockSize | New value. |
Definition at line 267 of file Alignment.cpp.
References Cleaner::blockSize, and Cleaning.
Referenced by trimAlManager::innerPerform().
Alignment::setKeepSequencesFlag | ( | bool | newFlagValue | ) |
Keep Sequences setter.
newFlagValue | New flag value |
Definition at line 274 of file Alignment.cpp.
References Cleaning, and Cleaner::keepSequences.
Referenced by trimAlManager::innerPerform().
Alignment::setWindowsSize | ( | int | ghWindow, |
int | shWindow | ||
) |
Windows setter.
ghWindow | Half the Gap Window. |
shWindow | Half the Similarity Window. |
Definition at line 259 of file Alignment.cpp.
References statistics::Manager::ghWindow, statistics::Manager::shWindow, and Statistics.
Referenced by trimAlManager::set_window_size().
bool Alignment::statSVG | ( | const char *const | destFile | ) |
Definition at line 1431 of file Alignment.cpp.
References statistics::Manager::consistency, filename, statistics::Manager::gaps, statistics::Gaps::getGapsWindow(), statistics::Similarity::getMdkWindowedVector(), statistics::Consistency::getValues(), numberOfResidues, originalNumberOfResidues, originalNumberOfSequences, utils::quicksort(), statistics::Manager::similarity, and Statistics.
Referenced by trimAlManager::svg_stats_out().
Alignment::updateSequencesAndResiduesNums | ( | bool | countSequences = true , |
bool | countResidues = true |
||
) |
Updates the sequence number and residue number based on saveResidues and saveSequences.
countSequences | Whether to count sequences |
countResidues | Whether to count residues |
Definition at line 2347 of file Alignment.cpp.
References numberOfResidues, numberOfSequences, originalNumberOfResidues, originalNumberOfSequences, saveResidues, and saveSequences.
Referenced by Cleaner::removeColumns(), Cleaner::removeOnlyTerminal(), and Cleaner::removeSequences().
Alignment::alignmentInfo |
String containing information of the alignment as a whole.
Definition at line 86 of file Alignment.h.
Referenced by Alignment().
Alignment::Cleaning = nullptr |
Trimming submodule.
It contains methods and variables related to trimming.
Definition at line 58 of file Alignment.h.
Referenced by Alignment(), Cleaner::calculateRepresentativeSeq(), Cleaner::cleanByCutValueFallBehind(), Cleaner::cleanByCutValueOverpass(), Cleaner::cleanByCutValueOverpassOrEquals(), Cleaner::cleanOverlapSeq(), trimAlManager::CleanResiduesAuto(), trimAlManager::CleanResiduesNonAuto(), trimAlManager::CleanSequences(), Cleaner::cleanStrict(), trimAlManager::innerPerform(), trimAlManager::postprocess_alignment(), printSeqIdentity(), Cleaner::removeColumns(), Cleaner::removeSequences(), setBlockSize(), setKeepSequencesFlag(), and ~Alignment().
|
private |
Flag variable that represents the type of the alignment.
Definition at line 51 of file Alignment.h.
Referenced by Alignment(), and getAlignmentType().
Alignment::filename |
Filename where this alignment was loaded from.
Definition at line 84 of file Alignment.h.
Referenced by Alignment(), alignmentSummaryHTML(), alignmentSummarySVG(), fillMatrices(), trimAlManager::output_reports(), trimAlManager::perform(), statistics::Mold::printAccumulative(), statistics::Mold::printByColumn(), statistics::Similarity::printConservationAcl(), statistics::Similarity::printConservationColumns(), statistics::Gaps::printGapsAcl(), statistics::Gaps::printGapsColumns(), statistics::Consistency::printStatisticsFileAcl(), statistics::Consistency::printStatisticsFileColumns(), FormatHandling::FormatManager::replaceINtag(), FormatHandling::FormatManager::splitAlignmentKeeping(), and statSVG().
Alignment::identities = nullptr |
2D Matrix that represents how much each sequence resembles other sequences in the same MSA
Identity score between two sequences is the ratio of identical residues by the total length (common and no-common residues) among them.
Definition at line 89 of file Alignment.h.
Referenced by Alignment(), Cleaner::calculateRelaxedSeqIdentity(), Cleaner::calculateRepresentativeSeq(), Cleaner::calculateSeqIdentity(), Cleaner::getCutPointClusters(), printSeqIdentity(), Cleaner::selectMethod(), and ~Alignment().
Alignment::isAligned = false |
Flag that indicates if all sequences on the alignment have the same length (Including gaps).
Definition at line 76 of file Alignment.h.
Referenced by Alignment(), alignmentSummaryHTML(), alignmentSummarySVG(), fillMatrices(), isFileAligned(), and main().
Alignment::numberOfResidues = 0 |
Number of residues present on the alignment if it is aligned.
Definition at line 74 of file Alignment.h.
Referenced by Alignment(), alignmentSummaryHTML(), alignmentSummarySVG(), calculateColIdentity(), Cleaner::calculateSeqIdentity(), calculateSeqOverlap(), statistics::Similarity::calculateVectors(), trimAlManager::check_select_cols_and_seqs_incompatibilities(), Cleaner::cleanByCutValueFallBehind(), Cleaner::cleanByCutValueOverpass(), Cleaner::cleanByCutValueOverpassOrEquals(), Cleaner::cleanStrict(), Cleaner::computeComplementaryAlig(), fillMatrices(), getAlignmentType(), getNumAminos(), getTranslationCDS(), printAlignmentInfo(), statistics::Consistency::printStatisticsFileAcl(), statistics::Consistency::printStatisticsFileColumns(), Cleaner::removeAllGapsSeqsAndCols(), Cleaner::removeSmallerBlocks(), FormatHandling::FormatManager::saveAlignments(), FormatHandling::FormatManager::splitAlignmentKeeping(), statSVG(), and updateSequencesAndResiduesNums().
Alignment::numberOfSequences = 0 |
Number of sequences present on the alignment.
Definition at line 70 of file Alignment.h.
Referenced by Alignment(), alignmentSummaryHTML(), alignmentSummarySVG(), statistics::Gaps::calcCutPoint(), statistics::Gaps::calcCutPoint2ndSlope(), calculateColIdentity(), calculateSeqOverlap(), trimAlManager::check_select_cols_and_seqs_incompatibilities(), checkCorrespondence(), Cleaner::computeComplementaryAlig(), fillMatrices(), getAlignmentType(), Cleaner::getClustering(), Cleaner::getCutPointClusters(), getNumSpecies(), getSequenceNameOrder(), getSequences(), getTranslationCDS(), prepareCodingSequence(), printAlignmentInfo(), printSeqIdentity(), printSeqOverlap(), Cleaner::removeAllGapsSeqsAndCols(), FormatHandling::FormatManager::saveAlignments(), Cleaner::selectMethod(), FormatHandling::FormatManager::splitAlignmentKeeping(), updateSequencesAndResiduesNums(), ~Alignment(), and statistics::Similarity::~Similarity().
Alignment::originalNumberOfResidues = 0 |
Number of residues the alignment had when it was loaded.
Definition at line 72 of file Alignment.h.
Referenced by Alignment(), alignmentSummaryHTML(), alignmentSummarySVG(), statistics::Gaps::applyWindow(), statistics::Similarity::applyWindow(), statistics::Gaps::calcCutPoint(), statistics::Similarity::calcCutPoint(), statistics::Gaps::calcCutPoint2ndSlope(), statistics::Gaps::calcCutPointMixSlope(), statistics::Similarity::calculateMatrixIdentity(), Cleaner::calculateRelaxedSeqIdentity(), Cleaner::calculateSpuriousVector(), statistics::Similarity::calculateVectors(), statistics::Gaps::CalculateVectors(), Cleaner::cleanByCutValueFallBehind(), Cleaner::cleanByCutValueOverpass(), Cleaner::cleanByCutValueOverpassOrEquals(), Cleaner::cleanCombMethods(), Cleaner::cleanCompareFile(), Cleaner::cleanStrict(), Cleaner::computeComplementaryAlig(), statistics::Gaps::Gaps(), getTranslationCDS(), statistics::Consistency::perform(), statistics::Mold::printAccumulative(), statistics::Mold::printByColumn(), statistics::Similarity::printConservationAcl(), statistics::Similarity::printConservationColumns(), statistics::Manager::printCorrespondence(), statistics::Gaps::printGapsAcl(), statistics::Gaps::printGapsColumns(), Cleaner::removeAllGapsSeqsAndCols(), Cleaner::removeOnlyTerminal(), Alignment::sequencesMatrix::sequencesMatrix(), statistics::Similarity::Similarity(), FormatHandling::FormatManager::splitAlignmentKeeping(), statSVG(), and updateSequencesAndResiduesNums().
Alignment::originalNumberOfSequences = 0 |
Number of sequences the alignment had when it was loaded.
Definition at line 68 of file Alignment.h.
Referenced by Alignment(), alignmentSummaryHTML(), alignmentSummarySVG(), statistics::Gaps::applyWindow(), statistics::Gaps::calcCutPoint(), statistics::Gaps::calcCutPointMixSlope(), statistics::Similarity::calculateMatrixIdentity(), Cleaner::calculateRelaxedSeqIdentity(), Cleaner::calculateRepresentativeSeq(), Cleaner::calculateSeqIdentity(), Cleaner::calculateSpuriousVector(), statistics::Similarity::calculateVectors(), statistics::Gaps::CalculateVectors(), trimAlManager::check_absolute_gap_theshold(), Cleaner::cleanNoAllGaps(), Cleaner::cleanOverlapSeq(), Cleaner::cleanSpuriousSeq(), Cleaner::computeComplementaryAlig(), statistics::Gaps::Gaps(), Cleaner::getClustering(), Cleaner::getCutPointClusters(), getTranslationCDS(), trimAlManager::perform(), statistics::Gaps::printGapsAcl(), statistics::Gaps::printGapsColumns(), printSeqIdentity(), Cleaner::removeAllGapsSeqsAndCols(), Cleaner::removeDuplicates(), Alignment::sequencesMatrix::sequencesMatrix(), FormatHandling::FormatManager::splitAlignmentKeeping(), statSVG(), and updateSequencesAndResiduesNums().
Alignment::overlaps = nullptr |
Matrix that stores the overlap between sequences of an alignment. See Alignment::calculateSeqOverlap.
Definition at line 92 of file Alignment.h.
Referenced by calculateSeqOverlap(), and printSeqOverlap().
Alignment::saveResidues = nullptr |
Vector containing which residues will be kept/removed.
-1 Indicates a 'removed/toRemove'. Otherwise residue/column will be kept.
Definition at line 95 of file Alignment.h.
Referenced by Alignment(), alignmentSummaryHTML(), alignmentSummarySVG(), Cleaner::calculateRelaxedSeqIdentity(), Cleaner::calculateSeqIdentity(), Cleaner::cleanByCutValueFallBehind(), Cleaner::cleanByCutValueOverpass(), Cleaner::cleanByCutValueOverpassOrEquals(), Cleaner::cleanCombMethods(), Cleaner::cleanStrict(), Cleaner::computeComplementaryAlig(), fillMatrices(), getTranslationCDS(), statistics::Manager::printCorrespondence(), Cleaner::removeAllGapsSeqsAndCols(), Cleaner::removeColumns(), Cleaner::removeOnlyTerminal(), Cleaner::removeSmallerBlocks(), updateSequencesAndResiduesNums(), and ~Alignment().
Alignment::saveSequences = nullptr |
Vector containing which sequences will be kept/removed.
-1 Indicates a 'removed/toRemove'. Otherwise sequence will be kept.
Definition at line 96 of file Alignment.h.
Referenced by Alignment(), alignmentSummaryHTML(), alignmentSummarySVG(), Cleaner::calculateRelaxedSeqIdentity(), Cleaner::calculateRepresentativeSeq(), Cleaner::calculateSeqIdentity(), statistics::Gaps::CalculateVectors(), Cleaner::cleanOverlapSeq(), Cleaner::computeComplementaryAlig(), fillMatrices(), Cleaner::getClustering(), Cleaner::getCutPointClusters(), getTranslationCDS(), Cleaner::removeAllGapsSeqsAndCols(), Cleaner::removeDuplicates(), Cleaner::removeSequences(), updateSequencesAndResiduesNums(), and ~Alignment().
Alignment::SeqRef = nullptr |
Pointer to keep a count of how many alignments depend on the same shared-sequences.
Increased on each trim step, and automatically reduced on each alignment destroy.
If it arrives to 0, the last alignment being destroyed will also take care of releasing memory for shared information.
Definition at line 66 of file Alignment.h.
Referenced by Alignment(), and ~Alignment().
Alignment::seqsInfo = nullptr |
Vector containing extra information about each sequence.
Definition at line 82 of file Alignment.h.
Referenced by Alignment(), getTranslationCDS(), and ~Alignment().
Alignment::seqsName = nullptr |
String vector containing the sequences names.
Definition at line 80 of file Alignment.h.
Referenced by Alignment(), alignmentSummaryHTML(), alignmentSummarySVG(), checkCorrespondence(), fillMatrices(), getSequenceNameOrder(), getSequences(), getTranslationCDS(), trimAlManager::perform(), prepareCodingSequence(), printAlignmentInfo(), printSeqIdentity(), printSeqOverlap(), Cleaner::removeAllGapsSeqsAndCols(), Cleaner::removeDuplicates(), Alignment::sequencesMatrix::sequencesMatrix(), FormatHandling::FormatManager::splitAlignmentKeeping(), and ~Alignment().
Alignment::sequences = nullptr |
Vector containing the sequences.
Definition at line 78 of file Alignment.h.
Referenced by Alignment(), alignmentSummaryHTML(), alignmentSummarySVG(), calculateColIdentity(), statistics::Manager::calculateGapStats(), statistics::Similarity::calculateMatrixIdentity(), Cleaner::calculateRelaxedSeqIdentity(), Cleaner::calculateRepresentativeSeq(), Cleaner::calculateSeqIdentity(), calculateSeqOverlap(), Cleaner::calculateSpuriousVector(), statistics::Similarity::calculateVectors(), statistics::Gaps::CalculateVectors(), checkCorrespondence(), fillMatrices(), getAlignmentType(), Cleaner::getCutPointClusters(), getSequences(), getTranslationCDS(), prepareCodingSequence(), printAlignmentInfo(), Cleaner::removeAllGapsSeqsAndCols(), Cleaner::removeDuplicates(), Alignment::sequencesMatrix::sequencesMatrix(), FormatHandling::FormatManager::splitAlignmentKeeping(), and ~Alignment().
Alignment::SequencesMatrix = nullptr |
SequencesMatrix submodule
This helper module stores the sequences without gaps. Useful when comparing MSA.
Definition at line 62 of file Alignment.h.
Referenced by Alignment(), statistics::Consistency::compareAndChoose(), statistics::Consistency::forceComparison(), statistics::Consistency::perform(), and ~Alignment().
Alignment::Statistics = nullptr |
Statistics submodule.
It contains methods and variables related to statistics calculation and reporting.
Definition at line 60 of file Alignment.h.
Referenced by Alignment(), alignmentSummaryHTML(), alignmentSummarySVG(), statistics::Manager::calculateConservationStats(), statistics::Similarity::calculateVectors(), Cleaner::clean(), Cleaner::clean2ndSlope(), Cleaner::cleanCombMethods(), Cleaner::cleanConservation(), Cleaner::cleanGaps(), Cleaner::cleanNoAllGaps(), trimAlManager::CleanResiduesNonAuto(), trimAlManager::CleanSequences(), trimAlManager::create_or_use_similarity_matrix(), trimAlManager::innerPerform(), statistics::Consistency::perform(), trimAlManager::print_statistics(), statistics::Manager::printStatisticsConservationColumns(), statistics::Manager::printStatisticsConservationTotal(), statistics::Manager::printStatisticsGapsColumns(), statistics::Manager::printStatisticsGapsTotal(), Cleaner::removeOnlyTerminal(), statistics::Manager::setSimilarityMatrix(), setWindowsSize(), statSVG(), and ~Alignment().