NCBI C++ ToolKit
Classes | Typedefs | Functions
blast_psi.h File Reference

High level definitions and declarations for the PSSM engine of PSI-BLAST. More...

#include <algo/blast/core/ncbi_std.h>
#include <algo/blast/core/blast_export.h>
#include <algo/blast/core/blast_options.h>
#include <algo/blast/core/blast_stat.h>
+ Include dependency graph for blast_psi.h:
+ This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Go to the SVN repository for this file.

Classes

struct  PSIMsaCell
 Structure to describe the characteristics of a position in the multiple sequence alignment data structure. More...
 
struct  PSIMsaDimensions
 Structure representing the dimensions of the multiple sequence alignment data structure. More...
 
struct  PSIMsa
 Multiple sequence alignment (msa) data structure containing the raw data needed by the PSSM engine to create a PSSM. More...
 
struct  PSICdMsaCellData
 Data needed for PSSM computation stored in MSA cell for single column in CD aligned to a position in the query. More...
 
struct  PSICdMsaCell
 Alignment cell that represents one column of CD aligned to a position in the query. More...
 
struct  PSICdMsa
 Data structure representing multiple alignemnt of CDs and query sequence along with data needed for PSSM computation. More...
 
struct  PSIMatrix
 This is the main return value from the PSSM engine. More...
 
struct  PSIDiagnosticsRequest
 Structure to allow requesting various diagnostics data to be collected by PSSM engine. More...
 
struct  PSIDiagnosticsResponse
 This structure contains the diagnostics information requested using the PSIDiagnosticsRequest structure. More...
 

Typedefs

typedef struct PSIMsaCell PSIMsaCell
 Structure to describe the characteristics of a position in the multiple sequence alignment data structure. More...
 
typedef struct PSIMsaDimensions PSIMsaDimensions
 Structure representing the dimensions of the multiple sequence alignment data structure. More...
 
typedef struct PSIMsa PSIMsa
 Multiple sequence alignment (msa) data structure containing the raw data needed by the PSSM engine to create a PSSM. More...
 
typedef struct PSICdMsaCellData PSICdMsaCellData
 Data needed for PSSM computation stored in MSA cell for single column in CD aligned to a position in the query. More...
 
typedef struct PSICdMsaCell PSICdMsaCell
 Alignment cell that represents one column of CD aligned to a position in the query. More...
 
typedef struct PSICdMsa PSICdMsa
 Data structure representing multiple alignemnt of CDs and query sequence along with data needed for PSSM computation. More...
 
typedef struct PSIMatrix PSIMatrix
 This is the main return value from the PSSM engine. More...
 
typedef struct PSIDiagnosticsRequest PSIDiagnosticsRequest
 Structure to allow requesting various diagnostics data to be collected by PSSM engine. More...
 
typedef struct PSIDiagnosticsResponse PSIDiagnosticsResponse
 This structure contains the diagnostics information requested using the PSIDiagnosticsRequest structure. More...
 

Functions

PSIMsaPSIMsaNew (const PSIMsaDimensions *dimensions)
 Allocates and initializes the multiple sequence alignment data structure for use as input to the PSSM engine. More...
 
PSIMsaPSIMsaFree (PSIMsa *msa)
 Deallocates the PSIMsa structure. More...
 
PSIMatrixPSIMatrixNew (Uint4 query_length, Uint4 alphabet_size)
 Allocates a new PSIMatrix structure. More...
 
PSIMatrixPSIMatrixFree (PSIMatrix *matrix)
 Deallocates the PSIMatrix structure passed in. More...
 
PSIDiagnosticsRequestPSIDiagnosticsRequestNew (void)
 Allocates a PSIDiagnosticsRequest structure, setting all fields to false. More...
 
PSIDiagnosticsRequestPSIDiagnosticsRequestNewEx (Boolean save_ascii_pssm)
 Allocates a PSIDiagnosticsRequest structure, setting fields to their default values for their use in the context of the PSI-BLAST application. More...
 
PSIDiagnosticsRequestPSIDiagnosticsRequestFree (PSIDiagnosticsRequest *diags_request)
 Deallocates the PSIDiagnosticsRequest structure passed in. More...
 
PSIDiagnosticsResponsePSIDiagnosticsResponseNew (Uint4 query_length, Uint4 alphabet_size, const PSIDiagnosticsRequest *request)
 Allocates a new PSI-BLAST diagnostics structure based on which fields of the PSIDiagnosticsRequest structure are TRUE. More...
 
PSIDiagnosticsResponsePSIDiagnosticsResponseFree (PSIDiagnosticsResponse *diags)
 Deallocates the PSIDiagnosticsResponse structure passed in. More...
 
int PSICreatePssm (const PSIMsa *msap, const PSIBlastOptions *options, BlastScoreBlk *sbp, PSIMatrix **pssm)
 Main entry point to core PSSM engine to calculate the PSSM. More...
 
int PSICreatePssmWithDiagnostics (const PSIMsa *msap, const PSIBlastOptions *options, BlastScoreBlk *sbp, const PSIDiagnosticsRequest *request, PSIMatrix **pssm, PSIDiagnosticsResponse **diagnostics)
 Main entry point to core PSSM engine which allows to request diagnostics information. More...
 
int PSICreatePssmFromCDD (const PSICdMsa *cd_msa, const PSIBlastOptions *options, BlastScoreBlk *sbp, const PSIDiagnosticsRequest *request, PSIMatrix **pssm, PSIDiagnosticsResponse **diagnostics)
 Main entry point to core PSSM engine for computing CDD-based PSSMs. More...
 
int PSICreatePssmFromFrequencyRatios (const Uint1 *query, Uint4 query_length, BlastScoreBlk *sbp, double **freq_ratios, double impala_scaling_factor, PSIMatrix **pssm)
 Top-level function to create a PSSM given a matrix of frequency ratios and perform scaling on the resulting PSSM (i.e. More...
 

Detailed Description

High level definitions and declarations for the PSSM engine of PSI-BLAST.

Definition in file blast_psi.h.

Typedef Documentation

◆ PSICdMsa

typedef struct PSICdMsa PSICdMsa

Data structure representing multiple alignemnt of CDs and query sequence along with data needed for PSSM computation.

◆ PSICdMsaCell

typedef struct PSICdMsaCell PSICdMsaCell

Alignment cell that represents one column of CD aligned to a position in the query.

◆ PSICdMsaCellData

Data needed for PSSM computation stored in MSA cell for single column in CD aligned to a position in the query.

◆ PSIDiagnosticsRequest

Structure to allow requesting various diagnostics data to be collected by PSSM engine.

◆ PSIDiagnosticsResponse

This structure contains the diagnostics information requested using the PSIDiagnosticsRequest structure.

◆ PSIMatrix

typedef struct PSIMatrix PSIMatrix

This is the main return value from the PSSM engine.

◆ PSIMsa

typedef struct PSIMsa PSIMsa

Multiple sequence alignment (msa) data structure containing the raw data needed by the PSSM engine to create a PSSM.

By convention, the first row of the data field contains the query sequence

◆ PSIMsaCell

typedef struct PSIMsaCell PSIMsaCell

Structure to describe the characteristics of a position in the multiple sequence alignment data structure.

◆ PSIMsaDimensions

Structure representing the dimensions of the multiple sequence alignment data structure.

Function Documentation

◆ PSICreatePssm()

int PSICreatePssm ( const PSIMsa msap,
const PSIBlastOptions options,
BlastScoreBlk sbp,
PSIMatrix **  pssm 
)

Main entry point to core PSSM engine to calculate the PSSM.

Parameters
msapmultiple sequence alignment data structure [in]
optionsoptions to the PSSM engine [in]
sbpBLAST score block structure [in|out]
pssmPSSM and statistical information (the latter is also returned in the sbp->kbp_gap_psi[0])
Returns
PSI_SUCCESS on success, otherwise one of the PSIERR_* constants

Definition at line 95 of file blast_psi.c.

References NULL, and PSICreatePssmWithDiagnostics().

◆ PSICreatePssmFromCDD()

int PSICreatePssmFromCDD ( const PSICdMsa cd_msa,
const PSIBlastOptions options,
BlastScoreBlk sbp,
const PSIDiagnosticsRequest request,
PSIMatrix **  pssm,
PSIDiagnosticsResponse **  diagnostics 
)

Main entry point to core PSSM engine for computing CDD-based PSSMs.

Parameters
cd_msainformation about CDs that match to query sequence [in]
optionsoptions to PSSM engine [in]
sbpBLAST score block structure [in|out]
requestdiagnostics information request [in]
pssmPSSM [out]
diagnosticsdiagnostics information response, expects a pointer to uninitialized structure [in|out]

Definition at line 229 of file blast_psi.c.

References _PSIComputeFreqRatiosFromCDs(), _PSIComputeFrequenciesFromCDs(), _PSICreateAndScalePssmFromFrequencyRatios(), _PSIInternalPssmDataNew(), _PSISaveCDDiagnostics(), _PSISequenceWeightsNew(), _PSIValidateCdMSA(), BlastScoreBlk::alphabet_size, PSICdMsa::dimensions, PSIBlastOptions::impala_scaling_factor, NULL, PSIBlastOptions::pseudo_count, PSI_SUCCESS, PSIDiagnosticsResponseFree(), PSIDiagnosticsResponseNew(), PSIERR_BADPARAM, PSIERR_OUTOFMEM, PSIMatrixNew(), PSICdMsa::query, PSIMsaDimensions::query_length, s_PSICreatePssmCleanUp(), s_PSISavePssm(), and _PSISequenceWeights::std_prob.

Referenced by CPssmEngine::x_CreatePssmFromCDD().

◆ PSICreatePssmFromFrequencyRatios()

int PSICreatePssmFromFrequencyRatios ( const Uint1 query,
Uint4  query_length,
BlastScoreBlk sbp,
double **  freq_ratios,
double  impala_scaling_factor,
PSIMatrix **  pssm 
)

Top-level function to create a PSSM given a matrix of frequency ratios and perform scaling on the resulting PSSM (i.e.

: performs the last two stages of the algorithm) Note that no diagnostics can be returned as those are calculated in earlier stages of the algorithm.

Parameters
queryquery sequence in ncbistdaa format, no sentinels needed [in]
query_lengthlength of the query sequence [in]
sbpBLAST score block structure [in|out]
freq_ratiosmatrix of frequency ratios, dimensions are query_length by BLASTAA_SIZE [in]
impala_scaling_factorscaling factor used in IMPALA-style scaling if its value is NOT kPSSM_NoImpalaScaling (otherwise it performs standard PSI-BLAST scaling) [in]
pssmPSSM and statistical information [in|out]
Returns
PSI_SUCCESS on success, otherwise one of the PSIERR_* constants

Definition at line 344 of file blast_psi.c.

References _PSICopyMatrix_double(), _PSICreateAndScalePssmFromFrequencyRatios(), _PSIInternalPssmDataNew(), BlastScoreBlk::alphabet_size, BLAST_GetStandardAaProbabilities(), _PSIInternalPssmData::freq_ratios, _PSIInternalPssmData::ncols, _PSIInternalPssmData::nrows, NULL, PSI_SUCCESS, PSIERR_OUTOFMEM, PSIMatrixNew(), query, s_PSICreatePssmFromFrequencyRatiosCleanUp(), and s_PSISavePssm().

Referenced by CPssmEngine::x_CreatePssmFromFreqRatios().

◆ PSICreatePssmWithDiagnostics()

int PSICreatePssmWithDiagnostics ( const PSIMsa msap,
const PSIBlastOptions options,
BlastScoreBlk sbp,
const PSIDiagnosticsRequest request,
PSIMatrix **  pssm,
PSIDiagnosticsResponse **  diagnostics 
)

Main entry point to core PSSM engine which allows to request diagnostics information.

Parameters
msapmultiple sequence alignment data structure [in]
optionsoptions to the PSSM engine [in]
sbpBLAST score block structure [in|out]
requestdiagnostics information request [in]
pssmPSSM and statistical information (the latter is also returned in the sbp->kbp_gap_psi[0]) [out]
diagnosticsdiagnostics information response, expects a pointer to an uninitialized structure which will be populated with data requested in requests [in|out]
Returns
PSI_SUCCESS on success, otherwise one of the PSIERR_* constants

Definition at line 105 of file blast_psi.c.

References _PSIAlignedBlockNew(), _PSIComputeAlignmentBlocks(), _PSIComputeFreqRatios(), _PSIComputeSequenceWeights(), _PSICreateAndScalePssmFromFrequencyRatios(), _PSIInternalPssmDataNew(), _PSIMsaNew(), _PSIPackedMsaFree(), _PSIPackedMsaNew(), _PSIPurgeBiasedSegments(), _PSISaveDiagnostics(), _PSISequenceWeightsNew(), _PSIStructureGroupCustomization(), _PSIValidateMSA(), _PSIValidateMSA_StructureGroup(), BlastScoreBlk::alphabet_size, _PSIMsa::dimensions, PSIBlastOptions::ignore_unaligned_positions, PSIBlastOptions::impala_scaling_factor, PSIBlastOptions::nsg_compatibility_mode, NULL, PSIBlastOptions::pseudo_count, PSI_SUCCESS, PSIDiagnosticsResponseFree(), PSIDiagnosticsResponseNew(), PSIERR_BADPARAM, PSIERR_OUTOFMEM, PSIMatrixNew(), _PSIMsa::query, PSIMsaDimensions::query_length, s_PSICreatePssmCleanUp(), s_PSISavePssm(), and _PSISequenceWeights::std_prob.

Referenced by PSICreatePssm(), and CPssmEngine::x_CreatePssmFromMsa().

◆ PSIDiagnosticsRequestFree()

PSIDiagnosticsRequest* PSIDiagnosticsRequestFree ( PSIDiagnosticsRequest diags_request)

Deallocates the PSIDiagnosticsRequest structure passed in.

Parameters
diags_requeststructure to deallocate [in]
Returns
NULL

Definition at line 611 of file blast_psi.c.

References NULL, and sfree.

Referenced by CPsiBlastInputClustalW::~CPsiBlastInputClustalW().

◆ PSIDiagnosticsRequestNew()

PSIDiagnosticsRequest* PSIDiagnosticsRequestNew ( void  )

Allocates a PSIDiagnosticsRequest structure, setting all fields to false.

Returns
newly allocated structure or NULL in case of memory allocation failure

Definition at line 585 of file blast_psi.c.

References calloc().

Referenced by CPsiBlastInputClustalW::CPsiBlastInputClustalW(), PSIDiagnosticsRequestNewEx(), and CPsiBlastTestFixture::x_ComputePssmForNextIteration().

◆ PSIDiagnosticsRequestNewEx()

PSIDiagnosticsRequest* PSIDiagnosticsRequestNewEx ( Boolean  save_ascii_pssm)

Allocates a PSIDiagnosticsRequest structure, setting fields to their default values for their use in the context of the PSI-BLAST application.

Parameters
save_ascii_pssmcorresponds to the command line argument to save the PSSM in ASCII format [in]
Returns
newly allocated structure or NULL in case of memory allocation failure

Definition at line 591 of file blast_psi.c.

References PSIDiagnosticsRequest::frequency_ratios, PSIDiagnosticsRequest::gapless_column_weights, PSIDiagnosticsRequest::information_content, PSIDiagnosticsRequest::interval_sizes, NULL, PSIDiagnosticsRequest::num_matching_seqs, PSIDiagnosticsRequestNew(), PSIDiagnosticsRequest::sigma, TRUE, and PSIDiagnosticsRequest::weighted_residue_frequencies.

Referenced by BOOST_AUTO_TEST_CASE(), CPsiBlastApp::ComputePssmForNextIteration(), CDeltaBlastApp::ComputePssmForNextPsiBlastIteration(), CDeltaBlast::Run(), and CPsiBlastArgs::x_CreatePssmFromMsa().

◆ PSIDiagnosticsResponseFree()

PSIDiagnosticsResponse* PSIDiagnosticsResponseFree ( PSIDiagnosticsResponse diags)

◆ PSIDiagnosticsResponseNew()

PSIDiagnosticsResponse* PSIDiagnosticsResponseNew ( Uint4  query_length,
Uint4  alphabet_size,
const PSIDiagnosticsRequest request 
)

Allocates a new PSI-BLAST diagnostics structure based on which fields of the PSIDiagnosticsRequest structure are TRUE.

Note: this is declared here for consistency - this does not need to be called by client code of this API, it is called in the PSICreatePssm* functions to allocate the diagnostics response structure.

Parameters
query_lengthlength of the query sequence [in]
alphabet_sizelength of the alphabet [in]
requestdiagnostics to retrieve from PSSM engine [in]
Returns
pointer to allocated PSIDiagnosticsResponse or NULL if dimensions or request are NULL

Definition at line 618 of file blast_psi.c.

References _PSIAllocateMatrix(), PSIDiagnosticsResponse::alphabet_size, calloc(), PSIDiagnosticsRequest::frequency_ratios, PSIDiagnosticsResponse::frequency_ratios, PSIDiagnosticsRequest::gapless_column_weights, PSIDiagnosticsResponse::gapless_column_weights, PSIDiagnosticsRequest::independent_observations, PSIDiagnosticsResponse::independent_observations, PSIDiagnosticsRequest::information_content, PSIDiagnosticsResponse::information_content, PSIDiagnosticsRequest::interval_sizes, PSIDiagnosticsResponse::interval_sizes, NULL, PSIDiagnosticsRequest::num_matching_seqs, PSIDiagnosticsResponse::num_matching_seqs, PSIDiagnosticsResponseFree(), PSIDiagnosticsResponse::query_length, PSIDiagnosticsResponse::residue_freqs, PSIDiagnosticsRequest::residue_frequencies, PSIDiagnosticsRequest::sigma, PSIDiagnosticsResponse::sigma, PSIDiagnosticsResponse::weighted_residue_freqs, and PSIDiagnosticsRequest::weighted_residue_frequencies.

Referenced by PSICreatePssmFromCDD(), and PSICreatePssmWithDiagnostics().

◆ PSIMatrixFree()

PSIMatrix* PSIMatrixFree ( PSIMatrix matrix)

Deallocates the PSIMatrix structure passed in.

Parameters
matrixstructure to deallocate [in]
Returns
NULL

Definition at line 569 of file blast_psi.c.

References _PSIDeallocateMatrix(), PSIMatrix::ncols, NULL, PSIMatrix::pssm, and sfree.

Referenced by PSIMatrixNew(), s_PSICreatePssmCleanUp(), and s_PSICreatePssmFromFrequencyRatiosCleanUp().

◆ PSIMatrixNew()

PSIMatrix* PSIMatrixNew ( Uint4  query_length,
Uint4  alphabet_size 
)

Allocates a new PSIMatrix structure.

Parameters
query_lengthnumber of columns allocated for the PSSM [in]
alphabet_sizenumber of rows allocated for the PSSM [in]
Returns
pointer to allocated PSIMatrix structure or NULL if out of memory

Definition at line 541 of file blast_psi.c.

References _PSIAllocateMatrix(), PSIMatrix::h, PSIMatrix::kappa, PSIMatrix::lambda, malloc(), PSIMatrix::ncols, PSIMatrix::nrows, NULL, PSIMatrixFree(), PSIMatrix::pssm, PSIMatrix::ung_h, PSIMatrix::ung_kappa, and PSIMatrix::ung_lambda.

Referenced by PSICreatePssmFromCDD(), PSICreatePssmFromFrequencyRatios(), and PSICreatePssmWithDiagnostics().

◆ PSIMsaFree()

PSIMsa* PSIMsaFree ( PSIMsa msa)

◆ PSIMsaNew()

PSIMsa* PSIMsaNew ( const PSIMsaDimensions dimensions)
Modified on Fri Sep 20 14:57:10 2024 by modify_doxy.py rev. 669887