table of contents
other versions
- jessie 2.1.2-1
- jessie-backports 2.1.2+20150609-1~bpo8+1
- stretch 2.1.2+20160831-5
- testing 3.1.2+dfsg-3
- unstable 3.1.2+dfsg-5
libhmsbeagle/beagle.h(3) | BEAGLE | libhmsbeagle/beagle.h(3) |
NAME¶
libhmsbeagle/beagle.h - This file documents the API as well as header for the Broad-platform Evolutionary Analysis General Likelihood Evaluator.SYNOPSIS¶
#include 'libhmsbeagle/platform.h'Classes¶
struct BeagleInstanceDetails
Enumerations¶
enum BeagleReturnCodes { BEAGLE_SUCCESS = 0, BEAGLE_ERROR_GENERAL = -1, BEAGLE_ERROR_OUT_OF_MEMORY = -2, BEAGLE_ERROR_UNIDENTIFIED_EXCEPTION = -3, BEAGLE_ERROR_UNINITIALIZED_INSTANCE = -4, BEAGLE_ERROR_OUT_OF_RANGE = -5, BEAGLE_ERROR_NO_RESOURCE = -6, BEAGLE_ERROR_NO_IMPLEMENTATION = -7, BEAGLE_ERROR_FLOATING_POINT = -8 }
Functions¶
BEAGLE_DLLEXPORT const char * beagleGetVersion (void)
Detailed Description¶
This file documents the API as well as header for the Broad-platform Evolutionary Analysis General Likelihood Evaluator. Copyright 2009-2013 Phylogenetic Likelihood Working Group This file is part of BEAGLE. BEAGLE is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. BEAGLE is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details. You should have received a copy of the GNU Lesser General Public License along with BEAGLE. If not, see http://www.gnu.org/licenses/. KEY CONCEPTS The key to BEAGLE performance lies in delivering fine-scale parallelization while minimizing data transfer and memory copy overhead. To accomplish this, the library lacks the concept of data structure for a tree, in spite of the intended use for phylogenetic analysis. Instead, BEAGLE acts directly on flexibly indexed data storage (called buffers) for observed character states and partial likelihoods. The client program can set the input buffers to reflect the data and can calculate the likelihood of a particular phylogeny by invoking likelihood calculations on the appropriate input and output buffers in the correct order. Because of this design simplicity, the library can support many different tree inference algorithms and likelihood calculation on a variety of models. Arbitrary numbers of states can be used, as can nonreversible substitution matrices via complex eigen decompositions, and mixture models with multiple rate categories and/or multiple eigen decompositions. Finally, BEAGLE application programming interface (API) calls can be asynchronous, allowing the calling program to implement other coarse-scale parallelization schemes such as evaluating independent genes or running concurrent Markov chains. USAGE To use the library, a client program first creates an instance of BEAGLE by calling beagleCreateInstance; multiple instances per client are possible and encouraged. All additional functions are called with a reference to this instance. The client program can optionally request that an instance run on certain hardware (e.g., a GPU) or have particular features (e.g., double-precision math). Next, the client program must specify the data dimensions and specify key aspects of the phylogenetic model. Character state data are then loaded and can be in the form of discrete observed states or partial likelihoods for ambiguous characters. The observed data are usually unchanging and loaded only once at the start to minimize memory copy overhead. The character data can be compressed into unique “site patterns” and associated weights for each. The parameters of the substitution process can then be specified, including the equilibrium state frequencies, the rates for one or more substitution rate categories and their weights, and finally, the eigen decomposition for the substitution process. In order to calculate the likelihood of a particular tree, the client program then specifies a series of integration operations that correspond to steps in Felsenstein’s algorithm. Finite-time transition probabilities for each edge are loaded directly if considering a nondiagonalizable model or calculated in parallel from the eigen decomposition and edge lengths specified. This is performed within BEAGLE’s memory space to minimize data transfers. A single function call will then request one or more integration operations to calculate partial likelihoods over some or all nodes. The operations are performed in the order they are provided, typically dictated by a postorder traversal of the tree topology. The client needs only specify nodes for which the partial likelihoods need updating, but it is up to the calling software to keep track of these dependencies. The final step in evaluating the phylogenetic model is done using an API call that yields a single log likelihood for the model given the data. Aspects of the BEAGLE API design support both maximum likelihood (ML) and Bayesian phylogenetic tree inference. For ML inference, API calls can calculate first and second derivatives of the likelihood with respect to the lengths of edges (branches). In both cases, BEAGLE provides the ability to cache and reuse previously computed partial likelihood results, which can yield a tremendous speedup over recomputing the entire likelihood every time a new phylogenetic model is evaluated. Author:Likelihood API Working Group
Daniel Ayres
Peter Beerli
Michael Cummings
Aaron Darling
Mark Holder
John Huelsenbeck
Paul Lewis
Michael Ott
Andrew Rambaut
Fredrik Ronquist
Marc Suchard
David Swofford
Derrick Zwickl
Enumeration Type Documentation¶
enum BeagleFlags¶
Hardware and implementation capability flags. This enumerates all possible hardware and implementation capability flags. Each capability is a bit in a 'long' Enumerator- BEAGLE_FLAG_PRECISION_SINGLE
- Single precision computation
- BEAGLE_FLAG_PRECISION_DOUBLE
- Double precision computation
- BEAGLE_FLAG_COMPUTATION_SYNCH
- Synchronous computation (blocking)
- BEAGLE_FLAG_COMPUTATION_ASYNCH
- Asynchronous computation (non-blocking)
- BEAGLE_FLAG_EIGEN_REAL
- Real eigenvalue computation
- BEAGLE_FLAG_EIGEN_COMPLEX
- Complex eigenvalue computation
- BEAGLE_FLAG_SCALING_MANUAL
- Manual scaling
- BEAGLE_FLAG_SCALING_AUTO
- Auto-scaling on (deprecated, may not work correctly)
- BEAGLE_FLAG_SCALING_ALWAYS
- Scale at every updatePartials (deprecated, may not work correctly)
- BEAGLE_FLAG_SCALING_DYNAMIC
- Manual scaling with dynamic checking (deprecated, may not work correctly)
- BEAGLE_FLAG_SCALERS_RAW
- Save raw scalers
- BEAGLE_FLAG_SCALERS_LOG
- Save log scalers
- BEAGLE_FLAG_INVEVEC_STANDARD
- Inverse eigen vectors passed to BEAGLE have not been transposed
- BEAGLE_FLAG_INVEVEC_TRANSPOSED
- Inverse eigen vectors passed to BEAGLE have been transposed
- BEAGLE_FLAG_VECTOR_SSE
- SSE computation
- BEAGLE_FLAG_VECTOR_AVX
- AVX computation
- BEAGLE_FLAG_VECTOR_NONE
- No vector computation
- BEAGLE_FLAG_THREADING_OPENMP
- OpenMP threading
- BEAGLE_FLAG_THREADING_NONE
- No threading
- BEAGLE_FLAG_PROCESSOR_CPU
- Use CPU as main processor
- BEAGLE_FLAG_PROCESSOR_GPU
- Use GPU as main processor
- BEAGLE_FLAG_PROCESSOR_FPGA
- Use FPGA as main processor
- BEAGLE_FLAG_PROCESSOR_CELL
- Use Cell as main processor
- BEAGLE_FLAG_PROCESSOR_PHI
- Use Intel Phi as main processor
- BEAGLE_FLAG_PROCESSOR_OTHER
- Use other type of processor
- BEAGLE_FLAG_FRAMEWORK_CUDA
- Use CUDA implementation with GPU resources
- BEAGLE_FLAG_FRAMEWORK_OPENCL
- Use OpenCL implementation with GPU resources
- BEAGLE_FLAG_FRAMEWORK_CPU
- Use CPU implementation
enum BeagleOpCodes¶
Operation codes. This enumerates all possible BEAGLE operation codes. Enumerator- BEAGLE_OP_COUNT
- Total number of integers per beagleUpdatePartials operation
- BEAGLE_OP_NONE
- Specify no use for indexed buffer
enum BeagleReturnCodes¶
Error return codes. This enumerates all possible BEAGLE return codes. Error codes are always negative. Enumerator- BEAGLE_SUCCESS
- Success
- BEAGLE_ERROR_GENERAL
- Unspecified error
- BEAGLE_ERROR_OUT_OF_MEMORY
- Not enough memory could be allocated
- BEAGLE_ERROR_UNIDENTIFIED_EXCEPTION
- Unspecified exception
- BEAGLE_ERROR_UNINITIALIZED_INSTANCE
- The instance index is out of range, or the instance has not been created
- BEAGLE_ERROR_OUT_OF_RANGE
- One of the indices specified exceeded the range of the array
- BEAGLE_ERROR_NO_RESOURCE
- No resource matches requirements
- BEAGLE_ERROR_NO_IMPLEMENTATION
- No implementation matches requirements
- BEAGLE_ERROR_FLOATING_POINT
- Floating-point range exceeded
Function Documentation¶
BEAGLE_DLLEXPORT int beagleAccumulateScaleFactors (intinstance, const int *scaleIndices, intcount, intcumulativeScaleIndex)¶
Accumulate scale factors. This function adds (log) scale factors from a list of scaleBuffers to a cumulative scale buffer. It is used to calculate the marginal scaling at a specific node for each site. Parameters:instance Instance number (input)
scaleIndices List of scaleBuffers to add (input)
count Number of scaleBuffers in list (input)
cumulativeScaleIndex Index number of scaleBuffer to accumulate factors
into (input)
BEAGLE_DLLEXPORT int beagleCalculateEdgeLogLikelihoods (intinstance, const int *parentBufferIndices, const int *childBufferIndices, const int *probabilityIndices, const int *firstDerivativeIndices, const int *secondDerivativeIndices, const int *categoryWeightsIndices, const int *stateFrequenciesIndices, const int *cumulativeScaleIndices, intcount, double *outSumLogLikelihood, double *outSumFirstDerivative, double *outSumSecondDerivative)¶
Calculate site log likelihoods and derivatives along an edge. This function integrates a list of partials at a parent and child node with respect to a set of partials-weights and state frequencies to return the log likelihood and first and second derivative sums Parameters:instance Instance number (input)
parentBufferIndices List of indices of parent partialsBuffers (input)
childBufferIndices List of indices of child partialsBuffers (input)
probabilityIndices List indices of transition probability matrices for
this edge (input)
firstDerivativeIndices List indices of first derivative matrices (input)
secondDerivativeIndices List indices of second derivative matrices
(input)
categoryWeightsIndices List of weights to apply to each partialsBuffer
(input)
stateFrequenciesIndices List of state frequencies for each partialsBuffer
(input). There should be one set for each of parentBufferIndices
cumulativeScaleIndices List of scaleBuffers containing accumulated
factors to apply to each partialsBuffer (input). There should be one index for
each of parentBufferIndices
count Number of partialsBuffers (input)
outSumLogLikelihood Pointer to destination for resulting log likelihood
(output)
outSumFirstDerivative Pointer to destination for resulting first
derivative (output)
outSumSecondDerivative Pointer to destination for resulting second
derivative (output)
Returns:
error code
BEAGLE_DLLEXPORT int beagleCalculateRootLogLikelihoods (intinstance, const int *bufferIndices, const int *categoryWeightsIndices, const int *stateFrequenciesIndices, const int *cumulativeScaleIndices, intcount, double *outSumLogLikelihood)¶
Calculate site log likelihoods at a root node. This function integrates a list of partials at a node with respect to a set of partials-weights and state frequencies to return the log likelihood sum Parameters:instance Instance number (input)
bufferIndices List of partialsBuffer indices to integrate (input)
categoryWeightsIndices List of weights to apply to each partialsBuffer
(input). There should be one categoryCount sized set for each of
parentBufferIndices
stateFrequenciesIndices List of state frequencies for each partialsBuffer
(input). There should be one set for each of parentBufferIndices
cumulativeScaleIndices List of scaleBuffers containing accumulated
factors to apply to each partialsBuffer (input). There should be one index for
each of parentBufferIndices
count Number of partialsBuffer to integrate (input)
outSumLogLikelihood Pointer to destination for resulting log likelihood
(output)
Returns:
error code
BEAGLE_DLLEXPORT int beagleConvolveTransitionMatrices (intinstance, const int *firstIndices, const int *secondIndices, const int *resultIndices, intmatrixCount)¶
Convolve lists of transition probability matrices. This function convolves two lists of transition probability matrices. Parameters:instance Instance number (input)
firstIndices List of indices of the first transition probability matrices
to convolve (input)
secondIndices List of indices of the second transition probability
matrices to convolve (input)
resultIndices List of indices of resulting transition probability
matrices (input)
matrixCount Length of lists
BEAGLE_DLLEXPORT int beagleCopyScaleFactors (intinstance, intdestScalingIndex, intsrcScalingIndex)¶
Copy scale factors. This function copies scale factors from one buffer to another. Parameters:instance Instance number (input)
destScalingIndex Destination scaleBuffer (input)
srcScalingIndex Source scaleBuffer (input)
BEAGLE_DLLEXPORT int beagleCreateInstance (inttipCount, intpartialsBufferCount, intcompactBufferCount, intstateCount, intpatternCount, inteigenBufferCount, intmatrixBufferCount, intcategoryCount, intscaleBufferCount, int *resourceList, intresourceCount, longpreferenceFlags, longrequirementFlags, BeagleInstanceDetails *returnInfo)¶
Create a single instance. This function creates a single instance of the BEAGLE library and can be called multiple times to create multiple data partition instances each returning a unique identifier. Parameters:tipCount Number of tip data elements (input)
partialsBufferCount Number of partials buffers to create (input)
compactBufferCount Number of compact state representation buffers to
create (input)
stateCount Number of states in the continuous-time Markov chain (input)
patternCount Number of site patterns to be handled by the instance
(input)
eigenBufferCount Number of rate matrix eigen-decomposition, category
weight, and state frequency buffers to allocate (input)
matrixBufferCount Number of transition probability matrix buffers (input)
categoryCount Number of rate categories (input)
scaleBufferCount Number of scale buffers to create, ignored for auto
scale or always scale (input)
resourceList List of potential resources on which this instance is
allowed (input, NULL implies no restriction)
resourceCount Length of resourceList list (input)
preferenceFlags Bit-flags indicating preferred implementation
characteristics, see BeagleFlags (input)
requirementFlags Bit-flags indicating required implementation
characteristics, see BeagleFlags (input)
returnInfo Pointer to return implementation and resource details
Returns:
the unique instance identifier (<0 if failed, see
BeagleReturnCodes)
BEAGLE_DLLEXPORT int beagleFinalize (void)¶
Finalize the library. This function finalizes the library and releases all allocated memory. This function is automatically called under GNU C via attribute ((destructor)). Returns:error code
BEAGLE_DLLEXPORT int beagleFinalizeInstance (intinstance)¶
Finalize this instance. This function finalizes the instance by releasing allocated memory Parameters:instance Instance number
Returns:
error code
BEAGLE_DLLEXPORT const char* beagleGetCitation (void)¶
Get citation. This function returns a pointer to a string describing the version of the library and how to cite it. Returns:A string describing the version of the library and how to
cite it
BEAGLE_DLLEXPORT int beagleGetPartials (intinstance, intbufferIndex, intscaleIndex, double *outPartials)¶
Get partials from an instance buffer. This function copies an instance buffer into the array outPartials. The outPartials array should be stateCount * patternCount * categoryCount in length. Parameters:instance Instance number from which to get
partialsBuffer (input)
bufferIndex Index of source partialsBuffer (input)
scaleIndex Index of scaleBuffer to apply to partialsBuffer (input)
outPartials Pointer to which to receive partialsBuffer (output)
Returns:
error code
BEAGLE_DLLEXPORT BeagleResourceList* beagleGetResourceList (void)¶
Get list of hardware resources. This function returns a pointer to a BeagleResourceList struct, which includes a BeagleResource array describing the available hardware resources. Returns:A list of hardware resources available to the library as
a BeagleResourceList
BEAGLE_DLLEXPORT int beagleGetScaleFactors (intinstance, intsrcScalingIndex, double *outScaleFactors)¶
Get scale factors. This function retrieves a buffer of scale factors. Parameters:instance Instance number (input)
srcScalingIndex Source scaleBuffer (input)
outScaleFactors Pointer to which to receive scaleFactors (output)
BEAGLE_DLLEXPORT int beagleGetSiteDerivatives (intinstance, double *outFirstDerivatives, double *outSecondDerivatives)¶
Get site derivatives for last beagleCalculateEdgeLogLikelihoods call. This function returns the derivatives for each site Parameters:instance Instance number (input)
outFirstDerivatives Pointer to destination for resulting first
derivatives (output)
outSecondDerivatives Pointer to destination for resulting second
derivatives (output)
Returns:
error code
BEAGLE_DLLEXPORT int beagleGetSiteLogLikelihoods (intinstance, double *outLogLikelihoods)¶
Get site log likelihoods for last beagleCalculateRootLogLikelihoods or beagleCalculateEdgeLogLikelihoods call. This function returns the log likelihoods for each site Parameters:instance Instance number (input)
outLogLikelihoods Pointer to destination for resulting log likelihoods
(output)
Returns:
error code
BEAGLE_DLLEXPORT int beagleGetTransitionMatrix (intinstance, intmatrixIndex, double *outMatrix)¶
Get a finite-time transition probability matrix. This function copies a finite-time transition matrix buffer into the array outMatrix. The outMatrix array should be of size stateCount * stateCount * categoryCount and will be filled with one matrix for each rate category. Parameters:instance Instance number (input)
matrixIndex Index of matrix buffer (input)
outMatrix Pointer to destination transition probability matrix
(output)
Returns:
error code
BEAGLE_DLLEXPORT const char* beagleGetVersion (void)¶
Get version. This function returns a pointer to a string with the library version number. Returns:A string with the version number
BEAGLE_DLLEXPORT int beagleRemoveScaleFactors (intinstance, const int *scaleIndices, intcount, intcumulativeScaleIndex)¶
Remove scale factors. This function removes (log) scale factors from a cumulative scale buffer. The scale factors to be removed are indicated in a list of scaleBuffers. Parameters:instance Instance number (input)
scaleIndices List of scaleBuffers to remove (input)
count Number of scaleBuffers in list (input)
cumulativeScaleIndex Index number of scaleBuffer containing accumulated
factors (input)
BEAGLE_DLLEXPORT int beagleResetScaleFactors (intinstance, intcumulativeScaleIndex)¶
Reset scalefactors. This function resets a cumulative scale buffer. Parameters:instance Instance number (input)
cumulativeScaleIndex Index number of cumulative scaleBuffer (input)
BEAGLE_DLLEXPORT int beagleSetCategoryRates (intinstance, const double *inCategoryRates)¶
Set category rates. This function sets the vector of category rates for an instance. Parameters:instance Instance number (input)
inCategoryRates Array containing categoryCount rate scalers (input)
Returns:
error code
BEAGLE_DLLEXPORT int beagleSetCategoryWeights (intinstance, intcategoryWeightsIndex, const double *inCategoryWeights)¶
Set a category weights buffer. This function copies a category weights array into an instance buffer. Parameters:instance Instance number (input)
categoryWeightsIndex Index of category weights buffer (input)
inCategoryWeights Category weights array (categoryCount) (input)
Returns:
error code
BEAGLE_DLLEXPORT int beagleSetEigenDecomposition (intinstance, inteigenIndex, const double *inEigenVectors, const double *inInverseEigenVectors, const double *inEigenValues)¶
Set an eigen-decomposition buffer. This function copies an eigen-decomposition into an instance buffer. Parameters:instance Instance number (input)
eigenIndex Index of eigen-decomposition buffer (input)
inEigenVectors Flattened matrix (stateCount x stateCount) of
eigen-vectors (input)
inInverseEigenVectors Flattened matrix (stateCount x stateCount) of
inverse-eigen- vectors (input)
inEigenValues Vector of eigenvalues
Returns:
error code
BEAGLE_DLLEXPORT int beagleSetPartials (intinstance, intbufferIndex, const double *inPartials)¶
Set an instance partials buffer. This function copies an array of partials into an instance buffer. The inPartials array should be stateCount * patternCount * categoryCount in length. Parameters:instance Instance number in which to set a
partialsBuffer (input)
bufferIndex Index of destination partialsBuffer (input)
inPartials Pointer to partials values to set (input)
Returns:
error code
BEAGLE_DLLEXPORT int beagleSetPatternWeights (intinstance, const double *inPatternWeights)¶
Set pattern weights. This function sets the vector of pattern weights for an instance. Parameters:instance Instance number (input)
inPatternWeights Array containing patternCount weights (input)
Returns:
error code
BEAGLE_DLLEXPORT int beagleSetStateFrequencies (intinstance, intstateFrequenciesIndex, const double *inStateFrequencies)¶
Set a state frequency buffer. This function copies a state frequency array into an instance buffer. Parameters:instance Instance number (input)
stateFrequenciesIndex Index of state frequencies buffer (input)
inStateFrequencies State frequencies array (stateCount) (input)
Returns:
error code
BEAGLE_DLLEXPORT int beagleSetTipPartials (intinstance, inttipIndex, const double *inPartials)¶
Set an instance partials buffer for tip node. This function copies an array of partials into an instance buffer. The inPartials array should be stateCount * patternCount in length. For most applications this will be used to set the partial likelihoods for the observed states. Internally, the partials will be copied categoryCount times. Parameters:instance Instance number in which to set a
partialsBuffer (input)
tipIndex Index of destination partialsBuffer (input)
inPartials Pointer to partials values to set (input)
Returns:
error code
BEAGLE_DLLEXPORT int beagleSetTipStates (intinstance, inttipIndex, const int *inStates)¶
Set the compact state representation for tip node. This function copies a compact state representation into an instance buffer. Compact state representation is an array of states: 0 to stateCount - 1 (missing = stateCount). The inStates array should be patternCount in length (replication across categoryCount is not required). Parameters:instance Instance number (input)
tipIndex Index of destination compactBuffer (input)
inStates Pointer to compact states (input)
Returns:
error code
BEAGLE_DLLEXPORT int beagleSetTransitionMatrices (intinstance, const int *matrixIndices, const double *inMatrices, const double *paddedValues, intcount)¶
Set multiple transition matrices. This function copies multiple transition matrices into matrix buffers. This function is used when the application wishes to explicitly set the transition matrices rather than using the beagleSetEigenDecomposition and beagleUpdateTransitionMatrices functions. The inMatrices array should be of size stateCount * stateCount * categoryCount * count. Parameters:instance Instance number (input)
matrixIndices Indices of matrix buffers (input)
inMatrices Pointer to source transition matrices (input)
paddedValues Values to be used for padding for ambiguous states (e.g. 1
for probability matrices, 0 for derivative matrices) (input)
count Number of transition matrices (input)
Returns:
error code
BEAGLE_DLLEXPORT int beagleSetTransitionMatrix (intinstance, intmatrixIndex, const double *inMatrix, doublepaddedValue)¶
Set a finite-time transition probability matrix. This function copies a finite-time transition probability matrix into a matrix buffer. This function is used when the application wishes to explicitly set the transition probability matrix rather than using the beagleSetEigenDecomposition and beagleUpdateTransitionMatrices functions. The inMatrix array should be of size stateCount * stateCount * categoryCount and will contain one matrix for each rate category. Parameters:instance Instance number (input)
matrixIndex Index of matrix buffer (input)
inMatrix Pointer to source transition probability matrix (input)
paddedValue Value to be used for padding for ambiguous states (e.g. 1 for
probability matrices, 0 for derivative matrices) (input)
Returns:
error code
BEAGLE_DLLEXPORT int beagleUpdatePartials (const intinstance, const BeagleOperation *operations, intoperationCount, intcumulativeScaleIndex)¶
Calculate or queue for calculation partials using a list of operations. This function either calculates or queues for calculation a list partials. Implementations supporting ASYNCH may queue these calculations while other implementations perform these operations immediately and in order. Parameters:instance Instance number (input)
operations BeagleOperation list specifying operations (input)
operationCount Number of operations (input)
cumulativeScaleIndex Index number of scaleBuffer to store accumulated
factors (input)
Returns:
error code
BEAGLE_DLLEXPORT int beagleUpdateTransitionMatrices (intinstance, inteigenIndex, const int *probabilityIndices, const int *firstDerivativeIndices, const int *secondDerivativeIndices, const double *edgeLengths, intcount)¶
Calculate a list of transition probability matrices. This function calculates a list of transition probabilities matrices and their first and second derivatives (if requested). Parameters:instance Instance number (input)
eigenIndex Index of eigen-decomposition buffer (input)
probabilityIndices List of indices of transition probability matrices to
update (input)
firstDerivativeIndices List of indices of first derivative matrices to
update (input, NULL implies no calculation)
secondDerivativeIndices List of indices of second derivative matrices to
update (input, NULL implies no calculation)
edgeLengths List of edge lengths with which to perform calculations
(input)
count Length of lists
Returns:
error code
BEAGLE_DLLEXPORT int beagleWaitForPartials (const intinstance, const int *destinationPartials, intdestinationPartialsCount)¶
Block until all calculations that write to the specified partials have completed. This function is optional and only has to be called by clients that 'recycle' partials. If used, this function must be called after an beagleUpdatePartials call and must refer to indices of 'destinationPartials' that were used in a previous beagleUpdatePartials call. The library will block until those partials have been calculated. Parameters:instance Instance number (input)
destinationPartials List of the indices of destinationPartials that must
be calculated before the function returns
destinationPartialsCount Number of destinationPartials (input)
Returns:
error code
Author¶
Generated automatically by Doxygen for BEAGLE from the source code.Fri Jul 24 2015 | Version 2.1.2 |