Legofit
infers population history from nucleotide site patterns.
Data Structures | Functions | Variables
gptree.c File Reference

Methods for simulating gene genealogies within a given tree of populations, and allowing populations to mix and also to split. More...

#include "gene.h"
#include "gptree.h"
#include "lblndx.h"
#include "parse.h"
#include "param.h"
#include "parstore.h"
#include "ptrptrmap.h"
#include "sampndx.h"
#include <errno.h>
#include <gsl/gsl_rng.h>
#include <string.h>
#include <unistd.h>
#include <pthread.h>

Data Structures

struct  GPTree
 GPTree stands for Gene-Population tree. More...
 

Functions

int GPTree_equals (const GPTree *lhs, const GPTree *rhs)
 Return 1 if two GPTree objects are equal, 0 if they differ. More...
 
void GPTree_print (void *vself, FILE *fp)
 
void GPTree_initStateVec (void *gpt, int ndx, int n, double x[n], gsl_rng *rng)
 Initialize vector x. More...
 
void GPTree_printParStore (void *vself, FILE *fp)
 Print a description of parameters.
 
void GPTree_printParStoreFree (void *vself, FILE *fp)
 Print a description of free parameters.
 
const char * GPTree_getNameFree (void *self, int i)
 Return pointer to name of i'th free parameter.
 
void GPTree_randomize (void *vself, gsl_rng *rng)
 Randomly perturb all free parameters while maintaining inequality constraints.
 
int GPTree_setParams (void *vself, int n, double x[n])
 Set free parameters from an array. More...
 
void GPTree_getParams (void *self, int n, double x[n])
 Copy free parameters from GPTree into an array. More...
 
int GPTree_nFree (const void *self)
 Return number of free parameters.
 
void GPTree_brlen (void *vself, BranchTab *branchtab, gsl_rng *rng, unsigned long nreps, int doSing, long unsigned *event_counter)
 Use coalescent simulation to estimate the mean branch length associated with each site pattern. More...
 
void * GPTree_new (const char *fname, Bounds bnd)
 GPTree constructor.
 
void GPTree_free (void *vself)
 GPTree destructor.
 
void * GPTree_dup (const void *vold)
 Duplicate a GPTree object.
 
void GPTree_sanityCheck (void *vself, const char *file, int line)
 
LblNdx GPTree_getLblNdx (void *vself)
 Get the LblNdx object from a GPTree.
 
unsigned GPTree_nsamples (const void *vself)
 
int GPTree_feasible (const void *vself, int verbose)
 Are parameters within the feasible region?
 
unsigned GPTree_nSamples (void *vself)
 

Variables

pthread_mutex_t outputLock
 

Detailed Description

Methods for simulating gene genealogies within a given tree of populations, and allowing populations to mix and also to split.

Author
Alan R. Rogers

Function Documentation

◆ GPTree_brlen()

void GPTree_brlen ( void *  vself,
BranchTab branchtab,
gsl_rng *  rng,
unsigned long  nreps,
int  doSing,
long unsigned *  event_counter 
)

Use coalescent simulation to estimate the mean branch length associated with each site pattern.

Parameters
selfGPTree object
[out]branchtabBranchTab object, which will tabulate branch lengths from this simulations.
[in,out]rngGSL random number generator
[in]nrepsnumber of replicate gene trees to simulate
[in]doSingif doSing is non-zero, singleton site patterns will be tabulated.

References PopNode_clear().

◆ GPTree_equals()

int GPTree_equals ( const GPTree lhs,
const GPTree rhs 
)

Return 1 if two GPTree objects are equal, 0 if they differ.

Abort with an error if the GPTree pointers are different but one or more of the internal pointers is the same. Does not access rootPop or rootGene. Used for testing GPTree_dup.

◆ GPTree_getParams()

void GPTree_getParams ( void *  self,
int  n,
double  x[n] 
)

Copy free parameters from GPTree into an array.

Parameters
[out]nnumber of parameters in array, which should equal the number of free parameters in the GPTree.
[out]xarray into which parameters will be copied

References ParStore_getFreeParams().

Referenced by GPTree_initStateVec(), and Network_init().

◆ GPTree_initStateVec()

void GPTree_initStateVec ( void *  gpt,
int  ndx,
int  n,
double  x[n],
gsl_rng *  rng 
)

Initialize vector x.

If ndx==0, simply copy the parameter vector from the GPTree object. Otherwise, randomize the GPTree first. This ensures that differential evolution starts with a set of points, one of which is the same as the values in the input file. This allows you to improve on existing estimates without starting from scratch each time.

References GPTree_dup(), GPTree_free(), GPTree_getParams(), and GPTree_randomize().

◆ GPTree_setParams()

int GPTree_setParams ( void *  vself,
int  n,
double  x[n] 
)

Set free parameters from an array.

Parameters
[in]nnumber of parameters in array, which should equal the number of free parameters in the GPTree.
[in]xarray of parameter values.
Returns
0 on success, 1 if values of x violate boundary constraints.

References ParStore_nFree(), and ParStore_setFreeParams().