de.jstacs.sequenceScores.statisticalModels.trainable.hmm.states.emissions.discrete
Class PhyloDiscreteEmission

java.lang.Object
  extended by de.jstacs.sequenceScores.statisticalModels.trainable.hmm.states.emissions.discrete.AbstractConditionalDiscreteEmission
      extended by de.jstacs.sequenceScores.statisticalModels.trainable.hmm.states.emissions.discrete.DiscreteEmission
          extended by de.jstacs.sequenceScores.statisticalModels.trainable.hmm.states.emissions.discrete.PhyloDiscreteEmission
All Implemented Interfaces:
SamplingComponent, SamplingFromStatistic, DifferentiableEmission, Emission, SamplingEmission, Storable, Cloneable

public class PhyloDiscreteEmission
extends DiscreteEmission
implements SamplingEmission

Phylogenetic discrete emission This class uses a phylogenetic tree to describe multidimensional data It implements Felsensteins model for nucleotide substitution (F81)

Author:
Michael Scharfe

Field Summary
 
Fields inherited from class de.jstacs.sequenceScores.statisticalModels.trainable.hmm.states.emissions.discrete.AbstractConditionalDiscreteEmission
con, counter, ess, grad, hyperParams, logNorm, offset, params, paramsFile, probs, reader, samplingIndex, statistic, writer
 
Constructor Summary
PhyloDiscreteEmission(AlphabetContainer con, double[] hyperParams, PhyloTree t)
          This is a simple constructor for a DiscreteEmission defining the individual hyper parameters.
PhyloDiscreteEmission(AlphabetContainer con, double ess, PhyloTree t)
          This is a simple constructor for a PhyloDiscreteEmission based on the equivalent sample size.
 
Method Summary
 void addGradientOfLogPriorTerm(double[] gradient, int offset)
          This method computes the gradient of Emission.getLogPriorTerm() for each parameter of this model.
 void addToStatistic(boolean forward, int startPos, int endPos, double weight, Sequence seq)
          This method adds the weight to the internal sufficient statistic.
protected  void appendFurtherInformation(StringBuffer xml)
          This method appends further information to the XML representation.
 PhyloDiscreteEmission clone()
           
 void estimateFromStatistic()
          This method estimates the parameters from the internal sufficient statistic.
protected  void extractFurtherInformation(StringBuffer xml)
          This method extracts further information from the XML representation.
 double getLogGammaScoreFromStatistic()
          This method calculates a score for the current statistics, which is independent from the current parameters In general the gamma-score is a product of gamma-functions parameterized with the current statistics
 double getLogPosteriorFromStatistic()
          This method calculates the a-posteriori probability for the current statistics
 double getLogProbAndPartialDerivationFor(boolean forward, int startPos, int endPos, IntList indices, DoubleList partDer, Sequence seq)
          Returns the logarithmic score for a Sequence beginning at position start in the Sequence and fills lists with the indices and the partial derivations.
 double getLogProbFor(boolean forward, int startPos, int endPos, Sequence seq)
          This method computes the logarithm of the likelihood.
 double getLogProposalPosteriorFromStatistic()
          Returns the log posterior of the proposal distribution for the current statistic
 void resetStatistic()
          This method resets the internal sufficient statistic.
 
Methods inherited from class de.jstacs.sequenceScores.statisticalModels.trainable.hmm.states.emissions.discrete.DiscreteEmission
getConditionIndex, toString
 
Methods inherited from class de.jstacs.sequenceScores.statisticalModels.trainable.hmm.states.emissions.discrete.AbstractConditionalDiscreteEmission
acceptParameters, drawParameters, drawParametersFromStatistic, extendSampling, fillCurrentParameter, fillSamplingGroups, finalize, fromXML, getAlphabetContainer, getHyperParams, getLogPriorTerm, getNodeLabel, getNodeShape, getNumberOfParameters, getSizeOfEventSpace, initForSampling, initializeFunctionRandomly, isInSamplingMode, joinStatistics, parseNextParameterSet, parseParameterSet, precompute, samplingStopped, setLinear, setParameter, setParameterOffset, setParameters, setShape, toXML
 
Methods inherited from class java.lang.Object
equals, getClass, hashCode, notify, notifyAll, toString, wait, wait, wait
 
Methods inherited from interface de.jstacs.sampling.SamplingFromStatistic
drawParametersFromStatistic
 
Methods inherited from interface de.jstacs.sampling.SamplingComponent
acceptParameters, extendSampling, initForSampling, isInSamplingMode, parseNextParameterSet, parseParameterSet, samplingStopped
 
Methods inherited from interface de.jstacs.sequenceScores.statisticalModels.trainable.hmm.states.emissions.Emission
getAlphabetContainer, getLogPriorTerm, getNodeLabel, getNodeShape, initializeFunctionRandomly, joinStatistics, setParameters, toString
 
Methods inherited from interface de.jstacs.Storable
toXML
 

Constructor Detail

PhyloDiscreteEmission

public PhyloDiscreteEmission(AlphabetContainer con,
                             double ess,
                             PhyloTree t)
This is a simple constructor for a PhyloDiscreteEmission based on the equivalent sample size.

Parameters:
con - the AlphabetContainer of this emission
ess - the equivalent sample size (ess) of this emission that is equally distributed over all parameters
t - the phylogenetic tree PhyloTree

PhyloDiscreteEmission

public PhyloDiscreteEmission(AlphabetContainer con,
                             double[] hyperParams,
                             PhyloTree t)
This is a simple constructor for a DiscreteEmission defining the individual hyper parameters.

Parameters:
con - the AlphabetContainer of this emission
hyperParams - the individual hyper parameters for each parameter
t - the phylogenetic tree PhyloTree
Method Detail

clone

public PhyloDiscreteEmission clone()
                            throws CloneNotSupportedException
Overrides:
clone in class AbstractConditionalDiscreteEmission
Throws:
CloneNotSupportedException

getLogProbFor

public double getLogProbFor(boolean forward,
                            int startPos,
                            int endPos,
                            Sequence seq)
                     throws OperationNotSupportedException
Description copied from interface: Emission
This method computes the logarithm of the likelihood.

Specified by:
getLogProbFor in interface Emission
Overrides:
getLogProbFor in class AbstractConditionalDiscreteEmission
Parameters:
forward - whether to use the forward or the reverse strand
startPos - the start position
endPos - the end position
seq - the sequence
Returns:
the logarithm of the probability
Throws:
OperationNotSupportedException - if forward=false and the reverse complement of the sequence seq is not defined

getLogProbAndPartialDerivationFor

public double getLogProbAndPartialDerivationFor(boolean forward,
                                                int startPos,
                                                int endPos,
                                                IntList indices,
                                                DoubleList partDer,
                                                Sequence seq)
                                         throws OperationNotSupportedException
Description copied from interface: DifferentiableEmission
Returns the logarithmic score for a Sequence beginning at position start in the Sequence and fills lists with the indices and the partial derivations.

Specified by:
getLogProbAndPartialDerivationFor in interface DifferentiableEmission
Overrides:
getLogProbAndPartialDerivationFor in class AbstractConditionalDiscreteEmission
Parameters:
forward - a switch whether to use the forward or the reverse complementary strand of the sequence
startPos - the start position in the Sequence
endPos - the end position in the Sequence
indices - an IntList of indices, after method invocation the list should contain the indices i where $\frac{\partial \log score(seq)}{\partial \lambda_i}$ is not zero
partDer - a DoubleList of partial derivations, after method invocation the list should contain the corresponding $\frac{\partial \log score(seq)}{\partial \lambda_i}$ that are not zero
seq - the Sequence
Returns:
the logarithmic score for the Sequence
Throws:
OperationNotSupportedException - if forward==false and the reverse complement of the sequence can not be computed

addGradientOfLogPriorTerm

public void addGradientOfLogPriorTerm(double[] gradient,
                                      int offset)
Description copied from interface: DifferentiableEmission
This method computes the gradient of Emission.getLogPriorTerm() for each parameter of this model. The results are added to the array grad beginning at index (offset + internal offset).

Specified by:
addGradientOfLogPriorTerm in interface DifferentiableEmission
Overrides:
addGradientOfLogPriorTerm in class AbstractConditionalDiscreteEmission
Parameters:
gradient - the array of gradients
offset - the start index of the HMM in the grad array, where the partial derivations for the parameters of the HMM shall be entered
See Also:
Emission.getLogPriorTerm(), DifferentiableEmission.setParameterOffset(int)

estimateFromStatistic

public void estimateFromStatistic()
Description copied from interface: Emission
This method estimates the parameters from the internal sufficient statistic.

Specified by:
estimateFromStatistic in interface Emission
Overrides:
estimateFromStatistic in class AbstractConditionalDiscreteEmission

appendFurtherInformation

protected void appendFurtherInformation(StringBuffer xml)
Description copied from class: AbstractConditionalDiscreteEmission
This method appends further information to the XML representation. It allows subclasses to save further parameters that are not defined in the superclass.

Overrides:
appendFurtherInformation in class AbstractConditionalDiscreteEmission
Parameters:
xml - the XML representation

extractFurtherInformation

protected void extractFurtherInformation(StringBuffer xml)
                                  throws NonParsableException
Description copied from class: AbstractConditionalDiscreteEmission
This method extracts further information from the XML representation. It allows subclasses to cast further parameters that are not defined in the superclass.

Overrides:
extractFurtherInformation in class AbstractConditionalDiscreteEmission
Parameters:
xml - the XML representation
Throws:
NonParsableException - if the information could not be reconstructed out of the StringBuffer xml

getLogProposalPosteriorFromStatistic

public double getLogProposalPosteriorFromStatistic()
Returns the log posterior of the proposal distribution for the current statistic

Returns:
the log posterior

getLogPosteriorFromStatistic

public double getLogPosteriorFromStatistic()
Description copied from interface: SamplingFromStatistic
This method calculates the a-posteriori probability for the current statistics

Specified by:
getLogPosteriorFromStatistic in interface SamplingFromStatistic
Overrides:
getLogPosteriorFromStatistic in class AbstractConditionalDiscreteEmission
Returns:
the logarithm of the a-posteriori probability

resetStatistic

public void resetStatistic()
Description copied from interface: Emission
This method resets the internal sufficient statistic.

Specified by:
resetStatistic in interface Emission
Overrides:
resetStatistic in class AbstractConditionalDiscreteEmission

addToStatistic

public void addToStatistic(boolean forward,
                           int startPos,
                           int endPos,
                           double weight,
                           Sequence seq)
                    throws OperationNotSupportedException
Description copied from interface: Emission
This method adds the weight to the internal sufficient statistic.

Specified by:
addToStatistic in interface Emission
Overrides:
addToStatistic in class AbstractConditionalDiscreteEmission
Parameters:
forward - whether to use the forward or the reverse strand
startPos - the start position
endPos - the end position
weight - the weight of the sequence
seq - the sequence
Throws:
OperationNotSupportedException - if forward=false and the reverse complement of the sequence seq is not defined

getLogGammaScoreFromStatistic

public double getLogGammaScoreFromStatistic()
Description copied from interface: SamplingEmission
This method calculates a score for the current statistics, which is independent from the current parameters In general the gamma-score is a product of gamma-functions parameterized with the current statistics

Specified by:
getLogGammaScoreFromStatistic in interface SamplingEmission
Overrides:
getLogGammaScoreFromStatistic in class AbstractConditionalDiscreteEmission
Returns:
the logarithm of the gamma-score for the current statistics