Codebase list ghmm / HEAD ghmm / pmodel.h
HEAD

Tree @HEAD (Download .tar.gz)

pmodel.h @HEADraw · history · blame

/*******************************************************************************
*
*       This file is part of the General Hidden Markov Model Library,
*       GHMM version __VERSION__, see http://ghmm.org
*
*       Filename: sclass_change.c
*       Authors:  Matthias Heinig
*
*       Copyright (C) 1998-2004 Alexander Schliep
*       Copyright (C) 1998-2001 ZAIK/ZPR, Universitaet zu Koeln
*       Copyright (C) 2002-2004 Max-Planck-Institut fuer Molekulare Genetik,
*                               Berlin
*
*       Contact: schliep@ghmm.org
*
*       This library is free software; you can redistribute it and/or
*       modify it under the terms of the GNU Library General Public
*       License as published by the Free Software Foundation; either
*       version 2 of the License, or (at your option) any later version.
*
*       This library 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
*       Library General Public License for more details.
*
*       You should have received a copy of the GNU Library General Public
*       License along with this library; if not, write to the Free
*       Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
*
*
*       This file is version $Revision: 2277 $
*                       from $Date: 2009-04-28 08:44:31 -0400 (Tue, 28 Apr 2009) $
*             last change by $Author: grunau $.
*
*******************************************************************************/
#ifndef GHMM_PMODEL_H
#define GHMM_PMODEL_H

#ifdef __cplusplus
extern "C" {
#endif

#include "model.h" 
#include "psequence.h"

struct ghmm_dpmodel;

  /** Class change context for pair HMMs.
   *
   * Holds functions pointers for determining the class for switching HMMs.
   */
typedef struct ghmm_dpmodel_class_change_context{

    /* Names of class change module/function (for python callback) */
    char* python_module;
    char* python_function;
    
    /** pointer to class function called with seq X, Y and resp indices 
     in the void you can pass the user data */
    int (*get_class)(struct ghmm_dpmodel*, ghmm_dpseq*, ghmm_dpseq*, int, int,void*);
    
    /* space for any data necessary for class switch, USER is RESPONSIBLE */
    void* user_data;
} ghmm_dpmodel_class_change_context;

  /** State struct for pair HMMs.
   */
typedef struct ghmm_dpstate {
  /** Initial probability */ 
  double pi;
  /** Log of the initial probability */
  double log_pi;
  /** Output probability */
  double *b;
  int order;
  
  /** IDs of the following states */ 
  int *out_id;  
  /** IDs of the previous states */    
  int *in_id;

  /** transition probs to successor states. (for each transition class)*/
  double **out_a; 
  /** transition probs from predecessor states. (for each transition class)*/ 
  double **in_a;
  /** number of transition classes in this state **/
  int kclasses;
  /** pointer to class function   */
  ghmm_dpmodel_class_change_context *class_change; 
  /** int (*get_class)(int*,int); */
  
  /** Transition probability to a successor 
      double *out_a; */
  /** Transition probablity to a precursor 
      double *in_a;*/

  /** Number of successor states */     
  int out_states; 
  /** Number of precursor states */
  int in_states;  
  /** if fix == 1 --> b stays fix during the training */
  int fix;

  int label;
  
  /** EXTENSIONS for Pair HMMs **/
  /** read offset_x characters from sequence X **/
  int offset_x;
  /** read offset_y characters from sequence Y **/
  int offset_y;
  /** which emission alphabet **/
  int alphabet;
  /** x coordinate position for graph representation plotting **/
  int xPosition;
  /** y coordinate position for graph representation plotting **/
  int yPosition;
} ghmm_dpstate;

/** Model struct for pair HMMs.
 *
 * Contains all parameters, that define a HMM.
 */
typedef struct ghmm_dpmodel {
  /** Number of states */
  int N;
  /** Number of outputs */   
  int M;   
  /** Vector of the states */
  ghmm_dpstate *s; 
  /** The a priori probability for the model.
      A value of -1 indicates that no prior is defined. 
      Note: this is not to be confused with priors on emission
      distributions*/
  double prior;

  /* contains a arbitrary name for the model */
  char* name;
  
   /** Contains bit flags for varios model extensions such as
      kSilentStates, kTiedEmissions (see ghmm.h for a complete list)
  */
  int model_type;

  /** Flag variables for each state indicating whether it is emitting
      or not. 
      Note: silent != NULL iff (model_type & kSilentStates) == 1  */
  int* silent; /*AS*/

  /** Int variable for the maximum level of higher order emissions */
  int maxorder;
  /** saves the history of emissions as int, 
      the nth-last emission is (emission_history * |alphabet|^n+1) % |alphabet|
      see ...*/
  int emission_history;

  /** Flag variables for each state indicating whether the states emissions
      are tied to another state. Groups of tied states are represented
      by their tie group leader (the lowest numbered member of the group).
      
      tied_to[s] == kUntied  : s is not a tied state
      
      tied_to[s] == s        : s is a tie group leader

      tied_to[t] == s        : t is tied to state s

      Note: tied_to != NULL iff (model_type & kTiedEmissions) == 1  */
  int* tied_to; 
  
  /** Note: State store order information of the emissions.
      Classic HMMS have emission order 0, that is the emission probability
      is conditioned only on the state emitting the symbol.

      For higher order emissions, the emission are conditioned on the state s
      as well as the previous emission_order[s] observed symbols.

      The emissions are stored in the state's usual double* b. The order is
      set state.order.

      Note: state.order != NULL iff (model_type & kHigherOrderEmissions) == 1  */
  
  /** ghmm_dbackground is a pointer to a
      ghmm_dbackground structure, which holds (essentially) an
      array of background distributions (which are just vectors of floating
      point numbers like state.b).

      For each state the array background_id indicates which of the background
      distributions to use in parameter estimation. A value of kNoBackgroundDistribution
      indicates that none should be used.


      Note: background_id != NULL iff (model_type & kHasBackgroundDistributions) == 1  */
  int *background_id;
  ghmm_dbackground* bp;

  /** (WR) added these variables for topological ordering of silent states 
      Condition: topo_order != NULL iff (model_type & kSilentStates) == 1
   */
  int* topo_order; 
  int  topo_order_length;

  /** EXTENSIONS for Pair HMMs **/
  /** total number of alphabets **/
  int number_of_alphabets;
  /** list of sizes of the alphabets **/
  int * size_of_alphabet;
  /** number of double sequences to modify the transition classes */
  int number_of_d_seqs;
  /** maximal offset in sequence X (for the viterbi lookback matrix) **/
  int max_offset_x;
  /** maximal offset in sequence Y (for the viterbi lookback matrix) **/
  int max_offset_y;

  int debug;
} ghmm_dpmodel;

int ghmm_dpmodel_state_alloc(ghmm_dpstate *state, int M, int in_states, int out_states);

ghmm_dpmodel* ghmm_dpmodel_init(void);

ghmm_dpmodel_class_change_context* ghmm_dpmodel_init_class_change(void);

void ghmm_dpmodel_state_clean(ghmm_dpstate *my_state);

int ghmm_dpmodel_free(ghmm_dpmodel *mo);

ghmm_dpstate* get_pstateptr(ghmm_dpstate *ary, int index); 

/** functions dealing with the emission indices **/
int ghmm_dpmodel_pair(int symbol_x, int symbol_y, int alphabet_size, int off_x, int off_y);

int ghmm_dpmodel_emission_table_size(ghmm_dpmodel * mo, int state_index);
  
int ghmm_dpmodel_default_transition_class(ghmm_dpmodel * mo, ghmm_dpseq * X, ghmm_dpseq * Y, int index_x, int index_y, void * user_data) ;

void ghmm_dpmodel_set_to_default_transition_class(ghmm_dpmodel_class_change_context * pccc);


#ifdef __cplusplus
}
#endif

#endif /* GHMM_PMODEL_H */