Nie możesz wybrać więcej, niż 25 tematów Tematy muszą się zaczynać od litery lub cyfry, mogą zawierać myślniki ('-') i mogą mieć do 35 znaków.
 
 
 
 

139 wiersze
4.9 KiB

  1. /* The MIT License
  2. Copyright (c) 2014-2016 Genome Research Ltd.
  3. Author: Petr Danecek <pd3@sanger.ac.uk>
  4. Permission is hereby granted, free of charge, to any person obtaining a copy
  5. of this software and associated documentation files (the "Software"), to deal
  6. in the Software without restriction, including without limitation the rights
  7. to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
  8. copies of the Software, and to permit persons to whom the Software is
  9. furnished to do so, subject to the following conditions:
  10. The above copyright notice and this permission notice shall be included in
  11. all copies or substantial portions of the Software.
  12. THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
  13. IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
  14. FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
  15. AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
  16. LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
  17. OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
  18. THE SOFTWARE.
  19. */
  20. #ifndef __HMM_H__
  21. #define __HMM_H__
  22. #define MAT(matrix,ndim,i,j) (matrix)[(ndim)*(i)+(j)] // P(i|j), that is, transition j->i
  23. typedef struct _hmm_t hmm_t;
  24. typedef void (*set_tprob_f) (hmm_t *hmm, uint32_t prev_pos, uint32_t pos, void *data, double *tprob);
  25. /**
  26. * hmm_init() - initialize HMM
  27. * @nstates: number of states
  28. * @tprob: transition probabilities matrix (nstates x nstates), for elements ordering
  29. * see the MAT macro above.
  30. * @ntprob: number of precalculated tprob matrices or 0 for constant probs, independent
  31. * of distance
  32. */
  33. hmm_t *hmm_init(int nstates, double *tprob, int ntprob);
  34. void hmm_set_tprob(hmm_t *hmm, double *tprob, int ntprob);
  35. #define HMM_VIT 1
  36. #define HMM_FWD 2
  37. #define HMM_BWD 4
  38. /**
  39. * hmm_init_states() - initial state probabilities
  40. * @probs: initial state probabilities or NULL to reset to default
  41. *
  42. * If uncalled, all states are initialized with the same likelihood
  43. */
  44. void hmm_init_states(hmm_t *hmm, double *probs);
  45. /**
  46. * hmm_snapshot() - take the model's snapshot, intended for sliding HMM
  47. * @snapshot: NULL or snapshot returned by previous hmm_snapshot() call, must be free()-ed by the caller
  48. * @pos: take the snapshot at this position
  49. *
  50. * If both restore() and snapshot() are needed, restore() must be called first.
  51. */
  52. void *hmm_snapshot(hmm_t *hmm, void *snapshot, uint32_t pos);
  53. /**
  54. * hmm_restore() - restore model's snapshot, intended for sliding HMM
  55. * @snapshot: snapshot returned by hmm_snapshot() call or NULL to reset
  56. * @isite: take the snapshot at i-th step
  57. *
  58. * If both restore() and snapshot() are needed, restore() must be called first.
  59. */
  60. void hmm_restore(hmm_t *hmm, void *snapshot);
  61. void hmm_reset(hmm_t *hmm, void *snapshot);
  62. /**
  63. * hmm_get_tprob() - return the array of transition matrices, precalculated
  64. * to ntprob positions. The first matrix is the initial tprob matrix
  65. * set by hmm_init() or hmm_set_tprob()
  66. */
  67. double *hmm_get_tprob(hmm_t *hmm);
  68. int hmm_get_nstates(hmm_t *hmm);
  69. /**
  70. * hmm_set_tprob_func() - custom setter of transition probabilities
  71. */
  72. void hmm_set_tprob_func(hmm_t *hmm, set_tprob_f set_tprob, void *data);
  73. /**
  74. * hmm_run_viterbi() - run Viterbi algorithm
  75. * @nsites: number of sites
  76. * @eprob: emission probabilities for each site and state (nsites x nstates)
  77. * @sites: list of positions
  78. *
  79. * When done, hmm->vpath[] contains the calculated Viterbi path. The states
  80. * are indexed starting from 0, a state at i-th site can be accessed as
  81. * vpath[nstates*i].
  82. */
  83. void hmm_run_viterbi(hmm_t *hmm, int nsites, double *eprob, uint32_t *sites);
  84. /**
  85. * hmm_get_viterbi_path() - the viterbi path: state at ith site is the
  86. * (nstates*isite)-th element
  87. */
  88. uint8_t *hmm_get_viterbi_path(hmm_t *hmm);
  89. /**
  90. * hmm_run_fwd_bwd() - run the forward-backward algorithm
  91. * @nsites: number of sites
  92. * @eprob: emission probabilities for each site and state (nsites x nstates)
  93. * @sites: list of positions
  94. */
  95. void hmm_run_fwd_bwd(hmm_t *hmm, int nsites, double *eprob, uint32_t *sites);
  96. /**
  97. * hmm_get_fwd_bwd_prob() - the probability of i-th state at j-th site can
  98. * be accessed as fwd_bwd[j*nstates+i].
  99. */
  100. double *hmm_get_fwd_bwd_prob(hmm_t *hmm);
  101. /**
  102. * hmm_run_baum_welch() - run one iteration of Baum-Welch algorithm
  103. * @nsites: number of sites
  104. * @eprob: emission probabilities for each site and state (nsites x nstates)
  105. * @sites: list of positions
  106. *
  107. * Same as hmm_run_fwd_bwd, in addition a pointer to a matrix with the new
  108. * transition probabilities is returned. In this verison, emission
  109. * probabilities are not updated.
  110. */
  111. double *hmm_run_baum_welch(hmm_t *hmm, int nsites, double *eprob, uint32_t *sites);
  112. void hmm_destroy(hmm_t *hmm);
  113. #endif